TABLE OF CONTENTS


ABINIT/fxgkkphase [ Functions ]

[ Top ] [ Functions ]

NAME

 fxgkkphase

FUNCTION

   Set phase factors to eliminate gauge variable in gkk matrix elements
    (comes from phase factors in wavefunctions)

COPYRIGHT

 Copyright (C) 2004-2018 ABINIT group (MVer)
 This file is distributed under the terms of the
 GNU General Public Licence, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

   elph_ds = datastructure for elph data (dimensions and eventually data)
   gkk_flag = flags for presence of gkk matrix elements
   h1_mat_el = irreducible matrix elements to be completed
   iqptfull = qpoint number in full zone

OUTPUT

   h1_mat_el = changed on output

PARENTS

CHILDREN

SOURCE

 33 #if defined HAVE_CONFIG_H
 34 #include "config.h"
 35 #endif
 36 
 37 #include "abi_common.h"
 38 
 39 
 40 subroutine fxgkkphase(elph_ds,gkk_flag,h1_mat_el,iqptfull)
 41 
 42  use defs_basis
 43  use defs_elphon
 44  use m_profiling_abi
 45 
 46 !This section has been created automatically by the script Abilint (TD).
 47 !Do not modify the following lines by hand.
 48 #undef ABI_FUNC
 49 #define ABI_FUNC 'fxgkkphase'
 50 !End of the abilint section
 51 
 52  implicit none
 53 
 54 !Arguments ------------------------------------
 55 !scalars
 56  integer,intent(in) :: iqptfull
 57  type(elph_type),intent(in) :: elph_ds
 58 !arrays
 59  integer,intent(in) :: gkk_flag(elph_ds%nbranch,elph_ds%k_phon%nkpt,elph_ds%nqpt_full)
 60  real(dp),intent(inout) :: h1_mat_el(2,elph_ds%nFSband,elph_ds%nFSband,elph_ds%nbranch,elph_ds%k_phon%nkpt)
 61 
 62 !Local variables-------------------------------
 63 !scalars
 64  integer :: ikpt_phon,ib1,ib2,ibranch
 65  real(dp) :: imphase,norm,rephase,tmpim
 66  real(dp) :: tmpre
 67 !arrays
 68 
 69 ! *************************************************************************
 70 
 71  do ikpt_phon=1,elph_ds%k_phon%nkpt
 72    do ib1=1,elph_ds%nFSband
 73      do ib2=1,elph_ds%nFSband
 74 !      Determine phase for first perturbation
 75        rephase =  h1_mat_el(1,ib2,ib1,1,ikpt_phon)
 76        imphase = -h1_mat_el(2,ib2,ib1,1,ikpt_phon)
 77        norm = sqrt(rephase**2+imphase**2)
 78        rephase = rephase / norm
 79        imphase = imphase / norm
 80 
 81 !      DEBUG
 82 !      if (ikpt_phon == 1) then
 83 !      write(std_out,*) 'fxgkkphase : rephase,imphase = ', rephase,imphase
 84 !      end if
 85 !      ENDDEBUG
 86 
 87 !      apply same phase factor to all perturbations
 88 !      ----------------------------------------------------------
 89 !      Very important ! Otherwise the scalar product with the
 90 !      displacement vector will not be preserved.
 91 !      ----------------------------------------------------------
 92        do ibranch=1,elph_ds%nbranch
 93 !        if we already have data
 94          if (gkk_flag(ibranch,1,iqptfull) /= -1) then
 95            tmpre =    rephase * h1_mat_el(1,ib2,ib1,ibranch,ikpt_phon)&
 96 &           -imphase * h1_mat_el(2,ib2,ib1,ibranch,ikpt_phon)
 97            tmpim =    rephase * h1_mat_el(2,ib2,ib1,ibranch,ikpt_phon)&
 98 &           +imphase * h1_mat_el(1,ib2,ib1,ibranch,ikpt_phon)
 99            h1_mat_el(1,ib2,ib1,ibranch,ikpt_phon) = tmpre
100            h1_mat_el(2,ib2,ib1,ibranch,ikpt_phon) = tmpim
101          end if
102 !        end if
103        end do
104      end do
105    end do
106  end do
107 !end loop over FS kpt
108 
109 end subroutine fxgkkphase