@brief Form the ROHF Fock matrix in MO basis.
This implementation is based on the work of M. F. Guest and V. Saunders, as described in Mol. Phys. 28, 819 (1974).
@brief Form the ROHF Fock matrix in the MO basis using the Guest-Saunders method.
@detail This subroutine transforms the alpha and beta Fock matrices from the AO to the MO basis, constructs the ROHF Fock matrix using the Guest-Saunders method, and optionally applies a level shift to the diagonal elements of the virtual-virtual block.
@author Konstantin Komarov, 2023
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp), | intent(inout), | dimension(:) | :: | fock_a_ao | ||
real(kind=dp), | intent(inout), | dimension(:) | :: | fock_b_ao | ||
real(kind=dp), | intent(out), | dimension(:) | :: | fock_mo | ||
real(kind=dp), | intent(in), | dimension(:,:) | :: | MOs | ||
real(kind=dp), | intent(in), | dimension(:) | :: | overlap_tri | ||
real(kind=dp), | intent(out), | dimension(:) | :: | work | ||
integer, | intent(in) | :: | nocca | |||
integer, | intent(in) | :: | noccb | |||
integer, | intent(in) | :: | nbf | |||
real(kind=dp), | intent(in) | :: | vshift |
subroutine form_rohf_fock(fock_a_ao, fock_b_ao, fock_mo, & MOs, overlap_tri, work, & nocca, noccb, nbf, vshift) use mathlib, only: orthogonal_transform_sym, & orthogonal_transform2, & triangular_to_full, & unpack_matrix, & pack_matrix implicit none real(kind=dp), intent(inout), dimension(:) :: fock_a_ao real(kind=dp), intent(inout), dimension(:) :: fock_b_ao real(kind=dp), intent(out), dimension(:) :: fock_mo real(kind=dp), intent(in), dimension(:,:) :: MOs real(kind=dp), intent(in), dimension(:) :: overlap_tri real(kind=dp), intent(out), dimension(:) :: work integer, intent(in) :: nocca, noccb, nbf real(kind=dp), intent(in) :: vshift real(kind=dp), allocatable, dimension(:,:) :: & overlap, work_matrix, fock, fock_a, fock_b real(kind=dp) :: acc, aoo, avv, bcc, boo, bvv integer :: i, nbf_tri acc = 0.5_dp; aoo = 0.5_dp; avv = 0.5_dp bcc = 0.5_dp; boo = 0.5_dp; bvv = 0.5_dp nbf_tri = nbf * (nbf + 1) / 2 ! Allocate full matrices allocate(overlap(nbf, nbf), & work_matrix(nbf, nbf), & fock(nbf, nbf), & fock_a(nbf, nbf), & fock_b(nbf, nbf), & source=0.0_dp) ! Transform alpha and beta Fock matrices to MO basis call orthogonal_transform_sym(nbf, nbf, fock_a_ao, MOs, nbf, fock_mo) fock_a_ao(:nbf_tri) = fock_mo(:nbf_tri) call orthogonal_transform_sym(nbf, nbf, fock_b_ao, MOs, nbf, fock_mo) fock_b_ao(:nbf_tri) = fock_mo(:nbf_tri) ! Unpack triangular matrices to full matrices call unpack_matrix(fock_a_ao, fock_a) call unpack_matrix(fock_b_ao, fock_b) ! Construct ROHF Fock matrix in MO basis using Guest-Saunders method associate ( na => nocca & , nb => noccb & ) fock(1:nb, 1:nb) = acc * fock_a(1:nb, 1:nb) & + bcc * fock_b(1:nb, 1:nb) fock(nb+1:na, nb+1:na) = aoo * fock_a(nb+1:na, nb+1:na) & + boo * fock_b(nb+1:na, nb+1:na) fock(na+1:nbf, na+1:nbf) = avv * fock_a(na+1:nbf, na+1:nbf) & + bvv * fock_b(na+1:nbf, na+1:nbf) fock(1:nb, nb+1:na) = fock_b(1:nb, nb+1:na) fock(nb+1:na, 1:nb) = fock_b(nb+1:na, 1:nb) fock(1:nb, na+1:nbf) = 0.5_dp * (fock_a(1:nb, na+1:nbf) & + fock_b(1:nb, na+1:nbf)) fock(na+1:nbf, 1:nb) = 0.5_dp * (fock_a(na+1:nbf, 1:nb) & + fock_b(na+1:nbf, 1:nb)) fock(nb+1:na, na+1:nbf) = fock_a(nb+1:na, na+1:nbf) fock(na+1:nbf, nb+1:na) = fock_a(na+1:nbf, nb+1:na) ! Apply Vshift to the diagonal do i = nb+1, na fock(i,i) = fock(i,i) + vshift * 0.5_dp end do do i = na+1, nbf fock(i,i) = fock(i,i) + vshift end do end associate ! Pack fock to fock_mo, make it triangular call pack_matrix(fock, fock_mo) call triangular_to_full(fock_mo, nbf, 'u') ! Back-transform ROHF Fock matrix to AO basis call unpack_matrix(overlap_tri, overlap) call dsymm('l', 'u', nbf, nbf, & 1.0_dp, overlap, nbf, & MOs, nbf, & 0.0_dp, work, nbf) call orthogonal_transform2('t', nbf, nbf, work, nbf, fock, nbf, & work_matrix, nbf, overlap) call pack_matrix(work_matrix, fock_a_ao) end subroutine form_rohf_fock