TABLE OF CONTENTS
ABINIT/get_nv_fs_temp [ Functions ]
NAME
get_nv_fs_temp
FUNCTION
This routine calculates the fermi energy, FD smeared DOS(Ef) and Veloc_sq(Ef) at looped temperatures.
COPYRIGHT
Copyright (C) 2012-2018 ABINIT group (BXU) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
INPUTS
elph_ds elph_ds%nband = number of bands in ABINIT elph_ds%k_fine%nkptirr = Number of irreducible points for which there exist at least one band that crosses the Fermi level. elph_ds%nFSband = number of bands included in the FS integration elph_ds%k_fine%nkpt = number of k points for fine k-grid elph_ds%k_phon%nkpt = number of k points for coarse k-grid elph_ds%tempermin = minimum temperature at which resistivity etc are calculated (in K) elph_ds%temperinc = interval temperature grid on which resistivity etc are calculated (in K) elph_ds%ep_b_min= first band taken into account in FS integration (if telphint==2) elph_ds%ep_b_max= last band taken into account in FS integration (if telphint==2) elph_ds%telphint = flag for integration over the FS with 0=tetrahedra 1=gaussians elph_ds%elphsmear = smearing width for gaussian integration or buffer in energy for calculations with tetrahedra (telphint=0) eigenGS = Ground State eigenvalues gprimd = reciprocal lattice vectors (dimensionful) kptrlatt_fine = k-point grid vectors (if divided by determinant of present matrix) max_occ = maximal occupancy for a band
OUTPUT
elph_ds%fermie=Fermi level at input temperature elph_tr_ds%dos_n0=DOS(Ef) at looped temperatures elph_tr_ds%veloc_sq0=FS averaged velocity at Ef at looped temperatures
PARENTS
elphon
CHILDREN
ebands_update_occ,ep_fs_weights,get_veloc_tr,wrtout
SOURCE
49 #if defined HAVE_CONFIG_H 50 #include "config.h" 51 #endif 52 53 #include "abi_common.h" 54 55 subroutine get_nv_fs_temp(elph_ds,BSt,eigenGS,gprimd,max_occ,elph_tr_ds) 56 57 58 use defs_basis 59 use defs_datatypes 60 use defs_elphon 61 use m_profiling_abi 62 use m_io_tools 63 64 use m_ebands, only : ebands_update_occ 65 66 !This section has been created automatically by the script Abilint (TD). 67 !Do not modify the following lines by hand. 68 #undef ABI_FUNC 69 #define ABI_FUNC 'get_nv_fs_temp' 70 use interfaces_14_hidewrite 71 use interfaces_77_ddb, except_this_one => get_nv_fs_temp 72 !End of the abilint section 73 74 implicit none 75 76 !Arguments ------------------------------------ 77 78 !data_type 79 type(elph_type),intent(inout) :: elph_ds 80 type(ebands_t),intent(inout) :: BSt 81 type(elph_tr_type),intent(inout) :: elph_tr_ds 82 83 !Scalars 84 real(dp), intent(in) :: max_occ 85 86 ! arrays 87 real(dp), intent(in) :: gprimd(3,3) 88 real(dp), intent(in) :: eigenGS(elph_ds%nband,elph_ds%k_fine%nkptirr,elph_ds%nsppol) 89 90 !Local variables------------------------------- 91 92 integer :: isppol!, ie1 93 integer :: itemp, tmp_nenergy 94 95 character(len=500) :: message 96 97 real(dp) :: Temp, tmp_elphsmear, tmp_delta_e 98 ! real(dp) :: xtr, e1 99 ! real(dp),allocatable :: tmp_wtk(:,:) 100 101 ! ************************************************************************* 102 103 ABI_ALLOCATE(elph_tr_ds%dos_n0,(elph_ds%ntemper,elph_ds%nsppol)) 104 ABI_ALLOCATE(elph_tr_ds%veloc_sq0,(elph_ds%ntemper,3,elph_ds%nsppol)) 105 !if (elph_ds%use_k_fine == 1) then 106 !ABI_ALLOCATE(tmp_wtk,(elph_ds%nFSband,elph_ds%k_fine%nkpt)) 107 !else 108 !ABI_ALLOCATE(tmp_wtk,(elph_ds%nFSband,elph_ds%k_phon%nkpt)) 109 !end if 110 111 elph_tr_ds%dos_n0 = zero 112 elph_tr_ds%veloc_sq0 = zero 113 114 tmp_nenergy = 8 115 do itemp=1,elph_ds%ntemper ! runs over temperature in K 116 Temp=elph_ds%tempermin + elph_ds%temperinc*dble(itemp) 117 tmp_delta_e = kb_HaK*Temp 118 Bst%occopt = 3 119 Bst%tsmear = Temp*kb_HaK 120 tmp_elphsmear = Temp*kb_HaK 121 call ebands_update_occ(Bst,-99.99_dp) 122 write(message,'(a,f12.6,a,E20.12)')'At T=',Temp,' Fermi level is:',Bst%fermie 123 call wrtout(std_out,message,'COLL') 124 if (abs(elph_ds%fermie) < tol10) then 125 elph_ds%fermie = BSt%fermie 126 end if 127 128 ! FD smeared DOS and veloc 129 130 call ep_fs_weights(elph_ds%ep_b_min, elph_ds%ep_b_max, eigenGS, tmp_elphsmear, & 131 & elph_ds%fermie, gprimd, elph_ds%k_fine%irredtoGS, elph_ds%kptrlatt_fine,& 132 & max_occ, elph_ds%minFSband, elph_ds%nband, elph_ds%nFSband, & 133 & elph_ds%nsppol, elph_ds%telphint, elph_ds%k_fine) 134 135 do isppol=1,elph_ds%nsppol 136 elph_ds%n0(isppol) = sum(elph_ds%k_fine%wtk(:,:,isppol))/elph_ds%k_fine%nkpt 137 write(message,'(a,f12.6,a,f12.6)')'At T=',Temp,' The DOS at Ef is:', elph_ds%n0(isppol) 138 call wrtout(std_out,message,'COLL') 139 140 ! For the non-LOVA case, N(Ef) is not that important (canceled out eventually). 141 ! Should not be important for metal, comment out for now 142 ! tmp_wtk = zero 143 ! do ie1=-tmp_nenergy,tmp_nenergy ! use ie1 here, hope there is no confusion 144 ! e1=Bst%fermie+ie1*tmp_delta_e 145 ! xtr=(e1-Bst%fermie)/(2.0_dp*kb_HaK*Temp) 146 ! 147 ! call ep_fs_weights(elph_ds%ep_b_min, elph_ds%ep_b_max, eigenGS, elph_ds%elphsmear, & 148 ! & e1, gprimd, elph_ds%k_fine%irredtoGS, elph_ds%kptrlatt_fine, & 149 ! & max_occ, elph_ds%minFSband, elph_ds%nband, elph_ds%nFSband, & 150 ! & elph_ds%nsppol, elph_ds%telphint, elph_ds%k_fine) 151 ! 152 ! tmp_wtk(:,:) = tmp_wtk(:,:) + elph_ds%k_fine%wtk(:,:,isppol)* & 153 ! & tmp_delta_e/(4.0d0*kb_HaK*Temp)/(COSH(xtr)**2.0d0) 154 ! end do ! ie1 155 156 ! elph_ds%k_fine%wtk(:,:,isppol) = tmp_wtk(:,:) 157 elph_tr_ds%dos_n0(itemp,isppol) = sum(elph_ds%k_fine%wtk(:,:,isppol))/elph_ds%k_fine%nkpt 158 ! elph_ds%n0(isppol) = sum(elph_ds%k_fine%wtk(:,:,isppol))/elph_ds%k_fine%nkpt ! for get_veloc_tr 159 ! write(message,'(a,f12.6,a,f12.6)')'At T=',Temp,' The eff. DOS at Ef is:', elph_tr_ds%dos_n0(itemp,isppol) 160 ! call wrtout(std_out,message,'COLL') 161 end do ! isppol 162 call get_veloc_tr(elph_ds,elph_tr_ds) 163 elph_tr_ds%veloc_sq0(itemp,:,:) = elph_tr_ds%FSelecveloc_sq(:,:) 164 165 end do ! temperature 166 167 end subroutine get_nv_fs_temp