The Elk Code
gwchgk.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2018 P. Elliott, J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine gwchgk(ik,chgk)
7 use modmain
8 use modgw
9 use modomp
10 implicit none
11 ! arguments
12 integer, intent(in) :: ik
13 real(8), intent(out) :: chgk
14 ! local variables
15 integer ist,iw,i,nthd
16 real(8) e
17 complex(8) z1
18 ! automatic arrays
19 complex(8) gs(nstsv),g(nstsv,nstsv),ge(4,nstsv)
20 ! allocatable arrays
21 complex(8), allocatable :: se(:,:,:)
22 ! external functions
23 complex(8), external :: gwtails
24 ! read the self-energy from file
25 allocate(se(nstsv,nstsv,0:nwfm))
26 call getgwsefm(ik,se)
27 ! zero the charge
28 chgk=0.d0
29 ! loop over fermionic Matsubara frequencies
30 call holdthd(nwfm+1,nthd)
31 !$OMP PARALLEL DO DEFAULT(SHARED) &
32 !$OMP PRIVATE(gs,g,ist,e,z1,i) &
33 !$OMP REDUCTION(+:chgk) &
34 !$OMP SCHEDULE(DYNAMIC) &
35 !$OMP NUM_THREADS(nthd)
36 do iw=0,nwfm
37 ! compute the diagonal matrix Gₛ
38  do ist=1,nstsv
39  e=evalsv(ist,ik)-efermi
40  gs(ist)=1.d0/(wfm(iw)-e)
41  end do
42 ! compute 1 - Gₛ Σ
43  do ist=1,nstsv
44  z1=-gs(ist)
45  g(ist,1:nstsv)=z1*se(ist,1:nstsv,iw)
46  g(ist,ist)=g(ist,ist)+1.d0
47  end do
48 ! invert this matrix
49  call zminv(nstsv,g)
50 ! take the trace of G = (1 - Gₛ Σ)⁻¹ Gₛ
51  do ist=1,nstsv
52  g(ist,ist)=g(ist,ist)*gs(ist)
53  chgk=chgk+dble(g(ist,ist))
54  end do
55 ! store the Green's function at the end point frequencies
56  i=0
57  if (iw == 0) i=1
58  if (iw == 1) i=2
59  if (iw == nwfm-1) i=3
60  if (iw == nwfm) i=4
61  if (i /= 0) then
62  do ist=1,nstsv
63  ge(i,ist)=g(ist,ist)
64  end do
65  end if
66 end do
67 !$OMP END PARALLEL DO
68 call freethd(nthd)
69 ! add the Matsubara tails analytically
70 do ist=1,nstsv
71  chgk=chgk+dble(gwtails(ge(:,ist)))
72 end do
73 chgk=chgk*wkpt(ik)*occmax*kboltz*tempk
74 deallocate(se)
75 end subroutine
76 
real(8) efermi
Definition: modmain.f90:907
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:921
subroutine getgwsefm(ik, se)
Definition: getgwsefm.f90:7
subroutine gwchgk(ik, chgk)
Definition: gwchgk.f90:7
Definition: modomp.f90:6
real(8), parameter kboltz
Definition: modmain.f90:1263
real(8), dimension(:), allocatable wkpt
Definition: modmain.f90:475
real(8) occmax
Definition: modmain.f90:901
real(8) tempk
Definition: modmain.f90:684
Definition: modgw.f90:6
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer nwfm
Definition: modgw.f90:19
subroutine zminv(n, a)
Definition: zminv.f90:7
complex(8), dimension(:), allocatable wfm
Definition: modgw.f90:25