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

286 lines
7.9 KiB
Fortran
Raw Permalink Normal View History

subroutine CCD(dotest,maxSCF,thresh,max_diis,nBasin,nCin,nOin,nVin,nRin,ERI,ENuc,ERHF,eHF)
2019-03-19 10:13:33 +01:00
! CCD module
implicit none
! Input variables
logical,intent(in) :: dotest
2019-03-19 10:13:33 +01:00
integer,intent(in) :: maxSCF
integer,intent(in) :: max_diis
double precision,intent(in) :: thresh
2020-03-26 16:57:46 +01:00
integer,intent(in) :: nBasin
integer,intent(in) :: nCin
integer,intent(in) :: nOin
integer,intent(in) :: nVin
integer,intent(in) :: nRin
2019-03-19 10:13:33 +01:00
double precision,intent(in) :: ENuc,ERHF
2020-03-26 16:57:46 +01:00
double precision,intent(in) :: eHF(nBasin)
double precision,intent(in) :: ERI(nBasin,nBasin,nBasin,nBasin)
2019-03-19 10:13:33 +01:00
! Local variables
2020-03-26 16:57:46 +01:00
integer :: nBas
integer :: nC
2019-03-19 10:13:33 +01:00
integer :: nO
integer :: nV
2020-03-26 16:57:46 +01:00
integer :: nR
2019-03-19 10:13:33 +01:00
integer :: nSCF
double precision :: Conv
double precision :: EcMP2,EcMP3,EcMP4
double precision :: ECC,EcCC
2019-03-19 10:13:33 +01:00
double precision,allocatable :: seHF(:)
double precision,allocatable :: sERI(:,:,:,:)
double precision,allocatable :: dbERI(:,:,:,:)
double precision,allocatable :: eO(:)
double precision,allocatable :: eV(:)
double precision,allocatable :: delta_OOVV(:,:,:,:)
double precision,allocatable :: OOOO(:,:,:,:)
double precision,allocatable :: OOVV(:,:,:,:)
double precision,allocatable :: OVOV(:,:,:,:)
double precision,allocatable :: VVVV(:,:,:,:)
2022-09-27 16:43:24 +02:00
double precision,allocatable :: OVVO(:,:,:,:)
2019-03-19 10:13:33 +01:00
double precision,allocatable :: X1(:,:,:,:)
double precision,allocatable :: X2(:,:)
double precision,allocatable :: X3(:,:)
double precision,allocatable :: X4(:,:,:,:)
double precision,allocatable :: u(:,:,:,:)
double precision,allocatable :: v(:,:,:,:)
2022-09-27 14:07:31 +02:00
double precision,allocatable :: r(:,:,:,:)
double precision,allocatable :: t(:,:,:,:)
2019-03-19 10:13:33 +01:00
2023-07-30 22:14:28 +02:00
integer :: n_diis
2020-03-22 23:01:36 +01:00
double precision :: rcond
double precision,allocatable :: error_diis(:,:)
double precision,allocatable :: t_diis(:,:)
2022-09-28 16:09:09 +02:00
logical :: do_EE_EOM_CC_1h1p = .false.
2022-09-27 14:07:31 +02:00
logical :: do_EA_EOM_CC_1p = .false.
logical :: do_IP_EOM_CC_1h = .false.
2022-09-27 16:43:24 +02:00
logical :: do_DEA_EOM_CC_2p = .false.
logical :: do_DIP_EOM_CC_2h = .false.
2022-09-27 14:07:31 +02:00
2019-03-19 10:13:33 +01:00
! Hello world
write(*,*)
write(*,*)'**************************************'
write(*,*)'| CCD calculation |'
write(*,*)'**************************************'
write(*,*)
! Spatial to spin orbitals
2020-03-26 16:57:46 +01:00
nBas = 2*nBasin
nC = 2*nCin
nO = 2*nOin
nV = 2*nVin
nR = 2*nRin
2019-03-19 10:13:33 +01:00
2020-03-26 16:57:46 +01:00
allocate(seHF(nBas),sERI(nBas,nBas,nBas,nBas))
2019-03-19 10:13:33 +01:00
2020-03-26 16:57:46 +01:00
call spatial_to_spin_MO_energy(nBasin,eHF,nBas,seHF)
call spatial_to_spin_ERI(nBasin,ERI,nBas,sERI)
2019-03-19 10:13:33 +01:00
! Antysymmetrize ERIs
2020-03-26 16:57:46 +01:00
allocate(dbERI(nBas,nBas,nBas,nBas))
2019-03-19 10:13:33 +01:00
2023-07-19 11:03:55 +02:00
call antisymmetrize_ERI(2,nBas,sERI,dbERI)
2019-03-19 10:13:33 +01:00
deallocate(sERI)
! Form energy denominator
2020-10-19 14:20:35 +02:00
allocate(eO(nO-nC),eV(nV-nR))
allocate(delta_OOVV(nO-nC,nO-nC,nV-nR,nV-nR))
2019-03-19 10:13:33 +01:00
2020-10-19 14:20:35 +02:00
eO(:) = seHF(nC+1:nO)
eV(:) = seHF(nO+1:nBas-nR)
2019-03-19 10:13:33 +01:00
2020-03-26 16:57:46 +01:00
call form_delta_OOVV(nC,nO,nV,nR,eO,eV,delta_OOVV)
2019-03-19 10:13:33 +01:00
deallocate(seHF)
! Create integral batches
2020-10-19 14:20:35 +02:00
allocate(OOOO(nO-nC,nO-nC,nO-nC,nO-nC),OOVV(nO-nC,nO-nC,nV-nR,nV-nR), &
OVOV(nO-nC,nV-nR,nO-nC,nV-nR),VVVV(nV-nR,nV-nR,nV-nR,nV-nR))
2019-03-19 10:13:33 +01:00
2020-10-19 14:20:35 +02:00
OOOO(:,:,:,:) = dbERI(nC+1:nO ,nC+1:nO ,nC+1:nO ,nC+1:nO )
OOVV(:,:,:,:) = dbERI(nC+1:nO ,nC+1:nO ,nO+1:nBas-nR,nO+1:nBas-nR)
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)
2019-03-19 10:13:33 +01:00
2022-09-27 16:43:24 +02:00
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)
2019-03-19 10:13:33 +01:00
deallocate(dbERI)
! MP2 guess amplitudes
2022-09-27 14:07:31 +02:00
allocate(t(nO-nC,nO-nC,nV-nR,nV-nR))
2019-03-19 10:13:33 +01:00
2022-09-27 14:07:31 +02:00
t(:,:,:,:) = -OOVV(:,:,:,:)/delta_OOVV(:,:,:,:)
2019-03-19 10:13:33 +01:00
2022-09-27 14:07:31 +02:00
call CCD_correlation_energy(nC,nO,nV,nR,OOVV,t,EcMP2)
2019-03-19 10:13:33 +01:00
EcMP4 = 0d0
2020-03-22 23:01:36 +01:00
! Memory allocation for DIIS
2020-10-19 14:20:35 +02:00
allocate(error_diis((nO-nR)**2*(nV-nR)**2,max_diis),t_diis((nO-nR)**2*(nV-nR)**2,max_diis))
2020-03-22 23:01:36 +01:00
2019-03-19 10:13:33 +01:00
! Initialization
2022-09-27 14:07:31 +02:00
allocate(r(nO-nC,nO-nC,nV-nR,nV-nR),u(nO-nC,nO-nC,nV-nR,nV-nR),v(nO-nC,nO-nC,nV-nR,nV-nR))
2020-10-19 14:20:35 +02:00
allocate(X1(nO-nC,nO-nC,nO-nC,nO-nC),X2(nV-nR,nV-nR),X3(nO-nC,nO-nC),X4(nO-nC,nO-nC,nV-nR,nV-nR))
2019-03-19 10:13:33 +01:00
Conv = 1d0
nSCF = 0
2020-03-22 23:01:36 +01:00
n_diis = 0
t_diis(:,:) = 0d0
error_diis(:,:) = 0d0
2019-03-19 10:13:33 +01:00
!------------------------------------------------------------------------
! Main SCF loop
!------------------------------------------------------------------------
write(*,*)
write(*,*)'----------------------------------------------------'
write(*,*)'| CCD calculation |'
write(*,*)'----------------------------------------------------'
write(*,'(1X,A1,1X,A3,1X,A1,1X,A16,1X,A1,1X,A10,1X,A1,1X,A10,1X,A1,1X)') &
'|','#','|','E(CCD)','|','Ec(CCD)','|','Conv','|'
write(*,*)'----------------------------------------------------'
do while(Conv > thresh .and. nSCF < maxSCF)
! Increment
nSCF = nSCF + 1
! Form linear array
2022-09-27 14:07:31 +02:00
call form_u(nC,nO,nV,nR,OOOO,VVVV,OVOV,t,u)
2019-03-19 10:13:33 +01:00
! Form interemediate arrays
2022-09-27 14:07:31 +02:00
call form_X(nC,nO,nV,nR,OOVV,t,X1,X2,X3,X4)
2019-03-19 10:13:33 +01:00
! Form quadratic array
2022-09-27 14:07:31 +02:00
call form_v(nC,nO,nV,nR,X1,X2,X3,X4,t,v)
2019-03-19 10:13:33 +01:00
! Compute residual
2022-09-27 14:07:31 +02:00
r(:,:,:,:) = OOVV(:,:,:,:) + delta_OOVV(:,:,:,:)*t(:,:,:,:) + u(:,:,:,:) + v(:,:,:,:)
2019-03-19 10:13:33 +01:00
! Check convergence
2022-09-27 14:07:31 +02:00
Conv = maxval(abs(r(:,:,:,:)))
2019-03-19 10:13:33 +01:00
! Update amplitudes
2022-09-27 14:07:31 +02:00
t(:,:,:,:) = t(:,:,:,:) - r(:,:,:,:)/delta_OOVV(:,:,:,:)
2019-03-19 10:13:33 +01:00
! Compute correlation energy
call CCD_correlation_energy(nC,nO,nV,nR,OOVV,t,EcCC)
2019-03-19 10:13:33 +01:00
2022-09-27 14:07:31 +02:00
if(nSCF == 1) call MP3_correlation_energy(nC,nO,nV,nR,OOVV,t,v,delta_OOVV,EcMP3)
2019-03-19 10:13:33 +01:00
! Dump results
ECC = ERHF + EcCC
2019-03-19 10:13:33 +01:00
2020-03-22 23:01:36 +01:00
! DIIS extrapolation
n_diis = min(n_diis+1,max_diis)
2022-09-27 14:07:31 +02:00
call DIIS_extrapolation(rcond,(nO-nC)**2*(nV-nR)**2,(nO-nC)**2*(nV-nR)**2,n_diis,error_diis,t_diis,-r/delta_OOVV,t)
2020-03-22 23:01:36 +01:00
! Reset DIIS if required
if(abs(rcond) < 1d-15) n_diis = 0
2019-03-19 10:13:33 +01:00
write(*,'(1X,A1,1X,I3,1X,A1,1X,F16.10,1X,A1,1X,F10.6,1X,A1,1X,F10.6,1X,A1,1X)') &
'|',nSCF,'|',ECC+ENuc,'|',EcCC,'|',Conv,'|'
2019-03-19 10:13:33 +01:00
2023-12-03 18:47:30 +01:00
end do
2019-03-19 10:13:33 +01:00
write(*,*)'----------------------------------------------------'
!------------------------------------------------------------------------
! End of SCF loop
!------------------------------------------------------------------------
! Did it actually converge?
if(nSCF == maxSCF) then
write(*,*)
write(*,*)'!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!'
write(*,*)' Convergence failed '
write(*,*)'!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!'
write(*,*)
stop
2023-12-03 18:47:30 +01:00
end if
2019-03-19 10:13:33 +01:00
2021-11-10 09:42:30 +01:00
write(*,*)
write(*,*)'----------------------------------------------------'
write(*,*)' CCD energy '
write(*,*)'----------------------------------------------------'
write(*,'(1X,A30,1X,F15.10)')' E(CCD) = ',ECC
write(*,'(1X,A30,1X,F15.10)')' Ec(CCD) = ',EcCC
2021-11-10 09:42:30 +01:00
write(*,*)'----------------------------------------------------'
write(*,*)
2019-03-19 10:13:33 +01:00
! Moller-Plesset energies
write(*,*)
write(*,'(1X,A15,1X,F10.6)') 'Ec(MP2) = ',EcMP2
write(*,'(1X,A15,1X,F10.6)') 'Ec(MP3) = ',EcMP3
write(*,'(1X,A15,1X,F10.6)') 'Ec(MP4-SDQ) = ',EcMP4
write(*,*)
2022-09-27 14:07:31 +02:00
!------------------------------------------------------------------------
! EOM section
!------------------------------------------------------------------------
! EE-EOM-CCD (1h1p)
2022-09-27 16:43:24 +02:00
if(do_EE_EOM_CC_1h1p) call EE_EOM_CCD_1h1p(nC,nO,nV,nR,eO,eV,OOVV,OVVO,t)
2022-09-27 14:07:31 +02:00
! EA-EOM (1p)
! if(do_EA-EOM-CC_1p) call EA-EOM-CCD_1p()
! IP-EOM-CCD(1h)
! if(do_IP-EOM-CC_1h) call IP-EOM-CCD_1h()
! DEA-EOM (2p)
if(do_DEA_EOM_CC_2p) call DEA_EOM_CCD_2p(nC,nO,nV,nR,eV,OOVV,VVVV,t)
! DIP-EOM-CCD(2h)
if(do_DIP_EOM_CC_2h) call DIP_EOM_CCD_2h(nC,nO,nV,nR,eO,OOVV,OOOO,t)
! Testing zone
if(dotest) then
call dump_test_value('R','CCD correlation energy',EcCC)
end if
2023-07-27 11:38:38 +02:00
end subroutine