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

173 lines
5.0 KiB
Fortran
Raw Permalink Normal View History

subroutine RCC(dotest,doCCD,dopCCD,doDCD,doCCSD,doCCSDT,dodrCCD,dorCCD,docrCCD,dolCCD, &
2023-11-13 22:15:00 +01:00
maxSCF,thresh,max_diis,nBas,nC,nO,nV,nR,ERI,ENuc,ERHF,eHF)
2023-07-22 22:19:46 +02:00
! Coupled-cluster module
implicit none
include 'parameters.h'
! Input variables
logical :: dotest
2023-07-22 22:19:46 +02:00
logical :: doCCD
logical :: dopCCD
logical :: doDCD
logical :: doCCSD
logical :: doCCSDT
2023-10-02 21:37:25 +02:00
logical :: dodrCCD
logical :: dorCCD
logical :: docrCCD
logical :: dolCCD
2023-07-22 22:19:46 +02:00
integer,intent(in) :: maxSCF
integer,intent(in) :: max_diis
double precision,intent(in) :: thresh
integer,intent(in) :: nBas
2023-11-29 16:20:53 +01:00
integer,intent(in) :: nC
integer,intent(in) :: nO
integer,intent(in) :: nV
integer,intent(in) :: nR
2023-07-22 22:19:46 +02:00
double precision,intent(in) :: ENuc
2023-11-13 22:15:00 +01:00
double precision,intent(in) :: ERHF
double precision,intent(in) :: eHF(nBas)
2023-07-22 22:19:46 +02:00
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
! Local variables
double precision :: start_CC ,end_CC ,t_CC
!------------------------------------------------------------------------
! Perform CCD calculation
!------------------------------------------------------------------------
if(doCCD) then
2023-07-31 16:01:12 +02:00
call wall_time(start_CC)
2023-11-13 22:15:00 +01:00
call CCD(dotest,maxSCF,thresh,max_diis,nBas,nC,nO,nV,nR,ERI,ENuc,ERHF,eHF)
2023-07-31 16:01:12 +02:00
call wall_time(end_CC)
2023-07-22 22:19:46 +02:00
t_CC = end_CC - start_CC
write(*,'(A65,1X,F9.3,A8)') 'Total CPU time for CCD = ',t_CC,' seconds'
write(*,*)
end if
!------------------------------------------------------------------------
! Perform DCD calculation
!------------------------------------------------------------------------
if(doDCD) then
2023-07-31 16:01:12 +02:00
call wall_time(start_CC)
call DCD(dotest,maxSCF,thresh,max_diis,nBas,nC,nO,nV,nR, &
2023-11-13 22:15:00 +01:00
ERI,ENuc,ERHF,eHF)
2023-07-31 16:01:12 +02:00
call wall_time(end_CC)
2023-07-22 22:19:46 +02:00
t_CC = end_CC - start_CC
write(*,'(A65,1X,F9.3,A8)') 'Total CPU time for DCD = ',t_CC,' seconds'
write(*,*)
end if
!------------------------------------------------------------------------
! Perform CCSD or CCSD(T) calculation
!------------------------------------------------------------------------
if(doCCSDT) doCCSD = .true.
if(doCCSD) then
2023-07-31 16:01:12 +02:00
call wall_time(start_CC)
2023-11-13 22:15:00 +01:00
call CCSD(dotest,maxSCF,thresh,max_diis,doCCSDT,nBas,nC,nO,nV,nR,ERI,ENuc,ERHF,eHF)
2023-07-31 16:01:12 +02:00
call wall_time(end_CC)
2023-07-22 22:19:46 +02:00
t_CC = end_CC - start_CC
write(*,'(A65,1X,F9.3,A8)') 'Total CPU time for CCSD or CCSD(T)= ',t_CC,' seconds'
write(*,*)
end if
!------------------------------------------------------------------------
! Perform direct ring CCD calculation
!------------------------------------------------------------------------
2023-10-02 21:37:25 +02:00
if(dodrCCD) then
2023-07-22 22:19:46 +02:00
2023-07-31 16:01:12 +02:00
call wall_time(start_CC)
2023-11-13 22:15:00 +01:00
call drCCD(dotest,maxSCF,thresh,max_diis,nBas,nC,nO,nV,nR,ERI,ENuc,ERHF,eHF)
2023-07-31 16:01:12 +02:00
call wall_time(end_CC)
2023-07-22 22:19:46 +02:00
t_CC = end_CC - start_CC
write(*,'(A65,1X,F9.3,A8)') 'Total CPU time for direct ring CCD = ',t_CC,' seconds'
write(*,*)
end if
!------------------------------------------------------------------------
! Perform ring CCD calculation
!------------------------------------------------------------------------
2023-10-02 21:37:25 +02:00
if(dorCCD) then
2023-07-22 22:19:46 +02:00
2023-07-31 16:01:12 +02:00
call wall_time(start_CC)
2023-11-13 22:15:00 +01:00
call rCCD(dotest,maxSCF,thresh,max_diis,nBas,nC,nO,nV,nR,ERI,ENuc,ERHF,eHF)
2023-07-31 16:01:12 +02:00
call wall_time(end_CC)
2023-07-22 22:19:46 +02:00
t_CC = end_CC - start_CC
write(*,'(A65,1X,F9.3,A8)') 'Total CPU time for rCCD = ',t_CC,' seconds'
write(*,*)
end if
!------------------------------------------------------------------------
! Perform crossed-ring CCD calculation
!------------------------------------------------------------------------
2023-10-02 21:37:25 +02:00
if(docrCCD) then
2023-07-22 22:19:46 +02:00
2023-07-31 16:01:12 +02:00
call wall_time(start_CC)
2023-11-13 22:15:00 +01:00
call crCCD(dotest,maxSCF,thresh,max_diis,nBas,nC,nO,nV,nR,ERI,ENuc,ERHF,eHF)
2023-07-31 16:01:12 +02:00
call wall_time(end_CC)
2023-07-22 22:19:46 +02:00
t_CC = end_CC - start_CC
write(*,'(A65,1X,F9.3,A8)') 'Total CPU time for crossed-ring CCD = ',t_CC,' seconds'
write(*,*)
end if
!------------------------------------------------------------------------
! Perform ladder CCD calculation
!------------------------------------------------------------------------
2023-10-02 21:37:25 +02:00
if(dolCCD) then
2023-07-22 22:19:46 +02:00
2023-07-31 16:01:12 +02:00
call wall_time(start_CC)
2023-11-13 22:15:00 +01:00
call lCCD(dotest,maxSCF,thresh,max_diis,nBas,nC,nO,nV,nR,ERI,ENuc,ERHF,eHF)
2023-07-31 16:01:12 +02:00
call wall_time(end_CC)
2023-07-22 22:19:46 +02:00
t_CC = end_CC - start_CC
write(*,'(A65,1X,F9.3,A8)') 'Total CPU time for ladder CCD = ',t_CC,' seconds'
write(*,*)
end if
!------------------------------------------------------------------------
! Perform pair CCD calculation
!------------------------------------------------------------------------
if(dopCCD) then
2023-07-31 16:01:12 +02:00
call wall_time(start_CC)
2023-11-13 22:15:00 +01:00
call pCCD(dotest,maxSCF,thresh,max_diis,nBas,nC,nO,nV,nR,ERI,ENuc,ERHF,eHF)
2023-07-31 16:01:12 +02:00
call wall_time(end_CC)
2023-07-22 22:19:46 +02:00
t_CC = end_CC - start_CC
write(*,'(A65,1X,F9.3,A8)') 'Total CPU time for pair CCD = ',t_CC,' seconds'
write(*,*)
end if
end subroutine