TABLE OF CONTENTS
ABINIT/ep_ph_weights [ Functions ]
NAME
ep_ph_weights
FUNCTION
This routine calculates the phonon integration weights for the electron phonon routines, by different methods 1) Gaussian smearing 0) Tetrahedron method
COPYRIGHT
Copyright (C) 2010-2018 ABINIT group (MVer,BXu) 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
phfrq = phonon energies elphsmear = smearing width for Gaussian method omega = input phonon energy gprimd = Reciprocal lattice vectors (dimensionful) kptrlatt = k-point grid vectors (if divided by determinant of present matrix) telphint = option for FS integration: 0 Tetrahedron method 1 Gaussian smearing k_obj%nkpt = number of FS k-points k_obj%kpt = FS k-points k_obj%full2full = mapping of FS k-points in full grid under symops
OUTPUT
tmp_wtq = integration weights
TODO
weights should be recalculated on-the-fly! The present implementation is not flexible!
PARENTS
get_tau_k,mka2f,mka2f_tr
CHILDREN
destroy_tetra,get_tetra_weight,init_tetra,matr3inv
SOURCE
47 #if defined HAVE_CONFIG_H 48 #include "config.h" 49 #endif 50 51 #include "abi_common.h" 52 53 subroutine ep_ph_weights(phfrq,elphsmear,omega_min,omega_max,nomega,gprimd,kptrlatt,nbranch,telphint,k_obj,tmp_wtq) 54 55 use defs_basis 56 use defs_elphon 57 use m_profiling_abi 58 use m_errors 59 use m_tetrahedron 60 use m_xmpi 61 62 !This section has been created automatically by the script Abilint (TD). 63 !Do not modify the following lines by hand. 64 #undef ABI_FUNC 65 #define ABI_FUNC 'ep_ph_weights' 66 use interfaces_32_util 67 !End of the abilint section 68 69 implicit none 70 71 !Arguments ------------------------------------ 72 !scalars 73 type(elph_kgrid_type), intent(inout) :: k_obj 74 integer,intent(in) :: nbranch 75 real(dp), intent(in) :: elphsmear 76 real(dp), intent(in) :: omega_min,omega_max 77 real(dp), intent(in) :: gprimd(3,3) 78 integer, intent(in) :: kptrlatt(3,3) 79 integer, intent(in) :: nomega 80 integer, intent(in) :: telphint 81 82 ! arrays 83 real(dp), intent(in) :: phfrq(nbranch,k_obj%nkpt) 84 real(dp), intent(out) :: tmp_wtq(nbranch,k_obj%nkpt,nomega) 85 86 !Local variables------------------------------- 87 !scalars 88 integer,parameter :: bcorr0=0 89 integer :: ikpt, ib1, ibranch 90 integer :: ierr, iomega 91 real(dp) :: rcvol, max_occ 92 real(dp) :: smdeltaprefactor, smdeltafactor, xx, gaussmaxarg 93 real(dp) :: domega,omega 94 95 ! arrays 96 real(dp) :: rlatt(3,3), klatt(3,3) 97 real(dp), allocatable :: tweight(:,:), dtweightde(:,:) 98 character (len=80) :: errstr 99 type(t_tetrahedron) :: tetrahedra 100 101 ! ************************************************************************* 102 103 !write(std_out,*) 'ep_ph : nqpt ', k_obj%nkpt 104 !=================================== 105 !Set up integration weights for FS 106 !=================================== 107 max_occ = one 108 gaussmaxarg = sqrt(-log(1.d-100)) 109 domega = (omega_max - omega_min)/(nomega - 1) 110 111 if (telphint == 0) then 112 113 ! ========================= 114 ! Tetrahedron integration 115 ! ========================= 116 117 rlatt(:,:) = kptrlatt(:,:) 118 call matr3inv(rlatt,klatt) 119 120 call init_tetra(k_obj%full2full(1,1,:), gprimd,klatt,k_obj%kpt, k_obj%nkpt,& 121 & tetrahedra, ierr, errstr) 122 ABI_CHECK(ierr==0,errstr) 123 124 rcvol = abs (gprimd(1,1)*(gprimd(2,2)*gprimd(3,3)-gprimd(3,2)*gprimd(2,3)) & 125 & -gprimd(2,1)*(gprimd(1,2)*gprimd(3,3)-gprimd(3,2)*gprimd(1,3)) & 126 & +gprimd(3,1)*(gprimd(1,2)*gprimd(2,3)-gprimd(2,2)*gprimd(1,3))) 127 128 ! do all the omega points for tetrahedron weight calculation 129 130 ABI_ALLOCATE(tweight,(k_obj%nkpt,nomega)) 131 ABI_ALLOCATE(dtweightde,(k_obj%nkpt,nomega)) 132 133 do ibranch = 1,nbranch 134 call get_tetra_weight(phfrq(ibranch,:),omega_min,omega_max,& 135 & max_occ,nomega,k_obj%nkpt,tetrahedra,bcorr0,& 136 & tweight,dtweightde,xmpi_comm_self) 137 138 tmp_wtq(ibranch,:,:) = dtweightde(:,:)*k_obj%nkpt 139 end do 140 ABI_DEALLOCATE(tweight) 141 ABI_DEALLOCATE(dtweightde) 142 143 call destroy_tetra(tetrahedra) 144 145 else if (telphint == 1) then 146 147 ! ============================================================== 148 ! Gaussian or integration: 149 ! Each kpt contributes a gaussian of integrated weight 1 150 ! for each branch. The gaussian being centered at the input energy 151 ! =============================================================== 152 153 ! took out factor 1/k_obj%nkpt which intervenes only at integration time 154 155 ! gaussian smdeltaprefactor = sqrt(piinv)/elphsmear/k_obj%nkpt 156 smdeltaprefactor = max_occ*sqrt(piinv)/elphsmear 157 smdeltafactor = one/elphsmear 158 159 tmp_wtq = zero 160 omega = omega_min 161 do iomega = 1, nomega 162 omega = omega + domega 163 do ikpt=1, k_obj%nkpt 164 do ib1=1,nbranch 165 xx = smdeltafactor*(phfrq(ib1,ikpt)-omega) 166 if (abs(xx) < gaussmaxarg) then 167 tmp_wtq(ib1,ikpt,iomega) = exp(-xx*xx)*smdeltaprefactor 168 end if 169 end do 170 end do 171 end do 172 end if ! if telphint 173 174 175 end subroutine ep_ph_weights