4
1
mirror of https://github.com/pfloos/quack synced 2024-06-18 11:15:30 +02:00

EE-EOM-CCD

This commit is contained in:
Pierre-Francois Loos 2022-09-27 16:43:24 +02:00
parent 1f2f6634bb
commit 69bab74d81
4 changed files with 152 additions and 10 deletions

View File

@ -9,7 +9,7 @@
# CIS* CIS(D) CID CISD FCI
F F F F F
# RPA* RPAx* crRPA ppRPA
F F F T
F T F F
# G0F2* evGF2* qsGF2* G0F3 evGF3
F F F F F
# G0W0* evGW* qsGW* ufG0W0 ufGW

View File

@ -44,6 +44,8 @@ subroutine CCD(BSE,maxSCF,thresh,max_diis,nBasin,nCin,nOin,nVin,nRin,ERI,ENuc,ER
double precision,allocatable :: OVOV(:,:,:,:)
double precision,allocatable :: VVVV(:,:,:,:)
double precision,allocatable :: OVVO(:,:,:,:)
double precision,allocatable :: X1(:,:,:,:)
double precision,allocatable :: X2(:,:)
double precision,allocatable :: X3(:,:)
@ -60,11 +62,11 @@ subroutine CCD(BSE,maxSCF,thresh,max_diis,nBasin,nCin,nOin,nVin,nRin,ERI,ENuc,ER
double precision,allocatable :: error_diis(:,:)
double precision,allocatable :: t_diis(:,:)
logical :: do_EE_EOM_CC_1h1p = .false.
logical :: do_EE_EOM_CC_1h1p = .true.
logical :: do_EA_EOM_CC_1p = .false.
logical :: do_IP_EOM_CC_1h = .false.
logical :: do_DEA_EOM_CC_2p = .true.
logical :: do_DIP_EOM_CC_2h = .true.
logical :: do_DEA_EOM_CC_2p = .false.
logical :: do_DIP_EOM_CC_2h = .false.
! Hello world
@ -125,6 +127,9 @@ subroutine CCD(BSE,maxSCF,thresh,max_diis,nBasin,nCin,nOin,nVin,nRin,ERI,ENuc,ER
OVOV(:,:,:,:) = dbERI(nC+1:nO ,nO+1:nBas-nR,nC+1:nO ,nO+1:nBas-nR)
VVVV(:,:,:,:) = dbERI(nO+1:nBas-nR,nO+1:nBas-nR,nO+1:nBas-nR,nO+1:nBas-nR)
allocate(OVVO(nO-nC,nV-nR,nV-nR,nO-nC))
OVVO(:,:,:,:) = dbERI(nC+1:nO,nO+1:nBas-nR,nO+1:nBas-nR,nC+1:nO)
deallocate(dbERI)
! MP2 guess amplitudes
@ -258,7 +263,7 @@ subroutine CCD(BSE,maxSCF,thresh,max_diis,nBasin,nCin,nOin,nVin,nRin,ERI,ENuc,ER
! EE-EOM-CCD (1h1p)
! if(do_EE-EOM-CC_1h1p) call EE-EOM-CCD_1h1p()
if(do_EE_EOM_CC_1h1p) call EE_EOM_CCD_1h1p(nC,nO,nV,nR,eO,eV,OOVV,OVVO,t)
! EA-EOM (1p)

View File

@ -45,14 +45,14 @@ subroutine DIP_EOM_CCD_2h(nC,nO,nV,nR,eO,OOVV,OOOO,t)
! Form one-body terms
do i=1,nO-nR
do j=1,nO-nR
do i=1,nO-nC
do j=1,nO-nC
F(i,j) = eO(i)*Kronecker_delta(i,j)
do k=1,nO-nR
do a=1,nV-nC
do b=1,nV-nC
do k=1,nO-nC
do a=1,nV-nR
do b=1,nV-nR
F(i,j) = F(i,j) + 0.5d0*OOVV(i,k,a,b)*t(j,k,a,b)

137
src/CC/EE_EOM_CCD_1h1p.f90 Normal file
View File

@ -0,0 +1,137 @@
subroutine EE_EOM_CCD_1h1p(nC,nO,nV,nR,eO,eV,OOVV,OVVO,t)
! EE-EOM-CCD calculation up to 1h1p
implicit none
! Input variables
integer,intent(in) :: nC
integer,intent(in) :: nO
integer,intent(in) :: nV
integer,intent(in) :: nR
double precision,intent(in) :: eO(nO)
double precision,intent(in) :: eV(nV)
double precision,intent(in) :: OOVV(nO,nO,nV,nV)
double precision,intent(in) :: OVVO(nO,nV,nV,nO)
double precision,intent(in) :: t(nO,nO,nV,nV)
! Local variables
integer :: a,b,c,d
integer :: i,j,k,l
integer :: ia,jb
integer :: nS
double precision,external :: Kronecker_delta
double precision,allocatable :: Fvv(:,:)
double precision,allocatable :: Foo(:,:)
double precision,allocatable :: Wovvo(:,:,:,:)
double precision,allocatable :: H(:,:)
double precision,allocatable :: Om(:)
! Hello world
write(*,*)
write(*,*)'*********************'
write(*,*)'| EE-EOM-CCD (1h1p) |'
write(*,*)'*********************'
write(*,*)
! Size of the EOM Hamiltonian
nS = (nO-nC)*(nV-nR)
! Memory allocation
allocate(Foo(nO,nO),Fvv(nV,nV),Wovvo(nO,nV,nV,nO),H(nS,nS),Om(nS))
! Form one-body terms
do a=1,nV-nR
do b=1,nV-nR
Fvv(a,b) = eV(a)*Kronecker_delta(a,b)
do i=1,nO-nC
do j=1,nO-nC
do c=1,nV-nR
Fvv(a,b) = Fvv(a,b) - 0.5d0*OOVV(i,j,b,c)*t(i,j,a,c)
end do
end do
end do
end do
end do
do i=1,nO-nC
do j=1,nO-nC
Foo(i,j) = eO(i)*Kronecker_delta(i,j)
do k=1,nO-nC
do a=1,nV-nR
do b=1,nV-nR
Foo(i,j) = Foo(i,j) + 0.5d0*OOVV(i,k,a,b)*t(j,k,a,b)
end do
end do
end do
end do
end do
! Form two-body terms
do i=1,nO-nC
do b=1,nV-nR
do a=1,nV-nR
do j=1,nO-nC
Wovvo(i,b,a,j) = OVVO(i,b,a,j)
do k=1,nO-nC
do c=1,nV-nR
Wovvo(i,b,a,j) = Wovvo(i,b,a,j) + OOVV(i,k,a,c)*t(k,j,c,b)
end do
end do
end do
end do
end do
end do
! Form EOM Hamiltonian
ia = 0
do i=1,nO-nC
do a=1,nV-nR
ia = ia + 1
jb = 0
do j=1,nO-nC
do b=1,nV-nR
jb = jb + 1
H(ia,jb) = Fvv(a,b)*Kronecker_delta(i,j) - Kronecker_delta(a,b)*Foo(i,j) + Wovvo(i,b,a,j)
end do
end do
end do
end do
! Diagonalize EOM Hamiltonian
if(nS > 0) call diagonalize_matrix(nS,H,Om)
! Dump results
call print_excitation('EE-EOM-CCD ',3,nS,Om)
end subroutine EE_EOM_CCD_1h1p