4
1
mirror of https://github.com/pfloos/quack synced 2024-06-18 03:05:31 +02:00
quack/src/MBPT/G0F2.f90

123 lines
2.9 KiB
Fortran
Raw Normal View History

2021-03-05 22:34:48 +01:00
subroutine G0F2(BSE,TDA,dBSE,dTDA,evDyn,singlet,triplet,linearize,eta,nBas,nC,nO,nV,nR,nS,ENuc,ERHF,ERI,dipole_int,eHF)
2020-03-19 10:21:18 +01:00
! Perform a one-shot second-order Green function calculation
implicit none
include 'parameters.h'
! Input variables
2020-06-01 17:14:42 +02:00
logical,intent(in) :: BSE
logical,intent(in) :: TDA
2020-06-15 23:04:07 +02:00
logical,intent(in) :: dBSE
logical,intent(in) :: dTDA
2020-06-16 14:02:14 +02:00
logical,intent(in) :: evDyn
2021-03-05 22:34:48 +01:00
logical,intent(in) :: singlet
logical,intent(in) :: triplet
2020-03-19 10:21:18 +01:00
logical,intent(in) :: linearize
2020-06-01 17:14:42 +02:00
double precision,intent(in) :: eta
2020-03-19 10:21:18 +01:00
integer,intent(in) :: nBas
integer,intent(in) :: nO
integer,intent(in) :: nC
integer,intent(in) :: nV
integer,intent(in) :: nR
2020-06-01 17:14:42 +02:00
integer,intent(in) :: nS
double precision,intent(in) :: ENuc
double precision,intent(in) :: ERHF
double precision,intent(in) :: eHF(nBas)
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
double precision,intent(in) :: dipole_int(nBas,nBas,ncart)
2020-03-19 10:21:18 +01:00
! Local variables
double precision :: eps
2020-06-01 17:14:42 +02:00
double precision :: V
double precision :: EcBSE(nspin)
2020-03-19 10:21:18 +01:00
double precision,allocatable :: eGF2(:)
2020-03-21 16:31:39 +01:00
double precision,allocatable :: Sig(:)
2020-03-19 10:21:18 +01:00
double precision,allocatable :: Z(:)
integer :: i,j,a,b,p
! Hello world
write(*,*)
write(*,*)'************************************************'
write(*,*)'| One-shot second-order Green function |'
write(*,*)'************************************************'
write(*,*)
! Memory allocation
2020-03-21 16:31:39 +01:00
allocate(Sig(nBas),Z(nBas),eGF2(nBas))
2020-03-19 10:21:18 +01:00
2020-03-21 16:31:39 +01:00
if(linearize) then
write(*,*) '*** Quasiparticle equation will be linearized ***'
write(*,*)
2020-03-19 10:21:18 +01:00
2020-03-21 16:31:39 +01:00
end if
2020-03-19 10:21:18 +01:00
2020-03-21 16:31:39 +01:00
! Frequency-dependent second-order contribution
Sig(:) = 0d0
Z(:) = 0d0
2020-03-19 10:21:18 +01:00
do p=nC+1,nBas-nR
do i=nC+1,nO
do j=nC+1,nO
do a=nO+1,nBas-nR
2020-06-01 17:14:42 +02:00
eps = eHF(p) + eHF(a) - eHF(i) - eHF(j)
V = (2d0*ERI(p,a,i,j) - ERI(p,a,j,i))*ERI(p,a,i,j)
2020-06-03 12:24:38 +02:00
2020-06-03 12:06:16 +02:00
Sig(p) = Sig(p) + V*eps/(eps**2 + eta**2)
Z(p) = Z(p) - V*(eps**2 - eta**2)/(eps**2 + eta**2)**2
2020-03-19 10:21:18 +01:00
end do
end do
end do
end do
do p=nC+1,nBas-nR
do i=nC+1,nO
do a=nO+1,nBas-nR
do b=nO+1,nBas-nR
2020-06-01 17:14:42 +02:00
eps = eHF(p) + eHF(i) - eHF(a) - eHF(b)
V = (2d0*ERI(p,i,a,b) - ERI(p,i,b,a))*ERI(p,i,a,b)
2020-06-03 12:24:38 +02:00
2020-06-03 12:06:16 +02:00
Sig(p) = Sig(p) + V*eps/(eps**2 + eta**2)
Z(p) = Z(p) - V*(eps**2 - eta**2)/(eps**2 + eta**2)**2
2020-03-19 10:21:18 +01:00
end do
end do
end do
end do
2020-06-03 12:06:16 +02:00
Z(:) = 1d0/(1d0 - Z(:))
2020-03-19 10:21:18 +01:00
if(linearize) then
2020-06-01 17:14:42 +02:00
eGF2(:) = eHF(:) + Z(:)*Sig(:)
2020-03-19 10:21:18 +01:00
else
2020-06-01 17:14:42 +02:00
eGF2(:) = eHF(:) + Sig(:)
2020-03-19 10:21:18 +01:00
end if
2020-03-21 16:31:39 +01:00
2020-03-19 10:21:18 +01:00
! Print results
2020-06-01 17:14:42 +02:00
call print_G0F2(nBas,nO,eHF,Sig,eGF2,Z)
2020-03-19 10:21:18 +01:00
2020-06-01 17:14:42 +02:00
! Perform BSE2 calculation
if(BSE) then
2021-03-05 22:34:48 +01:00
call BSE2(TDA,dBSE,dTDA,evDyn,singlet,triplet,eta,nBas,nC,nO,nV,nR,nS,ERI,dipole_int,eHF,eGF2,EcBSE)
2020-06-01 17:14:42 +02:00
end if
2020-06-03 12:24:38 +02:00
2020-03-19 10:21:18 +01:00
end subroutine G0F2