TABLE OF CONTENTS
ABINIT/nmsq_gam_sumfs [ Functions ]
NAME
nmsq_gam_sumfs
FUNCTION
Calculate gamma matrices from original h1_mat_el_sq matrix elements averaging over bands near 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
displ_red = phonon mode displacement vectors, post-multiplied by gprim matrix (ie. turned to reduced coordinates) eigvec = eigenvectors of phonons (to turn to cartesian coord frame) elph_ds = datastructure with gkk matrix elements FSfullpqtofull = mapping of k+q to k kpt_phon = coordinates of kpoints near to FS h1_mat_el_sq = matrix elements $<psi_{k+q,m} | H^{1} | psi_{k,n}>$ squared iqptirred = index of present qpoint
OUTPUT
accum_mat = matrix for accumulating FS average of gkk (gamma matrix -> linewidths) accum_mat2 = matrix for accumulating FS average of gamma matrix with good prefactors
PARENTS
normsq_gkq
CHILDREN
gam_mult_displ,zgemm
SOURCE
40 #if defined HAVE_CONFIG_H 41 #include "config.h" 42 #endif 43 44 #include "abi_common.h" 45 46 subroutine nmsq_gam_sumFS(accum_mat,accum_mat2,displ_red,eigvec,elph_ds,FSfullpqtofull,& 47 & h1_mat_el_sq,iqptirred) 48 49 use defs_basis 50 use defs_elphon 51 use m_profiling_abi 52 use m_errors 53 54 !This section has been created automatically by the script Abilint (TD). 55 !Do not modify the following lines by hand. 56 #undef ABI_FUNC 57 #define ABI_FUNC 'nmsq_gam_sumFS' 58 use interfaces_77_ddb, except_this_one => nmsq_gam_sumFS 59 !End of the abilint section 60 61 implicit none 62 63 !Arguments ------------------------------------ 64 !scalars 65 integer,intent(in) :: iqptirred 66 type(elph_type),intent(inout) :: elph_ds 67 !arrays 68 integer,intent(in) :: FSfullpqtofull(elph_ds%k_phon%nkpt,elph_ds%nqpt_full) 69 real(dp),intent(in) :: displ_red(2,elph_ds%nbranch,elph_ds%nbranch) 70 real(dp),intent(in) :: eigvec(2,elph_ds%nbranch,elph_ds%nbranch) 71 real(dp),intent(inout) :: & 72 & h1_mat_el_sq(2,elph_ds%nFSband*elph_ds%nFSband,elph_ds%nbranch*elph_ds%nbranch,elph_ds%k_phon%my_nkpt,elph_ds%nsppol) 73 real(dp),intent(inout) :: accum_mat(2,elph_ds%nbranch,elph_ds%nbranch,elph_ds%nsppol) 74 real(dp),intent(inout) :: accum_mat2(2,elph_ds%nbranch,elph_ds%nbranch,elph_ds%nsppol) 75 76 !Local variables------------------------------- 77 !scalars 78 integer :: ikpt_phon,ikpt_phonq,ib1,ib2,ibeff,ibranch,ipert1,isppol,jbranch,iqpt_fullbz 79 integer :: ik_this_proc 80 real(dp) :: sd1,sd2 81 character(len=500) :: message 82 !arrays 83 real(dp) :: gkq_sum_bands(2,elph_ds%nbranch,elph_ds%nbranch) 84 real(dp) :: tmp_gkq_sum_bands(2,elph_ds%nbranch,elph_ds%nbranch) 85 real(dp) :: tmp_mat2(2,elph_ds%nbranch,elph_ds%nbranch) 86 real(dp),allocatable :: zgemm_tmp_mat(:,:,:) 87 88 ! ************************************************************************* 89 90 if (elph_ds%ep_keepbands /= 0) then 91 write (message,'(a,i0)')' elph_ds%ep_keepbands should be 0 in order to average over bands!',elph_ds%ep_keepbands 92 MSG_ERROR(message) 93 end if 94 95 iqpt_fullbz = elph_ds%qirredtofull(iqptirred) 96 97 98 !MG20060603 NOTE: 99 !accum_mat and accum_mat2 are real, the imaginary part is used for debugging purpose 100 !accum_mat2 is used to store the phonon-linewidhts before interpolation 101 102 ABI_ALLOCATE(zgemm_tmp_mat ,(2,elph_ds%nbranch,elph_ds%nbranch)) 103 104 do isppol=1,elph_ds%nsppol 105 do ik_this_proc =1, elph_ds%k_phon%my_nkpt 106 ikpt_phon = elph_ds%k_phon%my_ikpt(ik_this_proc) 107 108 ikpt_phonq = FSfullpqtofull(ikpt_phon,iqpt_fullbz) 109 110 gkq_sum_bands = zero 111 tmp_gkq_sum_bands = zero 112 113 114 do ib1=1,elph_ds%nFSband 115 ! weights for distance from the fermi surface 116 sd1 = elph_ds%k_phon%wtk(ib1,ikpt_phon,isppol) 117 118 do ib2=1,elph_ds%nFSband 119 ! weights for distance from the fermi surface 120 sd2 = elph_ds%k_phon%wtk(ib2,ikpt_phonq,isppol) 121 ibeff=ib2+(ib1-1)*elph_ds%nFSband 122 123 zgemm_tmp_mat = reshape(h1_mat_el_sq(:,ibeff,:,isppol,ik_this_proc),(/2,elph_ds%nbranch,elph_ds%nbranch/)) 124 125 call gam_mult_displ(elph_ds%nbranch, displ_red, zgemm_tmp_mat, tmp_mat2) 126 127 ! sum over bands in gkq_sum_bands 128 do ipert1=1,elph_ds%nbranch 129 gkq_sum_bands(1,ipert1,ipert1) = gkq_sum_bands(1,ipert1,ipert1) + sd1*sd2*tmp_mat2(1,ipert1,ipert1) 130 end do 131 132 133 134 end do 135 end do 136 ! END loop over bands 137 138 ! summing over k points, still diagonal in jbranch 139 accum_mat(:,:,:,isppol) = accum_mat(:,:,:,isppol) + gkq_sum_bands(:,:,:) 140 accum_mat2(:,:,:,isppol) = accum_mat2(:,:,:,isppol) + gkq_sum_bands(:,:,:) 141 142 ! summed over bands, now turn to cartesian coordinates 143 144 ! Final Gamma matrix (hermitian) = E * D_g * E^{+} 145 ! Where E^{+} is the hermitian conjugate of the eigenvector matrix E 146 ! And D_g is the diagonal matrix of values of gamma for this qpoint 147 148 ! Here gkq_sum_bands is indexed with real phonon modes (not atom+idir) 149 ! turn gkq_sum_bands to atom+cartesian coordinates (instead of normal coordinates for qpoint) 150 ! This is not a full matrix multiplication, just vector one, by 151 ! gkq_sum_bands(1,jbranch,jbranch) 152 tmp_mat2(:,:,:) = zero 153 do ibranch =1,elph_ds%nbranch 154 do jbranch =1,elph_ds%nbranch 155 tmp_mat2(1,ibranch,jbranch) = tmp_mat2(1,ibranch,jbranch) + & 156 & eigvec(1,ibranch,jbranch) * & 157 & gkq_sum_bands(1,jbranch,jbranch) 158 tmp_mat2(2,ibranch,jbranch) = tmp_mat2(2,ibranch,jbranch) + & 159 & eigvec(2,ibranch,jbranch) * & 160 & gkq_sum_bands(1,jbranch,jbranch) 161 end do 162 end do 163 164 ! here eigvec is transposed and complex conjugated. 165 zgemm_tmp_mat=zero 166 call zgemm('n','c',elph_ds%nbranch,elph_ds%nbranch,elph_ds%nbranch,cone,& 167 & tmp_mat2,elph_ds%nbranch,eigvec,elph_ds%nbranch,czero,zgemm_tmp_mat,elph_ds%nbranch) 168 169 gkq_sum_bands = zgemm_tmp_mat 170 171 ! ! gamma matrix contribution in cartesian coordinates (ie interpolatable form) 172 ! gamma matrix contribution in reduced coordinates (ie interpolatable form) 173 h1_mat_el_sq(:,1,:,ik_this_proc,isppol) = reshape(gkq_sum_bands(:,:,:),(/2,elph_ds%nbranch*elph_ds%nbranch/)) 174 175 ! accum_mat(:,:,:,isppol) = accum_mat(:,:,:,isppol) + gkq_sum_bands(:,:,:) 176 end do 177 ! END loop over kpt_phon 178 end do 179 !END loop over sppol 180 181 ABI_DEALLOCATE(zgemm_tmp_mat) 182 183 end subroutine nmsq_gam_sumFS