The Elk Code
gwdmatk.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2018 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 gwdmatk(ik)
7 use modmain
8 use modgw
9 use modomp
10 implicit none
11 ! arguments
12 integer, intent(in) :: ik
13 ! local variables
14 integer ist,jst,iw,i,nthd
15 real(8) e,t1
16 complex(8) z1
17 ! automatic arrays
18 complex(8) gs(nstsv),g(nstsv,nstsv),ge(4,nstsv,nstsv)
19 complex(8) evecsv(nstsv,nstsv),d(nstsv,nstsv),a(nstsv,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 density matrix
28 d(1:nstsv,1:nstsv)=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,jst,e,z1,i) &
33 !$OMP REDUCTION(+:d) &
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 ! compute G = (1 - Gₛ Σ)⁻¹ Gₛ
51  do jst=1,nstsv
52  z1=gs(jst)
53  g(1:nstsv,jst)=g(1:nstsv,jst)*z1
54  end do
55 ! add to the density matrix
56  d(1:nstsv,1:nstsv)=d(1:nstsv,1:nstsv)+g(1:nstsv,1:nstsv)
57 ! store the Green's function at the end point frequencies
58  i=0
59  if (iw == 0) i=1
60  if (iw == 1) i=2
61  if (iw == nwfm-1) i=3
62  if (iw == nwfm) i=4
63  if (i /= 0) ge(i,1:nstsv,1:nstsv)=g(1:nstsv,1:nstsv)
64 end do
65 !$OMP END PARALLEL DO
66 call freethd(nthd)
67 ! add the Matsubara tails analytically
68 do jst=1,nstsv
69  do ist=1,nstsv
70  d(ist,jst)=d(ist,jst)+gwtails(ge(:,ist,jst))
71  end do
72 end do
73 ! multiply by 1/β
74 t1=kboltz*tempk
75 d(1:nstsv,1:nstsv)=t1*d(1:nstsv,1:nstsv)
76 deallocate(se)
77 ! make density matrix Hermitian
78 do ist=1,nstsv
79  do jst=1,ist-1
80  z1=0.5d0*(d(ist,jst)+conjg(d(jst,ist)))
81  d(ist,jst)=z1
82  d(jst,ist)=conjg(z1)
83  end do
84  d(ist,ist)=dble(d(ist,ist))
85 end do
86 ! diagonalise the density matrix for the natural orbitals and occupation numbers
87 call eveqnzh(nstsv,nstsv,d,occsv(:,ik))
88 occsv(1:nstsv,ik)=occsv(1:nstsv,ik)*occmax
89 ! get the second-variational eigenvectors from file
90 call getevecsv(filext,ik,vkl(:,ik),evecsv)
91 ! apply unitary transformation to the third-variational states so that they
92 ! refer to the first-variational basis
93 call zgemm('N','N',nstsv,nstsv,nstsv,zone,evecsv,nstsv,d,nstsv,zzero,a,nstsv)
94 ! write the density matrix to file as second-variational eigenvectors
95 call putevecsv(filext,ik,a)
96 end subroutine
97 
real(8) efermi
Definition: modmain.f90:907
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
character(256) filext
Definition: modmain.f90:1301
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:921
subroutine getgwsefm(ik, se)
Definition: getgwsefm.f90:7
Definition: modomp.f90:6
real(8), parameter kboltz
Definition: modmain.f90:1263
complex(8), parameter zone
Definition: modmain.f90:1240
subroutine eveqnzh(n, ld, a, w)
Definition: eveqnzh.f90:7
real(8) occmax
Definition: modmain.f90:901
real(8) tempk
Definition: modmain.f90:684
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:905
complex(8), parameter zzero
Definition: modmain.f90:1240
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
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 putevecsv(fext, ik, evecsv)
Definition: putevecsv.f90:7
subroutine zminv(n, a)
Definition: zminv.f90:7
subroutine gwdmatk(ik)
Definition: gwdmatk.f90:7
complex(8), dimension(:), allocatable wfm
Definition: modgw.f90:25