TABLE OF CONTENTS


ABINIT/get_fs_bands [ Functions ]

[ Top ] [ Functions ]

NAME

 get_fs_bands

FUNCTION

 This routine determines the bands which contribute to the Fermi surface

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

  eigenGS = ground state eigenvalues
  hdr = header from input GS file
  ep_b_min, ep_b_max=A non-zero value is used to impose certain bands.
  fermie=Fermi level.
  eigenGS(hdr%nband(1),hdr%nkpt,hdr%nsppol)=Energies.

OUTPUT

  minFSband,maxFSband=Minimun and maximum index for the bands that cross the Fermi level
  nkptirr=Number of irreducible points for which there exist at least one band that crosses the Fermi level.

TODO

  1) Indeces and dimensions should should be spin dependent.
  2) In the present status of the code, all the k-points in the IBZ are used!

PARENTS

      elphon

CHILDREN

      wrtout

SOURCE

 40 #if defined HAVE_CONFIG_H
 41 #include "config.h"
 42 #endif
 43 
 44 #include "abi_common.h"
 45 
 46 
 47 subroutine get_fs_bands(eigenGS,hdr,fermie,ep_b_min,ep_b_max,minFSband,maxFSband,nkptirr)
 48 
 49  use defs_basis
 50  use defs_abitypes
 51  use m_profiling_abi
 52 
 53 !This section has been created automatically by the script Abilint (TD).
 54 !Do not modify the following lines by hand.
 55 #undef ABI_FUNC
 56 #define ABI_FUNC 'get_fs_bands'
 57  use interfaces_14_hidewrite
 58 !End of the abilint section
 59 
 60  implicit none
 61 
 62 !Arguments ------------------------------------
 63 !scalars
 64  integer, intent(in) :: ep_b_min, ep_b_max
 65  integer,intent(out) :: minFSband,maxFSband,nkptirr
 66  real(dp),intent(in) :: fermie
 67  type(hdr_type),intent(in) :: hdr
 68 !arrays
 69  real(dp),intent(in) :: eigenGS(hdr%nband(1),hdr%nkpt,hdr%nsppol)
 70 
 71 !Local variables-------------------------------
 72 !scalars
 73  integer :: iband,ikpt,isppol,nband
 74  real(dp) :: epsFS,gausstol,gaussig
 75  character(len=500) :: message
 76  integer :: kpt_phonflag(hdr%nkpt)
 77 
 78 ! *************************************************************************
 79 
 80 !supposes nband is equal for all kpts
 81  nband = hdr%nband(1)
 82 
 83 !gausstol = minimum weight value for integration weights on FS
 84 !should be set to reproduce DOS at Ef (Ref. PRB 34, 5065 p. 5067)
 85  gausstol = 1.0d-10
 86 
 87 !use same band indices in both spin channels
 88  maxFSband=1
 89  minFSband=nband
 90 
 91 !window of states around fermi Energy is contained in +/- epsFS
 92 !should be adjusted to take into account a minimal but sufficient
 93 !fraction of the kpoints: see the loop below.
 94 !The 1000 is purely empirical!!!
 95 !Should also take into account the density of kpoints.
 96 !gaussig = width of gaussian energy window around fermi energy
 97 !needed to get a good fraction of kpoints contributing to the FS
 98  
 99  gaussig = (maxval(eigenGS)-minval(eigenGS))/1000.0_dp
100  
101  write (message,'(a,f11.8,2a)')' get_fs_bands : initial energy window = ',gaussig,ch10,&
102 & ' The window energy will be increased until the full k-grid is inside the range'
103  call wrtout(std_out,message,'COLL')
104 
105 !NOTE: could loop back to here and change gaussig until we have
106 !a certain fraction of the kpoints in the FS region...
107  nkptirr = 0
108  
109 !Do not use restricted fermi surface: include all kpts -> one
110  do while (nkptirr < hdr%nkpt)
111    gaussig = gaussig*1.05_dp
112 
113 !  we must take into account kpoints with states within epsFS:
114    epsFS = gaussig*sqrt(log(one/(gaussig*sqrt(pi)*gausstol)))
115    
116 !  check if there are eigenvalues close to the Fermi surface
117 !  (less than epsFS from it)
118    kpt_phonflag(:) = 0
119 
120 !  do for each sppol channel
121    do isppol=1,hdr%nsppol
122      do ikpt=1,hdr%nkpt
123        do iband=1,nband
124          if (abs(eigenGS(iband,ikpt,isppol) - fermie) < epsFS) then
125            kpt_phonflag(ikpt) = 1
126            if (iband > maxFSband) maxFSband = iband
127            if (iband < minFSband) minFSband = iband
128          end if
129        end do
130      end do
131    end do ! isppol
132    
133 !  if user imposed certain bands for e-p, make sure they are kept
134    if (ep_b_min /= 0 .and. ep_b_min < minFSband) then
135      minFSband = ep_b_min
136    end if
137    if (ep_b_max /= 0 .and. ep_b_max > maxFSband) then
138      maxFSband = ep_b_max
139    end if
140    
141 !  number of irreducible kpoints (by all sym) contributing to the Fermi surface (to be completed by symops).
142    nkptirr = sum(kpt_phonflag(:))
143  end do
144 
145  write(std_out,*) ' Energy window around Fermi level= ',epsFS,' nkptirr= ',nkptirr
146 
147 end subroutine get_fs_bands