The Elk Code
 
Loading...
Searching...
No Matches
dysonr.f90
Go to the documentation of this file.
1
2! Copyright (C) 2016 A. Davydov, A. Sanna, J. K. Dewhurst, S. Sharma and
3! E. K. U. Gross. This file is distributed under the terms of the GNU General
4! Public License. See the file COPYING for license details.
5
6subroutine dysonr(ik,wr,sem,sf)
7use modmain
8use modgw
9use modomp
10implicit none
11! arguments
12integer, intent(in) :: ik
13real(8), intent(in) :: wr(nwplot)
14complex(8), intent(in) :: sem(nstsv,nstsv,0:nwfm)
15real(8), intent(out) :: sf(nwplot)
16! local variables
17integer ist,jst,iw,nthd
18real(8) w,e,sm,t1
19complex(8) z1
20! automatic arrays
21complex(8) gs(nstsv),g(nstsv,nstsv)
22! allocatable arrays
23complex(8), allocatable :: ser(:,:,:)
24allocate(ser(nstsv,nstsv,nwplot))
25ser(1:nstsv,1:nstsv,1:nwplot)=0.d0
26call holdthd(nstsv,nthd)
27!$OMP PARALLEL DO DEFAULT(SHARED) &
28!$OMP PRIVATE(ist) &
29!$OMP SCHEDULE(DYNAMIC) &
30!$OMP NUM_THREADS(nthd)
31do jst=1,nstsv
32 do ist=1,nstsv
33 if (tsediag.and.(ist /= jst)) cycle
34! perform numerical analytic continuation from the imaginary to the real axis
35 call acgwse(ist,jst,sem,wr,ser)
36 end do
37end do
38!$OMP END PARALLEL DO
39call freethd(nthd)
40! solve the Dyson equation for each frequency
41call holdthd(nwplot,nthd)
42!$OMP PARALLEL DO DEFAULT(SHARED) &
43!$OMP PRIVATE(gs,g,w,ist,jst) &
44!$OMP PRIVATE(e,t1,z1,sm) &
45!$OMP SCHEDULE(DYNAMIC) &
46!$OMP NUM_THREADS(nthd)
47do iw=1,nwplot
48 w=wr(iw)
49! compute the diagonal matrix Gₛ
50 do ist=1,nstsv
51 e=evalsv(ist,ik)-efermi
52 t1=sign(swidth,e)
53 gs(ist)=1.d0/cmplx(w-e,t1,8)
54 end do
55! compute 1 - Gₛ Σ
56 do ist=1,nstsv
57 z1=-gs(ist)
58 g(ist,1:nstsv)=z1*ser(ist,1:nstsv,iw)
59 g(ist,ist)=g(ist,ist)+1.d0
60 end do
61! invert this matrix
62 call zminv(nstsv,g)
63! compute G = (1 - Gₛ Σ)⁻¹ Gₛ
64 do jst=1,nstsv
65 z1=gs(jst)
66 g(1:nstsv,jst)=g(1:nstsv,jst)*z1
67 end do
68! determine the spectral function
69 sm=0.d0
70 do ist=1,nstsv
71 sm=sm+abs(aimag(g(ist,ist)))
72 end do
73 sf(iw)=sm*occmax/pi
74end do
75!$OMP END PARALLEL DO
76call freethd(nthd)
77deallocate(ser)
78end subroutine
79
subroutine acgwse(ist, jst, sem, wr, ser)
Definition acgwse.f90:7
subroutine dysonr(ik, wr, sem, sf)
Definition dysonr.f90:7
Definition modgw.f90:6
logical tsediag
Definition modgw.f90:27
real(8) efermi
Definition modmain.f90:904
real(8), parameter pi
Definition modmain.f90:1229
real(8) swidth
Definition modmain.f90:892
real(8) occmax
Definition modmain.f90:898
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine zminv(n, a)
Definition zminv.f90:7