TABLE OF CONTENTS


ABINIT/ep_el_weights [ Functions ]

[ Top ] [ Functions ]

NAME

 ep_el_weights

FUNCTION

 This routine calculates the Fermi Surface integration weights
  for the electron phonon routines, by different methods
    1) Gaussian smearing
    2) Tetrahedron method
    3) Window in bands for all k-points
    4) Fermi Dirac smearing, follows gaussian with a different smearing function

COPYRIGHT

 Copyright (C) 2010-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

   ep_b_min = minimal band to include in FS window integration
   ep_b_max = maximal band to include in FS window integration
   eigenGS = Ground State eigenvalues
   elphsmear = smearing width for Gaussian method
   fermie = Fermi level
   gprimd = Reciprocal lattice vectors (dimensionful)
   irredtoGS = mapping of elph k-points to ground state grid
   kptrlatt = k-point grid vectors (if divided by determinant of present matrix)
   max_occ = maximal occupancy for a band
   minFSband = minimal band included for Fermi Surface integration in Gaussian and Tetrahedron cases
   nFSband = number of bands in FS integration
   nsppol = number of spin polarizations
   telphint = option for FS integration:
      0 Tetrahedron method
      1 Gaussian smearing
      2 Window in bands for all k-points
      3 Fermi Dirac smearing
   k_obj%nkpt = number of FS k-points
   k_obj%kpt = FS k-points
   k_obj%full2irr = mapping of FS k-points from full grid to irred points
   k_obj%full2full = mapping of FS k-points in full grid under symops

OUTPUT

TODO

   weights should be recalculated on-the-fly! The present implementation is not flexible!

PARENTS

      get_nv_fs_en,get_tau_k

CHILDREN

      destroy_tetra,get_tetra_weight,init_tetra,matr3inv

SOURCE

 58 #if defined HAVE_CONFIG_H
 59 #include "config.h"
 60 #endif
 61 
 62 #include "abi_common.h"
 63 
 64 subroutine ep_el_weights(ep_b_min, ep_b_max, eigenGS, elphsmear, enemin, enemax, nene, gprimd, &
 65 &    irredtoGS, kptrlatt, max_occ, minFSband, nband, nFSband, nsppol, telphint, k_obj, tmp_wtk)
 66 
 67  use defs_basis
 68  use defs_elphon
 69  use m_profiling_abi
 70  use m_errors
 71  use m_tetrahedron
 72  use m_xmpi
 73 
 74 !This section has been created automatically by the script Abilint (TD).
 75 !Do not modify the following lines by hand.
 76 #undef ABI_FUNC
 77 #define ABI_FUNC 'ep_el_weights'
 78  use interfaces_32_util
 79 !End of the abilint section
 80 
 81  implicit none
 82 
 83 !Arguments ------------------------------------
 84 !scalars
 85  type(elph_kgrid_type), intent(in) :: k_obj
 86  integer, intent(in) :: ep_b_min
 87  integer, intent(in) :: ep_b_max
 88  integer,intent(in) :: nband,nene
 89  real(dp), intent(in) :: elphsmear
 90  real(dp), intent(in) :: enemin,enemax
 91  real(dp), intent(in) :: gprimd(3,3)
 92  integer, intent(in) :: kptrlatt(3,3)
 93  real(dp), intent(in) :: max_occ
 94  integer, intent(in) :: minFSband
 95  integer, intent(in) :: nFSband
 96  integer, intent(in) :: nsppol
 97  integer, intent(in) :: telphint
 98 
 99 ! arrays
