mirror of
https://github.com/pfloos/quack
synced 2025-05-06 23:24:58 +02:00
added linearization for complex G0F2, does not work yet
This commit is contained in:
parent
f4fd58f521
commit
d49b293c36
@ -84,8 +84,7 @@ subroutine complex_cRG0F2(dotest,dophBSE,doppBSE,TDA,dBSE,dTDA,singlet,triplet,l
|
||||
|
||||
write(*,*) ' *** Quasiparticle energies obtained by root search *** '
|
||||
write(*,*)
|
||||
write(*,*) "ONLY LINEARISATION IMPLEMENTED YET"
|
||||
!call cRGF2_QP_graph(eta,nOrb,nC,nO,nV,nR,eHF,e_cap,ERI,Re_eGFlin,Im_eGFlin,eHF,e_cap,Re_eGF,Im_eGF,Re_Z,Im_Z)
|
||||
call complex_cRGF2_QP_graph(eta,nOrb,nC,nO,nV,nR,Re_eHF,Im_eHF,ERI,Re_eGFlin,Im_eGFlin,Re_eHF,Im_eHF,Re_eGF,Im_eGF,Re_Z,Im_Z)
|
||||
|
||||
end if
|
||||
|
||||
|
99
src/GF/complex_cRGF2_QP_graph.f90
Normal file
99
src/GF/complex_cRGF2_QP_graph.f90
Normal file
@ -0,0 +1,99 @@
|
||||
subroutine complex_cRGF2_QP_graph(eta,nBas,nC,nO,nV,nR,Re_eHF,Im_eHF,&
|
||||
ERI,Re_eGFlin,Im_eGFlin,Re_eOld,Im_eold,Re_eGF,Im_eGF,Re_Z,Im_Z)
|
||||
|
||||
! 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
|
||||
integer,intent(in) :: nC
|
||||
integer,intent(in) :: nO
|
||||
integer,intent(in) :: nV
|
||||
integer,intent(in) :: nR
|
||||
double precision,intent(in) :: Re_eHF(nBas)
|
||||
double precision,intent(in) :: Im_eHF(nBas)
|
||||
double precision,intent(in) :: Re_eGFlin(nBas)
|
||||
double precision,intent(in) :: Im_eGFlin(nBas)
|
||||
double precision,intent(in) :: Re_eOld(nBas)
|
||||
double precision,intent(in) :: Im_eOld(nBas)
|
||||
complex*16,intent(in) :: ERI(nBas,nBas,nBas,nBas)
|
||||
|
||||
! Local variables
|
||||
|
||||
integer :: p
|
||||
integer :: nIt
|
||||
integer,parameter :: maxIt = 64
|
||||
double precision,parameter :: thresh = 1d-6
|
||||
double precision :: Re_SigC,Im_SigC,Re_dSigC,Im_dSigC
|
||||
double precision :: Re_f,Im_f,Re_df,Im_df
|
||||
double precision :: Re_w,Im_w
|
||||
|
||||
! Output variables
|
||||
|
||||
double precision,intent(out) :: Re_eGF(nBas),Im_eGF(nBas)
|
||||
double precision,intent(out) :: Re_Z(nBas),Im_Z(nBas)
|
||||
|
||||
! Run Newton's algorithm to find the root
|
||||
|
||||
write(*,*)'-----------------------------------------------------'
|
||||
write(*,'(A5,1X,A3,1X,A15,1X,A15,1X,A10)') 'Orb.','It.','Re(e_GFlin) (eV)','Re(e_GF) (eV)','Re(Z)'
|
||||
write(*,'(A5,1X,A3,1X,A15,1X,A15,1X,A10)') 'Orb.','It.','Im(e_GFlin) (eV)','Im(e_GF) (eV)','Im(Z)'
|
||||
write(*,*)'-----------------------------------------------------'
|
||||
|
||||
Re_SigC = 0d0
|
||||
Im_SigC = 0d0
|
||||
Re_dSigC = 0d0
|
||||
Im_dSigC = 0d0
|
||||
|
||||
do p=nC+1,nBas-nR
|
||||
|
||||
Re_w = Re_eGFlin(p)
|
||||
Im_w = Im_eGFlin(p)
|
||||
nIt = 0
|
||||
Re_f = 1d0
|
||||
Im_f = 0d0
|
||||
|
||||
do while (abs(cmplx(Re_f,Im_f,kind=8)) > thresh .and. nIt < maxIt)
|
||||
|
||||
nIt = nIt + 1
|
||||
|
||||
call complex_cRGF_SigC_dSigC(p,eta,nBas,nO,nV,nR,Re_w,Im_w,Re_eOld,Im_eOld,ERI,&
|
||||
Re_SigC,Im_SigC,Re_Z,Im_Z)
|
||||
|
||||
Re_f = Re_w - Re_eHF(p) - Re_SigC
|
||||
Im_f = Im_w - Im_eHF(p) - Im_SigC
|
||||
Re_df = (1d0 - Re_dSigC)/((1d0 - Re_dSigC)**2 + Im_dSigC**2)
|
||||
Im_df = Im_dSigC/((1d0 - Re_dSigC)**2 + Im_dSigC**2)
|
||||
|
||||
Re_w = Re_w - Re_df*Re_f + Im_df*Im_f
|
||||
Im_w = Im_w - Re_f*Im_df - Re_df*Im_f
|
||||
|
||||
end do
|
||||
|
||||
if(nIt == maxIt) then
|
||||
|
||||
Re_eGF(p) = Re_eGFlin(p)
|
||||
Im_eGF(p) = Im_eGFlin(p)
|
||||
write(*,'(I5,1X,I3,1X,F15.9,1X,F15.9,1X,F10.6,1X,A12)') p,nIt,Re_eGFlin(p)*HaToeV,Re_eGF(p)*HaToeV,Re_Z(p),'Cvg Failed!'
|
||||
write(*,'(I5,1X,I3,1X,F15.9,1X,F15.9,1X,F10.6,1X,A12)') p,nIt,Im_eGFlin(p)*HaToeV,Im_eGF(p)*HaToeV,Im_Z(p),'Cvg Failed!'
|
||||
|
||||
else
|
||||
|
||||
Re_eGF(p) = Re_w
|
||||
Im_eGF(p) = Im_w
|
||||
Re_Z(p) = Re_df
|
||||
Im_Z(p) = Im_df
|
||||
|
||||
write(*,'(I5,1X,I3,1X,F15.9,1X,F15.9,1X,F10.6)') p,nIt,Re_eGFlin(p)*HaToeV,Re_eGF(p)*HaToeV,Re_Z(p)
|
||||
write(*,'(I5,1X,I3,1X,F15.9,1X,F15.9,1X,F10.6)') p,nIt,Im_eGFlin(p)*HaToeV,Im_eGF(p)*HaToeV,Im_Z(p)
|
||||
|
||||
write(*,*)'-----------------------------------------------------'
|
||||
end if
|
||||
|
||||
end do
|
||||
|
||||
end subroutine
|
84
src/GF/complex_cRGF_SigC_dSigC.f90
Normal file
84
src/GF/complex_cRGF_SigC_dSigC.f90
Normal file
@ -0,0 +1,84 @@
|
||||
subroutine complex_cRGF_SigC_dSigC(p,eta,nBas,nC,nO,nV,nR,Re_w,Im_w,Re_e,Im_e,ERI,Re_SigC,Im_SigC,Re_DS,Im_DS)
|
||||
|
||||
|
||||
implicit none
|
||||
include 'parameters.h'
|
||||
|
||||
! Input variables
|
||||
|
||||
integer,intent(in) :: p
|
||||
double precision,intent(in) :: eta
|
||||
integer,intent(in) :: nBas
|
||||
integer,intent(in) :: nC
|
||||
integer,intent(in) :: nO
|
||||
integer,intent(in) :: nV
|
||||
integer,intent(in) :: nR
|
||||
double precision,intent(in) :: Re_e(nBas)
|
||||
double precision,intent(in) :: Im_e(nBas)
|
||||
double precision,intent(in) :: Re_w
|
||||
double precision,intent(in) :: Im_w
|
||||
complex*16,intent(in) :: ERI(nBas,nBas,nBas,nBas)
|
||||
|
||||
! Local variables
|
||||
|
||||
integer :: i,j,a,b
|
||||
double precision :: eps
|
||||
double precision :: eta_tilde
|
||||
complex*16 :: num
|
||||
complex*16 :: z_dummy
|
||||
|
||||
! Output variables
|
||||
|
||||
double precision,intent(out) :: Re_SigC
|
||||
double precision,intent(out) :: Im_SigC
|
||||
double precision,intent(out) :: Re_DS
|
||||
double precision,intent(out) :: Im_DS
|
||||
|
||||
! Initialize
|
||||
Re_SigC = 0d0
|
||||
Im_SigC = 0d0
|
||||
Re_DS = 0d0
|
||||
Im_DS = 0d0
|
||||
|
||||
|
||||
! Compute GF2 self-energy
|
||||
|
||||
do i=nC+1,nO
|
||||
do j=nC+1,nO
|
||||
do a=nO+1,nBas-nR
|
||||
|
||||
eps = Re_w + Re_e(a) - Re_e(i) - Re_e(j)
|
||||
eta_tilde = eta - Im_w + Im_e(i) + Im_e(a) - Im_e(j)
|
||||
num = (2d0*ERI(p,a,i,j) - ERI(p,a,j,i))*ERI(p,a,i,j)
|
||||
z_dummy = num*cmplx(eps/(eps**2 + eta_tilde**2),eta_tilde/(eps**2 + eta_tilde**2),kind=8)
|
||||
Re_SigC = Re_SigC + real(z_dummy)
|
||||
Im_SigC = Im_SigC + aimag(z_dummy)
|
||||
z_dummy = num*cmplx(-(eps**2 - eta_tilde**2)/(eps**2 + eta_tilde**2)**2,&
|
||||
-2*eta_tilde*eps/(eps**2 + eta_tilde**2)**2,kind=8)
|
||||
Re_DS = Re_DS + real(z_dummy)
|
||||
Im_DS = Im_DS + aimag(z_dummy)
|
||||
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
|
||||
do i=nC+1,nO
|
||||
do a=nO+1,nBas-nR
|
||||
do b=nO+1,nBas-nR
|
||||
|
||||
eps = Re_w + Re_e(i) - Re_e(a) - Re_e(b)
|
||||
eta_tilde = eta + Im_w - Im_e(a) - Im_e(b) + Im_e(i)
|
||||
num = (2d0*ERI(p,i,a,b) - ERI(p,i,b,a))*ERI(p,i,a,b)
|
||||
|
||||
z_dummy = num*cmplx(eps/(eps**2 + eta_tilde**2),-eta_tilde/(eps**2 + eta_tilde**2),kind=8)
|
||||
Re_SigC = Re_SigC + real(z_dummy)
|
||||
Im_SigC = Im_SigC + aimag(z_dummy)
|
||||
z_dummy = num*cmplx(-(eps**2 - eta_tilde**2)/(eps**2 + eta_tilde**2)**2,&
|
||||
2*eta_tilde*eps/(eps**2 + eta_tilde**2)**2,kind=8)
|
||||
Re_DS = Re_DS + real(z_dummy)
|
||||
Im_DS = Im_DS + aimag(z_dummy)
|
||||
|
||||
end do
|
||||
end do
|
||||
end do
|
||||
end subroutine
|
Loading…
x
Reference in New Issue
Block a user