4
1
mirror of https://github.com/pfloos/quack synced 2024-06-02 11:25:28 +02:00
quack/src/GT/GTeh_QP_graph.f90

87 lines
2.1 KiB
Fortran

subroutine GTeh_QP_graph(eta,nBas,nC,nO,nV,nR,nS,eHF,Om,rhoL,rhoR,eGTlin,eGT,Z)
! Compute the graphical solution of the QP equation
implicit none
include 'parameters.h'
! Input variables
integer,intent(in) :: nBas
integer,intent(in) :: nC
integer,intent(in) :: nO
integer,intent(in) :: nV
integer,intent(in) :: nR
integer,intent(in) :: nS
double precision,intent(in) :: eta
double precision,intent(in) :: eHF(nBas)
double precision,intent(in) :: Om(nS)
double precision,intent(in) :: rhoL(nBas,nBas,nS)
double precision,intent(in) :: rhoR(nBas,nBas,nS)
double precision,intent(in) :: eGTlin(nBas)
! Local variables
integer :: p
integer :: nIt
integer,parameter :: maxIt = 64
double precision,parameter :: thresh = 1d-6
double precision,external :: GTeh_SigC,GTeh_dSigC
double precision :: SigC,dSigC
double precision :: f,df
double precision :: w
! Output variables
double precision,intent(out) :: eGT(nBas)
double precision,intent(out) :: Z(nBas)
! Run Newton's algorithm to find the root
do p=nC+1,nBas-nR
write(*,*) '-----------------'
write(*,'(A10,I3)') 'Orbital ',p
write(*,*) '-----------------'
w = eGTlin(p)
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
SigC = GTeh_SigC(p,w,eta,nBas,nC,nO,nV,nR,nS,eGTlin,Om,rhoL,rhoR)
dSigC = GTeh_dSigC(p,w,eta,nBas,nC,nO,nV,nR,nS,eGTlin,Om,rhoL,rhoR)
f = w - eHF(p) - SigC
df = 1d0/(1d0 - dSigC)
w = w - df*f
write(*,'(A3,I3,A1,1X,3F15.9)') 'It.',nIt,':',w*HaToeV,df,f
end do
if(nIt == maxIt) then
eGT(p) = eGTlin(p)
write(*,*) 'Newton root search has not converged!'
else
eGT(p) = w
Z(p) = df
write(*,'(A32,F16.10)') 'Quasiparticle energy (eV) ',eGT(p)*HaToeV
write(*,*)
end if
end do
end subroutine