File: dysonr.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (94 lines) | stat: -rw-r--r-- 2,329 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94

! Copyright (C) 2016 A. Davydov, A. Sanna, J. K. Dewhurst, S. Sharma and
! E. K. U. Gross. This file is distributed under the terms of the GNU General
! Public License. See the file COPYING for license details.

subroutine dysonr(ik,wr,sem,sf)
use modmain
use modgw
use modomp
implicit none
! arguments
integer, intent(in) :: ik
real(8), intent(in) :: wr(nwplot)
complex(8), intent(in) :: sem(nstsv,nstsv,0:nwfm)
real(8), intent(out) :: sf(nwplot)
! local variables
integer ist,jst,iw
integer info,nthd
real(8) w,e,sum,t1
complex(8) z1
! allocatable arrays
integer, allocatable :: ipiv(:)
complex(8), allocatable :: ser(:,:,:),gs(:),g(:,:),work(:)
allocate(ser(nstsv,nstsv,nwplot))
ser(:,:,:)=0.d0
call omp_hold(nstsv,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(ist) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do jst=1,nstsv
  do ist=1,nstsv
    if (tsediag.and.(ist.ne.jst)) cycle
! perform analytic continuation from the imaginary to the real axis
    call acgwse(ist,jst,sem,wr,ser)
  end do
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
! solve the Dyson equation for each frequency
call omp_hold(nwplot,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(gs,g,ipiv,work) &
!$OMP PRIVATE(w,ist,jst,e,t1) &
!$OMP PRIVATE(z1,info,sum) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do iw=1,nwplot
  allocate(gs(nstsv),g(nstsv,nstsv))
  allocate(ipiv(nstsv),work(nstsv))
  w=wr(iw)
! compute the diagonal matrix G_s
  do ist=1,nstsv
    e=evalsv(ist,ik)-efermi
    t1=sign(swidth,e)
    gs(ist)=1.d0/cmplx(w-e,t1,8)
  end do
! compute 1 - G_s Sigma
  do ist=1,nstsv
    z1=-gs(ist)
    g(ist,:)=z1*ser(ist,:,iw)
    g(ist,ist)=g(ist,ist)+1.d0
  end do
! invert this matrix
  call zgetrf(nstsv,nstsv,g,nstsv,ipiv,info)
  if (info.eq.0) call zgetri(nstsv,g,nstsv,ipiv,work,nstsv,info)
  if (info.ne.0) then
    write(*,*)
    write(*,'("Error(dysonr): unable to solve the Dyson equation")')
    write(*,'(" for frequency ",I6)') iw
    write(*,*)
    stop
  end if
! compute G = (1 - G_s Sigma)^(-1) G_s
  do jst=1,nstsv
    z1=gs(jst)
    g(:,jst)=g(:,jst)*z1
  end do
! determine the spectral function
  sum=0.d0
  do ist=1,nstsv
    sum=sum+abs(aimag(g(ist,ist)))
  end do
  sf(iw)=sum*occmax/pi
  deallocate(gs,g,ipiv,work)
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
deallocate(ser)
return
end subroutine