100  real(dp), intent(in) :: eigenGS(nband,k_obj%nkptirr,nsppol)
101  real(dp), intent(out) :: tmp_wtk(nFSband,k_obj%nkpt,nsppol,nene)
102  integer, intent(in) :: irredtoGS(k_obj%nkptirr)
103 
104 !Local variables-------------------------------
105 !scalars
106  integer,parameter :: bcorr0=0
107  integer :: ikpt, ikptgs, ib1, iband
108  integer :: ierr, ie, isppol
109  real(dp) :: deltaene, rcvol, fermie
110  real(dp) :: smdeltaprefactor, smdeltafactor, xx
111 
112 ! arrays
113  real(dp) :: rlatt(3,3), klatt(3,3)
114  real(dp), allocatable :: tmp_eigen(:), tweight(:,:), dtweightde(:,:)
115  character (len=500) :: message
116  character (len=80) :: errstr
117  type(t_tetrahedron) :: tetrahedra
118 
119 ! *************************************************************************
120 
121  ! Initialize tmp_wtk with zeros
122  tmp_wtk = zero
123 
124  !write(std_out,*) 'ep_el : nkpt ', k_obj%nkpt
125 !===================================
126 !Set up integration weights for FS
127 !===================================
128  deltaene = (enemax-enemin)/dble(nene-1)
129 
130  if (telphint == 0) then
131 
132 !  =========================
133 !  Tetrahedron integration
134 !  =========================
135 
136    rlatt(:,:) = kptrlatt(:,:)
137    call matr3inv(rlatt,klatt)
138 
139    call init_tetra(k_obj%full2full(1,1,:), gprimd,klatt,k_obj%kpt, k_obj%nkpt,&
140 &   tetrahedra, ierr, errstr)
141    ABI_CHECK(ierr==0,errstr)
142 
143    rcvol = abs (gprimd(1,1)*(gprimd(2,2)*gprimd(3,3)-gprimd(3,2)*gprimd(2,3)) &
144 &   -gprimd(2,1)*(gprimd(1,2)*gprimd(3,3)-gprimd(3,2)*gprimd(1,3)) &
145 &   +gprimd(3,1)*(gprimd(1,2)*gprimd(2,3)-gprimd(2,2)*gprimd(1,3)))
146 
147 !  fix small window around fermie for tetrahedron weight calculation
148    deltaene = (enemax-enemin)/dble(nene-1)
149 
150    ABI_ALLOCATE(tmp_eigen,(k_obj%nkpt))
151    ABI_ALLOCATE(tweight,(k_obj%nkpt,nene))
152    ABI_ALLOCATE(dtweightde,(k_obj%nkpt,nene))
153 
154    do iband = 1,nFSband
155 !    for each spin pol
156      do isppol=1,nsppol
157 !    For this band get its contribution
158        tmp_eigen(:) = zero
159        do ikpt=1,k_obj%nkpt
160          ikptgs = irredtoGS(k_obj%full2irr(1,ikpt))
161          tmp_eigen(ikpt) = eigenGS(minFSband+iband-1,ikptgs,isppol)
162        end do
163 !      calculate general integration weights at each irred kpoint as in Blochl et al PRB 49 16223
164        call get_tetra_weight(tmp_eigen,enemin,enemax,&
165 &       max_occ,nene,k_obj%nkpt,tetrahedra,bcorr0,&
166 &       tweight,dtweightde,xmpi_comm_self)
167 
168        tmp_wtk(iband,:,isppol,:) = dtweightde(:,:)*k_obj%nkpt
169      end do
170    end do
171    ABI_DEALLOCATE(tmp_eigen)
172    ABI_DEALLOCATE(tweight)
173    ABI_DEALLOCATE(dtweightde)
174 
175    call destroy_tetra(tetrahedra)
176 
177  else if (telphint == 1) then
178 
179 !  ==============================================================
180 !  Gaussian or integration:
181 !  Each kpt contributes a gaussian of integrated weight 1 
182 !  for each band. The gaussian being centered at the Fermi level.
183 !  ===============================================================
184 
185 !  took out factor 1/k_obj%nkpt which intervenes only at integration time
186 
187 !  MJV 18/5/2008 does smdeltaprefactor need to contain max_occ?
188 
189 !  gaussian smdeltaprefactor = sqrt(piinv)/elphsmear/k_obj%nkpt
190    smdeltaprefactor = max_occ*sqrt(piinv)/elphsmear
191    smdeltafactor = one/elphsmear
192 
193 !  SPPOL loop on isppol as well to get 2 sets of weights
194    do isppol=1,nsppol
195      fermie = enemin
196      do ie = 1, nene
197        fermie = fermie + deltaene
198 !      fine grid
199        do ikpt=1, k_obj%nkpt
200          ikptgs = irredtoGS(k_obj%full2irr(1,ikpt))
201          do ib1=1,nFSband
202            xx = smdeltafactor*(eigenGS(minFSband-1+ib1,ikptgs,isppol) - fermie)
203            if (abs(xx) < 40._dp) then
204              tmp_wtk(ib1,ikpt,isppol,ie) = exp(-xx*xx)*smdeltaprefactor
205            end if
206          end do
207        end do
208      end do
209    end do
210    
211 
212  else if (telphint == 2) then ! range of bands occupied
213 
214 !  SPPOL eventually be able to specify bands for up and down separately
215    fermie = enemin
216    do ie = 1, nene
217      fermie = fermie + deltaene
218      do ikpt=1,k_obj%nkpt
219        do ib1=ep_b_min, ep_b_max
220 !        for the moment both spin channels same
221          tmp_wtk(ib1,ikpt,:,ie) = max_occ
222        end do
223      end do
224    end do
225    
226    write(std_out,*) ' ep_el_weights : DOS is calculated from states in bands ',ep_b_min,' to ',ep_b_max
227    
228  else if (telphint == 3) then
229 
230 !  ==============================================================
231 !  Fermi Dirac integration:
232 !  Each kpt contributes a Fermi Dirac smearing function of integrated weight 1 
233 !  for each band. The function being centered at the Fermi level.
234 !  ===============================================================
235 
236 !  took out factor 1/k_obj%nkpt which intervenes only at integration time
237 
238 !  MJV 18/5/2008 does smdeltaprefactor need to contain max_occ?
239 
240 !  gaussian smdeltaprefactor = sqrt(piinv)/elphsmear/k_obj%nkpt
241    smdeltaprefactor = half*max_occ/elphsmear
242    smdeltafactor = one/elphsmear
243 
244 !  SPPOL loop on isppol as well to get 2 sets of weights
245    do isppol=1,nsppol
246      fermie = enemin
247      do ie = 1, nene
248        fermie = fermie + deltaene
249 !      fine grid
250        do ikpt=1, k_obj%nkpt
251          ikptgs = irredtoGS(k_obj%full2irr(1,ikpt))
252          do ib1=1,nFSband
253            xx = smdeltafactor*(eigenGS(minFSband-1+ib1,ikptgs,isppol) - fermie)
254            tmp_wtk(ib1,ikpt,isppol,ie) = smdeltaprefactor / (one + cosh(xx))
255          end do
256        end do
257      end do
258    end do
259    
260 
261  else 
262    write (message,'(a,i0)')" telphint should be between 0 and 3, found: ",telphint
263    MSG_BUG(message)
264  end if ! if telphint
265 
266 
267 end subroutine ep_el_weights