TABLE OF CONTENTS


ABINIT/ep_fs_weights [ Functions ]

[ Top ] [ Functions ]

NAME

 ep_fs_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

   k_obj%wtk = integration weights

TODO

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

PARENTS

      elphon,get_nv_fs_en,get_nv_fs_temp

CHILDREN

      destroy_tetra,get_tetra_weight,init_tetra,matr3inv,wrtout

SOURCE

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