1
0
mirror of https://gitlab.com/scemama/qp_plugins_scemama.git synced 2024-11-07 14:43:41 +01:00
qp_plugins_scemama/devel/cc/curly_Fock.irp.f

102 lines
2.7 KiB
Fortran

BEGIN_PROVIDER [ double precision, c_spin_fock_matrix_mo_oo, (spin_occ_num,spin_occ_num) ]
implicit none
BEGIN_DOC
! Curly F in Occupied-Occupied block
END_DOC
double precision,external :: Kronecker_Delta
integer :: i,j,m,n
integer :: a,b,e,f
do m=1,spin_occ_num
do i=1,spin_occ_num
c_spin_fock_matrix_mo_oo(m,i) = (1d0 - Kronecker_delta(m,i))*spin_fock_matrix_mo_oo(m,i)
do e=1,spin_vir_num
c_spin_fock_matrix_mo_oo(m,i) = c_spin_fock_matrix_mo_oo(m,i) + 0.5d0*t1_cc(i,e)*spin_fock_matrix_mo_ov(m,e)
end do
do e=1,spin_vir_num
do n=1,spin_occ_num
c_spin_fock_matrix_mo_oo(m,i) = c_spin_fock_matrix_mo_oo(m,i) + t1_cc(n,e)*OOOV(m,n,i,e)
end do
end do
do e=1,spin_vir_num
do n=1,spin_occ_num
do f=1,spin_vir_num
c_spin_fock_matrix_mo_oo(m,i) = c_spin_fock_matrix_mo_oo(m,i) + 0.5d0*taus(i,n,e,f)*OOVV(m,n,e,f)
end do
end do
end do
end do
end do
END_PROVIDER
BEGIN_PROVIDER [ double precision, c_spin_fock_matrix_mo_ov, (spin_occ_num,spin_vir_num) ]
implicit none
BEGIN_DOC
! Curly F in Occupied-Virtual block
END_DOC
integer :: i,j,m,n
integer :: a,b,e,f
c_spin_fock_matrix_mo_ov(:,:) = spin_fock_matrix_mo_ov(:,:)
do m=1,spin_occ_num
do e=1,spin_vir_num
do n=1,spin_occ_num
do f=1,spin_vir_num
c_spin_fock_matrix_mo_ov(m,e) = c_spin_fock_matrix_mo_ov(m,e) + t1_cc(n,f)*OOVV(m,n,e,f)
end do
end do
end do
end do
END_PROVIDER
BEGIN_PROVIDER [ double precision, c_spin_fock_matrix_mo_vv, (spin_vir_num,spin_vir_num) ]
implicit none
BEGIN_DOC
! Curly F in Occupied-Virtual block
END_DOC
double precision,external :: Kronecker_Delta
integer :: i,j,m,n
integer :: a,b,e,f
do a=1,spin_vir_num
do e=1,spin_vir_num
c_spin_fock_matrix_mo_vv(a,e) = (1d0 - Kronecker_delta(a,e))*spin_fock_matrix_mo_vv(a,e)
do m=1,spin_occ_num
c_spin_fock_matrix_mo_vv(a,e) = c_spin_fock_matrix_mo_vv(a,e) - 0.5d0*t1_cc(m,a)*spin_fock_matrix_mo_ov(m,e)
end do
do m=1,spin_occ_num
do f=1,spin_vir_num
c_spin_fock_matrix_mo_vv(a,e) = c_spin_fock_matrix_mo_vv(a,e) + t1_cc(m,f)*OVVV(m,a,f,e)
end do
end do
do m=1,spin_occ_num
do n=1,spin_occ_num
do f=1,spin_vir_num
c_spin_fock_matrix_mo_vv(a,e) = c_spin_fock_matrix_mo_vv(a,e) - 0.5d0*taus(m,n,a,f)*OOVV(m,n,e,f)
end do
end do
end do
end do
end do
END_PROVIDER