TABLE OF CONTENTS
ABINIT/integrate_gamma [ Functions ]
NAME
integrate_gamma
FUNCTION
This routine integrates the electron phonon coupling matrix over the kpoints on the fermi surface. A dependency on qpoint remains for gamma_qpt
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
elph_ds = elphon datastructure with data and dimensions elph_ds%qpt_full = qpoint coordinates elph_ds%nqptirred = number of irred qpoints elph_ds%qirredtofull = indexing of the GKK qpoints found FSfullpqtofull = mapping of k+q to k
OUTPUT
elph_ds = modified elph_ds%gamma_qpt and created elph_ds%gamma_rpt
PARENTS
elphon
CHILDREN
get_rank_1kpt,wrtout,xmpi_sum
SOURCE
37 #if defined HAVE_CONFIG_H 38 #include "config.h" 39 #endif 40 41 #include "abi_common.h" 42 43 subroutine integrate_gamma(elph_ds,FSfullpqtofull) 44 45 use defs_basis 46 use defs_elphon 47 use m_profiling_abi 48 use m_kptrank 49 use m_errors 50 use m_xmpi 51 52 !This section has been created automatically by the script Abilint (TD). 53 !Do not modify the following lines by hand. 54 #undef ABI_FUNC 55 #define ABI_FUNC 'integrate_gamma' 56 use interfaces_14_hidewrite 57 !End of the abilint section 58 59 implicit none 60 61 !Arguments ------------------------------------ 62 !scalars 63 type(elph_type),intent(inout) :: elph_ds 64 !arrays 65 integer,intent(in) :: FSfullpqtofull(elph_ds%k_phon%nkpt,elph_ds%nqpt_full) 66 67 !Local variables------------------------------- 68 !scalars 69 integer :: comm,ikpt_phon,ikpt_phonq,ib1,ib2,ibeff,iqpt,iqpt_fullbz,isppol,ierr 70 integer :: irec, symrankkpt_phon,nbranch,nsppol,ngkkband, ik_this_proc 71 character(len=500) :: message 72 character(len=fnlen) :: fname 73 !arrays 74 real(dp),allocatable :: tmp_gkk(:,:,:,:) 75 76 ! ************************************************************************* 77 78 comm = xmpi_world 79 80 write (message,'(3a)')ch10,' entering integrate_gamma ',ch10 81 call wrtout(std_out,message,'COLL') 82 83 nsppol = elph_ds%nsppol 84 nbranch = elph_ds%nbranch 85 ngkkband = elph_ds%ngkkband 86 87 ABI_ALLOCATE(elph_ds%gamma_qpt,(2,nbranch**2,nsppol,elph_ds%nqpt_full)) 88 elph_ds%gamma_qpt = zero 89 90 ABI_ALLOCATE(tmp_gkk ,(2,ngkkband**2,nbranch**2,nsppol)) 91 92 if (elph_ds%gkqwrite == 0) then 93 call wrtout(std_out,' integrate_gamma : keeping gamma matrices in memory','COLL') 94 else if (elph_ds%gkqwrite == 1) then 95 fname=trim(elph_ds%elph_base_name) // '_GKKQ' 96 write (message,'(2a)')' integrate_gamma : reading gamma matrices from file ',trim(fname) 97 call wrtout(std_out,message,'COLL') 98 else 99 write (message,'(a,i0)')' Wrong value for gkqwrite = ',elph_ds%gkqwrite 100 MSG_BUG(message) 101 end if 102 103 104 105 do iqpt=1,elph_ds%nqptirred 106 iqpt_fullbz = elph_ds%qirredtofull(iqpt) 107 call get_rank_1kpt (elph_ds%k_phon%kpt(:,iqpt_fullbz),symrankkpt_phon, elph_ds%k_phon%kptrank_t) 108 write (std_out,*) ' iqpt_fullbz in qpt grid only, rank ', iqpt_fullbz, symrankkpt_phon 109 110 do ik_this_proc =1,elph_ds%k_phon%my_nkpt 111 ikpt_phon = elph_ds%k_phon%my_ikpt(ik_this_proc) 112 113 if (elph_ds%gkqwrite == 0) then 114 tmp_gkk = elph_ds%gkk_qpt(:,:,:,ik_this_proc,:,iqpt) 115 else if (elph_ds%gkqwrite == 1) then 116 irec = (iqpt-1)*elph_ds%k_phon%my_nkpt+ik_this_proc 117 if (ikpt_phon == 1) then 118 write (std_out,*) ' integrate_gamma read record ', irec 119 end if 120 read (elph_ds%unitgkq,REC=irec) tmp_gkk(:,:,:,:) 121 end if 122 123 do isppol=1,nsppol 124 ikpt_phonq = FSfullpqtofull(ikpt_phon,iqpt_fullbz) 125 ! 126 do ib1=1,ngkkband 127 do ib2=1,ngkkband 128 ibeff = ib2+(ib1-1)*ngkkband 129 elph_ds%gamma_qpt(:,:,isppol,iqpt_fullbz) = elph_ds%gamma_qpt(:,:,isppol,iqpt_fullbz) + & 130 & tmp_gkk(:,ibeff,:,isppol)& 131 & *elph_ds%gkk_intweight(ib1,ikpt_phon,isppol)*elph_ds%gkk_intweight(ib2,ikpt_phonq,isppol) 132 ! NOTE: if ngkkband==1 we are using trivial weights since average 133 ! over bands was done in normsq_gkk (nmsq_gam_sumFS or nmsq_pure_gkk) 134 end do ! ib2 135 end do ! ib1 136 end do ! isppol 137 end do ! ikpt_phon 138 end do ! iqpt 139 140 call xmpi_sum (elph_ds%gamma_qpt, comm, ierr) 141 142 ABI_DEALLOCATE(tmp_gkk) 143 144 !need prefactor of 1/nkpt for each integration over 1 kpoint index. NOT INCLUDED IN elph_ds%gkk_intweight 145 do iqpt=1,elph_ds%nqptirred 146 iqpt_fullbz = elph_ds%qirredtofull(iqpt) 147 ! elph_ds%gamma_qpt(:,:,:,iqpt_fullbz) = elph_ds%gamma_qpt(:,:,:,iqpt_fullbz) / elph_ds%k_phon%nkpt / n0(1) / n0(1) 148 ! elph_ds%gamma_qpt(:,:,:,iqpt_fullbz) = elph_ds%gamma_qpt(:,:,:,iqpt_fullbz) / elph_ds%k_phon%nkpt / elph_ds%k_phon%nkpt 149 elph_ds%gamma_qpt(:,:,:,iqpt_fullbz) = elph_ds%gamma_qpt(:,:,:,iqpt_fullbz) * elph_ds%occ_factor / elph_ds%k_phon%nkpt 150 end do 151 152 call wrtout(std_out,' integrate_gamma: gamma matrices have been calculated for recip space and irred qpoints ',"COLL") 153 154 end subroutine integrate_gamma