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/r2_cc.irp.f
2019-09-18 00:52:53 +02:00

171 lines
4.6 KiB
Fortran

BEGIN_PROVIDER [ double precision, r2_cc, (spin_occ_num,spin_occ_num,spin_vir_num,spin_vir_num) ]
implicit none
BEGIN_DOC
! t2 residues in non-canonical CCSD
END_DOC
integer :: i,j,m,n
integer :: a,b,e,f
double precision :: x
double precision, allocatable :: tmp_oo(:,:), tmp_oo2(:,:)
double precision, external :: u_dot_v
allocate ( tmp_oo(spin_occ_num,spin_occ_num) )
r2_cc(:,:,:,:) = OOVV(:,:,:,:)
do b=1,spin_vir_num
do a=1,spin_vir_num
if (a == b) cycle
do e=1,spin_vir_num
x = u_dot_v( t1_cc(1,b), c_spin_fock_matrix_mo_ov(1,e), spin_occ_num )
x = c_spin_fock_matrix_mo_vv(b,e) - x*0.5d0
if (x == 0.d0) cycle
tmp_oo(:,:) = t2_cc(:,:,a,e) * x
r2_cc(:,:,a,b) = r2_cc(:,:,a,b) + tmp_oo(:,:)
r2_cc(:,:,b,a) = r2_cc(:,:,b,a) - tmp_oo(:,:)
end do
do i=1,spin_occ_num
do m=1,spin_occ_num
tmp_oo(m,i) = t2_cc(i,m,a,b)
enddo
enddo
do j=1,spin_occ_num
do i=1,spin_occ_num
x = 0.d0
do m=1,spin_occ_num
do e=1,spin_vir_num
x = x + c_spin_fock_matrix_mo_ov(m,e) * &
( tmp_oo(m,j)*t1_cc(i,e) - &
tmp_oo(m,i)*t1_cc(j,e) )
end do
end do
x = x + u_dot_v( tau_cc(1,1,a,b), cWoooo(1,1,i,j), spin_occ_num*spin_occ_num )
x = x * 0.5d0 + &
u_dot_v( tmp_oo(1,j), c_spin_fock_matrix_mo_oo(1,i), spin_occ_num ) - &
u_dot_v( tmp_oo(1,i), c_spin_fock_matrix_mo_oo(1,j), spin_occ_num )
do e=1,spin_vir_num
x = x - sum(t2_cc(j,:,a,e)*cWovvo(:,b,e,i))
end do
r2_cc(i,j,a,b) = r2_cc(i,j,a,b) + x
end do
end do
tmp_oo(:,:) = 0.d0
do f=1,spin_vir_num
do e=1,spin_vir_num
x = 0.5d0*cWvvvv(a,b,e,f)
if (x == 0.d0) cycle
tmp_oo(:,:) = tmp_oo(:,:) + x*tau_cc(:,:,e,f)
end do
end do
r2_cc(:,:,a,b) = r2_cc(:,:,a,b) + tmp_oo(:,:)
do j=1,spin_occ_num
do e=1,spin_vir_num
do m=1,spin_occ_num
x = t1_cc(m,a)*OVVO(m,b,e,j)
do i=1,spin_occ_num
r2_cc(i,j,a,b) = r2_cc(i,j,a,b) + t2_cc(i,m,a,e)*cWovvo(m,b,e,j) - t1_cc(i,e)*x
enddo
end do
end do
end do
do e=1,spin_vir_num
do j=1,spin_occ_num
do i=1,spin_occ_num
do m=1,spin_occ_num
r2_cc(i,j,a,b) = r2_cc(i,j,a,b) + t1_cc(j,e)*t1_cc(m,a)*OVVO(m,b,e,i)
end do
end do
end do
end do
do e=1,spin_vir_num
do j=1,spin_occ_num
do m=1,spin_occ_num
do i=1,spin_occ_num
r2_cc(i,j,a,b) = r2_cc(i,j,a,b) - t2_cc(i,m,b,e)*cWovvo(m,a,e,j)
end do
end do
end do
end do
do j=1,spin_occ_num
do e=1,spin_vir_num
do i=1,spin_occ_num
do m=1,spin_occ_num
r2_cc(i,j,a,b) = r2_cc(i,j,a,b) + t1_cc(i,e)*t1_cc(m,b)*OVVO(m,a,e,j)
end do
end do
end do
end do
do e=1,spin_vir_num
do j=1,spin_occ_num
do i=1,spin_occ_num
do m=1,spin_occ_num
r2_cc(i,j,a,b) = r2_cc(i,j,a,b) + t2_cc(j,m,b,e)*cWovvo(m,a,e,i)
end do
end do
end do
end do
do e=1,spin_vir_num
do j=1,spin_occ_num
do i=1,spin_occ_num
do m=1,spin_occ_num
r2_cc(i,j,a,b) = r2_cc(i,j,a,b) - t1_cc(j,e)*t1_cc(m,b)*OVVO(m,a,e,i)
end do
end do
end do
end do
do j=1,spin_occ_num
do e=1,spin_vir_num
do i=1,spin_occ_num
r2_cc(i,j,a,b) = r2_cc(i,j,a,b) + t1_cc(i,e)*VVVO(a,b,e,j)
end do
end do
end do
do e=1,spin_vir_num
do i=1,spin_occ_num
do j=1,spin_occ_num
r2_cc(i,j,a,b) = r2_cc(i,j,a,b) - t1_cc(j,e)*VVVO(a,b,e,i)
end do
end do
end do
do j=1,spin_occ_num
do i=1,spin_occ_num
do m=1,spin_occ_num
r2_cc(i,j,a,b) = r2_cc(i,j,a,b) - t1_cc(m,a)*OVOO(m,b,i,j)
end do
end do
end do
do j=1,spin_occ_num
do i=1,spin_occ_num
do m=1,spin_occ_num
r2_cc(i,j,a,b) = r2_cc(i,j,a,b) + t1_cc(m,b)*OVOO(m,a,i,j)
end do
end do
end do
end do
end do
! Final expression of the t2 residue
r2_cc(:,:,:,:) = delta_oovv(:,:,:,:)*t2_cc(:,:,:,:) - r2_cc(:,:,:,:)
END_PROVIDER