10
1
mirror of https://github.com/pfloos/quack synced 2025-01-10 21:18:23 +01:00
QuAcK/src/GF/GF2_QP_graph.f90

76 lines
1.7 KiB
Fortran
Raw Normal View History

2023-07-27 17:43:18 +02:00
subroutine GF2_QP_graph(eta,nBas,nC,nO,nV,nR,nS,eHF,ERI,eGF)
2023-06-30 19:17:35 +02:00
! Compute the graphical solution of the GF2 QP equation
implicit none
include 'parameters.h'
! Input variables
double precision,intent(in) :: eta
integer,intent(in) :: nBas,nC,nO,nV,nR,nS
double precision,intent(in) :: eHF(nBas)
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
! Local variables
integer :: p
integer :: nIt
integer,parameter :: maxIt = 64
double precision,parameter :: thresh = 1d-6
2023-07-27 17:43:18 +02:00
double precision,external :: GF2_SigC,GF2_dSigC
2023-06-30 19:17:35 +02:00
double precision :: sigC,dsigC
double precision :: f,df
double precision :: w
! Output variables
2023-07-27 17:43:18 +02:00
double precision,intent(out) :: eGF(nBas)
2023-06-30 19:17:35 +02:00
! Run Newton's algorithm to find the root
do p=nC+1,nBas-nR
write(*,*) '-----------------'
write(*,'(A10,I3)') 'Orbital ',p
write(*,*) '-----------------'
2023-07-27 17:43:18 +02:00
w = eHF(p)
2023-06-30 19:17:35 +02:00
nIt = 0
f = 1d0
write(*,'(A3,I3,A1,1X,3F15.9)') 'It.',nIt,':',w*HaToeV,f
do while (abs(f) > thresh .and. nIt < maxIt)
nIt = nIt + 1
2023-07-27 17:43:18 +02:00
sigC = GF2_SigC(p,w,eta,nBas,nC,nO,nV,nR,nS,eHF,ERI)
dsigC = GF2_dSigC(p,w,eta,nBas,nC,nO,nV,nR,nS,eHF,ERI)
2023-06-30 19:17:35 +02:00
f = w - eHF(p) - sigC
df = 1d0 - dsigC
w = w - f/df
write(*,'(A3,I3,A1,1X,3F15.9)') 'It.',nIt,':',w*HaToeV,f,sigC
end do
if(nIt == maxIt) then
write(*,*) 'Newton root search has not converged!'
else
2023-07-27 17:43:18 +02:00
eGF(p) = w
2023-06-30 19:17:35 +02:00
2023-07-27 17:43:18 +02:00
write(*,'(A32,F16.10)') 'Quasiparticle energy (eV) ',eGF(p)*HaToeV
2023-06-30 19:17:35 +02:00
write(*,*)
end if
end do
end subroutine