TABLE OF CONTENTS


ABINIT/integrate_gamma_tr [ Functions ]

[ Top ] [ Functions ]

NAME

 integrate_gamma_tr

FUNCTION

 This routine integrates the TRANSPORT electron phonon coupling matrices
 over the kpoints on the fermi surface. A dependency on qpoint
 remains for gamma_qpt_in/out
 Copied from integrate_gamma

COPYRIGHT

 Copyright (C) 2004-2018 ABINIT group (BXu,MJV)
 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 = elphon datastructure with data and dimensions
      elph_ds%qpt_full = qpoint coordinates
   FSfullpqtofull = mapping of k+q to k
   veloc_sq1 = mean square electronic velocity on constant energy surface
   veloc_sq2 = mean square electronic velocity on constant energy surface

OUTPUT

   elph_tr_ds%gamma_qpt_tr and created elph_tr_ds%gamma_rpt_tr

PARENTS

      elphon

CHILDREN

      wrtout,xmpi_sum

SOURCE

 38 #if defined HAVE_CONFIG_H
 39 #include "config.h"
 40 #endif
 41 
 42 #include "abi_common.h"
 43 
 44 
 45 subroutine integrate_gamma_tr(elph_ds,FSfullpqtofull,s1,s2, &
 46 &                             veloc_sq1,veloc_sq2,elph_tr_ds)
 47 
 48  use defs_basis
 49  use defs_elphon
 50  use m_errors
 51  use m_profiling_abi
 52  use m_xmpi
 53 
 54 !This section has been created automatically by the script Abilint (TD).
 55 !Do not modify the following lines by hand.
 56 #undef ABI_FUNC
 57 #define ABI_FUNC 'integrate_gamma_tr'
 58  use interfaces_14_hidewrite
 59 !End of the abilint section
 60 
 61  implicit none
 62 
 63 !Arguments ------------------------------------
 64 !scalars
 65  integer,intent(in) :: s1,s2
 66  type(elph_tr_type), intent(inout) :: elph_tr_ds
 67  type(elph_type),intent(in) :: elph_ds
 68 !arrays
 69  integer,intent(in) :: FSfullpqtofull(elph_ds%k_phon%nkpt,elph_ds%nqpt_full)
 70  real(dp),intent(in) :: veloc_sq1(3,elph_ds%nsppol), veloc_sq2(3,elph_ds%nsppol)
 71 
 72 !Local variables-------------------------------
 73 !scalars
 74  integer :: ikpt_phon,ikpt_phonq,ib1,ib2,ibeff,ierr,iqpt,iqpt_fullbz,isppol
 75  integer :: itensor, icomp, jcomp,comm
 76  integer :: fib1, fib2
 77  integer :: ik_this_proc
 78 ! integer :: ikpttemp
 79  character(len=500) :: message
 80  real(dp) :: wtk, wtkpq, interm
 81  real(dp) :: veloc1_i, veloc1_j, veloc2_i, veloc2_j
 82 !arrays
 83  real(dp) :: elvelock(3), elvelockpq(3)
 84  real(dp) :: velocwtk(3), velocwtkpq(3)
 85  real(dp) :: vvelocwtk(3,3), vvelocwtkpq(3,3)
 86  real(dp),allocatable :: tmp_gkk(:,:,:,:)
 87 
 88 ! *************************************************************************
 89 
 90  comm = xmpi_world
 91 
 92 !information
 93  if (elph_ds%gkqwrite == 0) then
 94    write (message,'(a)')' integrate_gamma_tr : keeping gamma matrices in memory'
 95    call wrtout(std_out,message,'COLL')
 96  else if (elph_ds%gkqwrite == 1) then
 97    write (message,'(a)')' integrate_gamma_tr : reading gamma matrices from disk'
 98    call wrtout(std_out,message,'COLL')
 99  else
