1
0
mirror of https://gitlab.com/scemama/qp_plugins_scemama.git synced 2024-07-26 04:37:31 +02:00
qp_plugins_scemama/devel/cc/CCSD_Ec_nc.irp.f

53 lines
1.0 KiB
FortranFixed
Raw Normal View History

2019-09-14 14:15:25 +02:00
subroutine CCSD_Ec_nc(t1,t2,Fov,EcCCSD)
2019-09-11 17:09:30 +02:00
! Compute the CCSD correlatio energy in non-conanical form
implicit none
! Input variables
2019-09-14 14:15:25 +02:00
double precision,intent(in) :: t1(spin_occ_num,spin_vir_num)
double precision,intent(in) :: t2(spin_occ_num,spin_occ_num,spin_vir_num,spin_vir_num)
2019-09-11 17:09:30 +02:00
2019-09-14 14:15:25 +02:00
double precision,intent(in) :: Fov(spin_occ_num,spin_vir_num)
2019-09-11 17:09:30 +02:00
! Local variables
integer :: i,j,a,b
! Output variables
double precision,intent(out) :: EcCCSD
! Compute CCSD correlation energy
EcCCSD = 0d0
! Singles contribution
2019-09-14 14:15:25 +02:00
do i=1,spin_occ_num
do a=1,spin_vir_num
2019-09-11 17:09:30 +02:00
EcCCSD = EcCCSD + Fov(i,a)*t1(i,a)
end do
end do
! Doubles contribution
2019-09-14 14:15:25 +02:00
do i=1,spin_occ_num
do j=1,spin_occ_num
do a=1,spin_vir_num
do b=1,spin_vir_num
2019-09-11 17:09:30 +02:00
EcCCSD = EcCCSD &
+ 0.5d0*OOVV(i,j,a,b)*t1(i,a)*t1(j,b) &
+ 0.25d0*OOVV(i,j,a,b)*t2(i,j,a,b)
end do
end do
end do
end do
end subroutine CCSD_Ec_nc