TABLE OF CONTENTS
ABINIT/get_veloc_tr [ Functions ]
NAME
get_veloc_tr
FUNCTION
calculate the (in) and (out) velocity factors for transport
COPYRIGHT
Copyright (C) 2004-2018 ABINIT group (JPC) 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
elph_ds elph_ds%nFSband = number of bands included in the FS integration elph_ds%k_fine%nkpt = number of kpts included in the FS integration elph_ds%nFSband = number of bands included in the FS integration elph_ds%minFSband = index of the lowest FS band elph_ds%nqpt_full = number of Q pts elph_ds%nqptirred = number of irreducible Q pts to index the GS electronic states : kphon_full2irr = mapping of full FS kpts to irreducible ones FSfullpqtofull = mapping of k + q to k FSirredtoGS = mapping of irreducible kpoints to GS set
OUTPUT
elph_tr_ds%FSelecveloc_sq = avergae FS electronic velocity
PARENTS
elphon,get_nv_fs_en,get_nv_fs_temp
CHILDREN
SOURCE
39 #if defined HAVE_CONFIG_H 40 #include "config.h" 41 #endif 42 43 #include "abi_common.h" 44 45 46 subroutine get_veloc_tr(elph_ds,elph_tr_ds) 47 48 use defs_basis 49 use defs_elphon 50 use m_profiling_abi 51 use m_errors 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_veloc_tr' 57 !End of the abilint section 58 59 implicit none 60 61 62 !Arguments ------------------------------------ 63 !scalars 64 65 !arrays 66 type(elph_type),intent(in) :: elph_ds 67 type(elph_tr_type),intent(inout) :: elph_tr_ds 68 69 70 !Local variables------------------------------- 71 !scalars 72 integer :: ikpt_fine 73 integer :: ib1,fib1,isppol, ii 74 real(dp) :: eta2 75 !arrays 76 real(dp) :: elvelock(3) 77 78 ! ********************************************************************* 79 80 !DEBUG 81 !write(std_out,*)' get_veloc_tr : enter ' 82 !ENDDEBUG 83 84 ABI_CHECK(allocated(elph_tr_ds%FSelecveloc_sq),"FSele not associated") 85 86 87 !precalculate the Fermi speed modulus squared 88 elph_tr_ds%FSelecveloc_sq = zero 89 do isppol=1,elph_ds%nsppol 90 do ikpt_fine=1,elph_ds%k_fine%nkpt 91 do ib1=1,elph_ds%nFSband 92 fib1=ib1+elph_ds%minFSband-1 93 elvelock(:)=elph_tr_ds%el_veloc(ikpt_fine,fib1,:,isppol) 94 do ii=1, 3 95 eta2=elvelock(ii)*elvelock(ii) 96 elph_tr_ds%FSelecveloc_sq(ii, isppol)=elph_tr_ds%FSelecveloc_sq(ii, isppol)& 97 & +eta2*elph_ds%k_fine%wtk(ib1,ikpt_fine,isppol) 98 end do 99 end do 100 end do 101 elph_tr_ds%FSelecveloc_sq(:,isppol) = elph_tr_ds%FSelecveloc_sq(:,isppol)/elph_ds%k_fine%nkpt/elph_ds%n0(isppol) 102 ! for factor 1/elph_ds%n0(isppol) see eq 12 of Allen prb 17 3725: sum of v**2 over all k gives n0 times FSelecveloc_sq 103 end do ! end isppol 104 write (std_out,*) ' get_veloc_tr: FSelecveloc_sq ', elph_tr_ds%FSelecveloc_sq 105 106 write (std_out,*) 'out of get_veloc_tr' 107 108 end subroutine get_veloc_tr