100    write (message,'(3a,i3)')' integrate_gamma_tr : BUG-',ch10,&
101 &   ' Wrong value for gkqwrite = ',elph_ds%gkqwrite
102    MSG_BUG(message)
103  end if
104 
105 !allocate temp variables
106  ABI_STAT_ALLOCATE(tmp_gkk,(2,elph_ds%ngkkband**2,elph_ds%nbranch**2,elph_ds%nsppol), ierr)
107  ABI_CHECK(ierr==0, 'trying to allocate array tmp_gkkout')
108 
109  do iqpt=1,elph_ds%nqptirred
110    iqpt_fullbz = elph_ds%qirredtofull(iqpt)
111 !  write(std_out,*)'iqpt, iqptfullbz  ',iqpt, iqpt_fullbz
112 
113    do ik_this_proc =1,elph_ds%k_phon%my_nkpt
114      ikpt_phon = elph_ds%k_phon%my_ikpt(ik_this_proc)
115 
116      if (elph_ds%gkqwrite == 0) then
117        tmp_gkk = elph_ds%gkk_qpt(:,:,:,ik_this_proc,:,iqpt)
118      else if (elph_ds%gkqwrite == 1) then
119        read(elph_ds%unitgkq,REC=((iqpt-1)*elph_ds%k_phon%my_nkpt+ik_this_proc)) tmp_gkk
120      end if
121 
122      ikpt_phonq = FSfullpqtofull(ikpt_phon,iqpt_fullbz)
123 
124      do isppol=1,elph_ds%nsppol
125        do ib1=1,elph_ds%ngkkband !FS bands
126          fib1=ib1+elph_ds%minFSband-1 ! full bands
127          elvelock(:)=elph_tr_ds%el_veloc(ikpt_phon,fib1,:,isppol)
128          wtk=elph_tr_ds%tmp_gkk_intweight1(ib1,ikpt_phon,isppol)
129          velocwtk(:)=elph_tr_ds%tmp_velocwtk1(ib1,ikpt_phon,:,isppol)
130          vvelocwtk(:,:)=elph_tr_ds%tmp_vvelocwtk1(ib1,ikpt_phon,:,:,isppol)
131 
132          do ib2=1,elph_ds%ngkkband ! FS bands
133            ibeff=ib2+(ib1-1)*elph_ds%ngkkband ! full bands
134            fib2=ib2+elph_ds%minFSband-1
135            elvelockpq(:)= elph_tr_ds%el_veloc(ikpt_phonq,fib2,:,isppol)
136            wtkpq=elph_tr_ds%tmp_gkk_intweight2(ib2,ikpt_phonq,isppol)
137            velocwtkpq(:)=elph_tr_ds%tmp_velocwtk2(ib2,ikpt_phonq,:,isppol)
138            vvelocwtkpq(:,:)=elph_tr_ds%tmp_vvelocwtk2(ib2,ikpt_phonq,:,:,isppol)
139 
140 !          MJV 31/03/2009: Note that the following is valid for any geometry, not just cubic!
141 !          see eq 5 and 6 of prb 36 4103 (Al-Lehaibi et al 1987)
142 !          see also Allen PRB 17 3725
143 !          generalization to tensorial quantities is simple, by keeping the directional
144 !          references of velock and velockpq as indices.
145            do icomp = 1, 3
146              do jcomp = 1, 3
147                itensor = (icomp-1)*3+jcomp
148 !              FIXME: could use symmetry i <-> j
149 
150                veloc1_i = sqrt(veloc_sq1(icomp,isppol))
151                veloc1_j = sqrt(veloc_sq1(jcomp,isppol))
152                veloc2_i = sqrt(veloc_sq2(icomp,isppol))
153                veloc2_j = sqrt(veloc_sq2(jcomp,isppol))
154                if (elph_ds%use_k_fine == 1) then
155                  interm = vvelocwtk(icomp,jcomp)*wtkpq/veloc1_i/veloc1_j + &
156 &                 s1*s2*vvelocwtkpq(icomp,jcomp)*wtk/veloc2_i/veloc2_j - &
157 &                 s1*velocwtk(jcomp)*velocwtkpq(icomp)/veloc1_j/veloc2_i - &
158 &                 s2*velocwtk(icomp)*velocwtkpq(jcomp)/veloc1_i/veloc2_j
159 
160                  elph_tr_ds%gamma_qpt_tr(:,itensor,:,isppol,iqpt_fullbz) = &
161 &                 elph_tr_ds%gamma_qpt_tr(:,itensor,:,isppol,iqpt_fullbz) + &
162 &                 tmp_gkk(:,ibeff,:,isppol)*interm
163                else
164                  elph_tr_ds%gamma_qpt_tr(:,itensor,:,isppol,iqpt_fullbz) = &
165 &                 elph_tr_ds%gamma_qpt_tr(:,itensor,:,isppol,iqpt_fullbz) + &
166 &                 tmp_gkk(:,ibeff,:,isppol) &
167 &                 *(elvelock(icomp)/veloc1_i - s1*elvelockpq(icomp)/veloc2_i) &
168 &                 *(elvelock(jcomp)/veloc1_j - s2*elvelockpq(jcomp)/veloc2_j) &
169 &                 *wtk*wtkpq
170                end if
171              end do
172            end do
173 
174          end do
175        end do
176      end do ! isppol
177 
178    end do ! ik
179  end do ! iq
180 
181  call xmpi_sum (elph_tr_ds%gamma_qpt_tr, comm, ierr)
182 
183  ABI_DEALLOCATE(tmp_gkk)
184 
185 
186 !need prefactor of 1/nkpt for each integration over 1 kpoint index.
187 !NOT INCLUDED IN elph_ds%gkk_intweight
188 !Add a factor of 1/2 for the cross terms of (v-v')(v-v')
189  elph_tr_ds%gamma_qpt_tr = elph_tr_ds%gamma_qpt_tr* elph_ds%occ_factor*0.5_dp / elph_ds%k_phon%nkpt
190 
191  write (message,'(2a)')' integrate_gamma_tr : transport gamma matrices are calculated ',&
192 & ' in recip space and for irred qpoints'
193 !call wrtout(std_out,message,'COLL')
194 
195 end subroutine integrate_gamma_tr