1
!{\src2tex{textfont=tt}}
2
!!****f* ABINIT/integrate_gamma
8
!! This routine integrates the electron phonon coupling matrix
9
!! over the kpoints on the fermi surface. A dependency on qpoint
10
!! or irpt (real space) remains, for gamma_qpt and gamma_rpt resp.
13
!! Copyright (C) 2004-2007 ABINIT group (MVer)
14
!! This file is distributed under the terms of the
15
!! GNU General Public Licence, see ~abinit/COPYING
16
!! or http://www.gnu.org/copyleft/gpl.txt .
17
!! For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
20
!! elph_ds = elphon datastructure with data and dimensions
21
!! FSfullpqtofull = mapping of k+q to k
22
!! gprim = reciprocal space lattice vectors
23
!! n0 = DOS at fermi level
24
!! natom = number of atoms
25
!! nrpt = number of real space points for FT
26
!! rpt = coordinates of real space points for FT
27
!! spqpt = qpoint coordinates
28
!! wghatm = weights for FT of real-space points
31
!! elph_ds = modified elph_ds%gamma_qpt and created elph_ds%gamma_rpt
39
!! ftgam,leave_new,wrtout
43
#if defined HAVE_CONFIG_H
47
subroutine integrate_gamma(elph_ds,FSfullpqtofull,gprim,n0,natom,nrpt,rpt,spqpt,wghatm)
53
!This section has been created automatically by the script Abilint (TD). Do not modify these by hand.
54
#ifdef HAVE_FORTRAN_INTERFACES
55
use interfaces_01manage_mpi
56
use interfaces_17ddb, except_this_one => integrate_gamma
58
!End of the abilint section
62
!Arguments ------------------------------------
64
integer,intent(in) :: natom,nrpt
65
real(dp),intent(in) :: n0
66
type(elph_type),intent(inout) :: elph_ds
68
integer,intent(in) :: FSfullpqtofull(elph_ds%nFSkpt,elph_ds%nqpt)
69
real(dp),intent(in) :: gprim(3,3),rpt(3,nrpt),spqpt(3,elph_ds%nqpt)
70
real(dp),intent(in) :: wghatm(natom,natom,nrpt)
72
!Local variables-------------------------------
74
integer :: iFSkpt,iFSkptq,ib1,ib2,ierr,iqpt,irpt,qtor
75
character(len=500) :: message
76
character(len=fnlen) :: fname
78
real(dp),allocatable :: tmp_gkk(:,:,:,:,:,:)
80
! *************************************************************************
82
allocate(elph_ds%gamma_qpt(2,elph_ds%nbranch,elph_ds%nbranch,elph_ds%nqpt),stat=ierr)
84
write (message,'(3a)')' integrate_gamma : ERROR- ',ch10,&
85
& ' trying to allocate array elph_ds%gamma_qpt '
86
call wrtout(06,message,'COLL')
87
call leave_new('COLL')
89
elph_ds%gamma_qpt(:,:,:,:) = zero
91
allocate (elph_ds%gamma_rpt(2,elph_ds%nbranch,elph_ds%nbranch,nrpt),stat=ierr)
93
write (message,'(3a)')' integrate_gamma : ERROR- ',ch10,&
94
& ' trying to allocate array elph_ds%gamma_rpt '
95
call wrtout(06,message,'COLL')
96
call leave_new('COLL')
98
elph_ds%gamma_rpt(:,:,:,:) = zero
100
if (elph_ds%gkqwrite == 0) then
101
write (message,'(a)')' integrate_gamma : keeping gamma matrices in memory'
102
call wrtout(06,message,'COLL')
103
!NOTE: if ngkkband==1 we are using trivial weights since average
104
! over bands was done in normsq_gkk (nmsq_gam_sumFS or nmsq_pure_gkk)
106
do iqpt=1,elph_ds%nqpt
107
do iFSkpt=1,elph_ds%nFSkpt
108
iFSkptq = FSfullpqtofull(iFSkpt,iqpt)
110
do ib1=1,elph_ds%ngkkband
111
do ib2=1,elph_ds%ngkkband
112
elph_ds%gamma_qpt(:,:,:,iqpt) = elph_ds%gamma_qpt(:,:,:,iqpt) + &
113
& elph_ds%gkk_qpt(:,ib2,ib1,:,:,iFSkpt,iqpt)&
114
& *elph_ds%gkk_intweight(ib1,iFSkpt)*elph_ds%gkk_intweight(ib2,iFSkptq)
119
else if (elph_ds%gkqwrite == 1) then
121
fname=trim(elph_ds%elph_base_name) // '_GKKQ'
122
write (message,'(2a)')' integrate_gamma : reading gamma matrices from file ',trim(fname)
123
call wrtout(06,message,'COLL')
125
allocate (tmp_gkk (2,elph_ds%ngkkband,elph_ds%ngkkband,elph_ds%nbranch,elph_ds%nbranch,elph_ds%nFSkpt),stat=ierr)
127
write (message,'(3a)')' integrate_gamma : ERROR- ',ch10,&
128
& ' trying to allocate array tmp_gkk '
129
call wrtout(06,message,'COLL')
130
call leave_new('COLL')
133
do iqpt=1,elph_ds%nqpt
134
read (elph_ds%unitgkq,REC=iqpt) tmp_gkk
135
do iFSkpt=1,elph_ds%nFSkpt
136
iFSkptq = FSfullpqtofull(iFSkpt,iqpt)
138
do ib1=1,elph_ds%ngkkband
139
do ib2=1,elph_ds%ngkkband
140
elph_ds%gamma_qpt(:,:,:,iqpt) = elph_ds%gamma_qpt(:,:,:,iqpt) + &
141
& tmp_gkk(:,ib2,ib1,:,:,iFSkpt)&
142
& *elph_ds%gkk_intweight(ib1,iFSkpt)*elph_ds%gkk_intweight(ib2,iFSkptq)
147
! if (iqpt == 1) then
148
! write (102,*) ' tmp_gkk ==== in integrate_gamma '
149
! do iFSkpt=1,elph_ds%nFSkpt
150
! write (102,*) iFSkpt
151
! write (102,'(9(2E14.6,1x))') tmp_gkk(:,:,:,:,:,iFSkpt)
158
write (message,'(3a,i3)')' integrate_gamma : BUG-',ch10,&
159
& ' Wrong value for gkqwrite = ',elph_ds%gkqwrite
160
call wrtout(06,message,'COLL')
161
call leave_new('COLL')
164
! need prefactor of 1/nkpt for each integration over 1 kpoint index.
165
! NOT INCLUDED IN elph_ds%gkk_intweight
166
do iqpt=1,elph_ds%nqpt
167
!elph_ds%gamma_qpt(:,:,:,iqpt) = elph_ds%gamma_qpt(:,:,:,iqpt) / elph_ds%nFSkpt / n0 / n0
168
!elph_ds%gamma_qpt(:,:,:,iqpt) = elph_ds%gamma_qpt(:,:,:,iqpt) / elph_ds%nFSkpt / elph_ds%nFSkpt
169
elph_ds%gamma_qpt(:,:,:,iqpt) = elph_ds%gamma_qpt(:,:,:,iqpt) / elph_ds%nFSkpt
173
!write (100,*) ' gamma_qpt on qpts actually calculated ==== in integrate_gamma '
174
!do iqpt=1,elph_ds%nqpt
176
! write (100,'(9(2E14.6,1x))') elph_ds%gamma_qpt(:,:,:,iqpt)
180
! Now FT to real space too
181
write (message,'(a)')' integrate_gamma : Fourier transforming gamma matrices to real space'
182
call wrtout(06,message,'COLL')
185
call ftgam(wghatm,elph_ds%gamma_qpt,elph_ds%gamma_rpt,gprim,natom,&
186
& elph_ds%nqpt,nrpt,qtor,rpt,spqpt)
188
write (message,'(a)')' integrate_gamma : gamma matrices are in real space '
189
call wrtout(06,message,'COLL')
192
end subroutine integrate_gamma