ov_exact Subroutine

public subroutine ov_exact(temp1, i1, i2, ia1, ia2, s_mo, ilow, noc, itype)

Uses

  • proc~~ov_exact~~UsesGraph proc~ov_exact ov_exact module~precision precision proc~ov_exact->module~precision iso_fortran_env iso_fortran_env module~precision->iso_fortran_env

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(out) :: temp1
integer, intent(in) :: i1
integer, intent(in) :: i2
integer, intent(in) :: ia1
integer, intent(in) :: ia2
real(kind=dp), intent(in), dimension(:,:) :: s_mo
integer, intent(in) :: ilow
integer, intent(in) :: noc
integer, intent(in) :: itype

Calls

proc~~ov_exact~~CallsGraph proc~ov_exact ov_exact proc~comp_det comp_det proc~ov_exact->proc~comp_det

Called by

proc~~ov_exact~~CalledByGraph proc~ov_exact ov_exact proc~mrsf_tlf mrsf_tlf proc~mrsf_tlf->proc~ov_exact proc~compute_states_overlap compute_states_overlap proc~compute_states_overlap->proc~mrsf_tlf proc~get_states_overlap get_states_overlap proc~get_states_overlap->proc~compute_states_overlap proc~get_state_overlap_c get_state_overlap_C proc~get_state_overlap_c->proc~get_states_overlap

Source Code

  subroutine ov_exact(temp1, i1, i2, ia1, ia2, s_mo, ilow, noc, itype)

    use precision, only: dp

    implicit none

    real(kind=dp), intent(out) :: temp1
    integer, intent(in) :: i1, i2, ia1, ia2
    real(kind=dp), intent(in), dimension(:,:) :: s_mo
    integer, intent(in) :: ilow, noc, itype

    real(kind=dp), dimension(noc*noc) :: ddet
    integer :: i, iipp, imax, imin, ipp

    select case (itype)
    case (1)
       imin = min(i1,i2)
       imax = max(i1,i2)
    !  (1,1) block
       do i = 1, imin-1
          do ipp = 1, imin-1
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+ilow-1,ipp+ilow-1)
          end do
       end do
    !  (1,2) block
       do i = 1, imin-1
          do ipp = imin, imax-2
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+ilow-1,ipp+1+ilow-1)
          end do
       end do
    !  (1,3) block
       do i = 1, imin-1
          do ipp = imax-1, noc-1
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+ilow-1,ipp+2+ilow-1)
          end do
       end do
    !  (2,1) block
       do i = imin, imax-2
          do ipp = 1, imin-1
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+1+ilow-1,ipp+ilow-1)
          end do
       end do
    !  (2,2) block
       do i = imin, imax-2
          do ipp = imin, imax-2
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+1+ilow-1,ipp+1+ilow-1)
          end do
       end do
    !  (2,3) block
       do i = imin, imax-2
          do ipp = imax-1, noc-1
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+1+ilow-1,ipp+2+ilow-1)
          end do
       end do
    !  (3,1) block
       do i = imax-1, noc-1
          do ipp = 1, imin-1
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+2+ilow-1,ipp+ilow-1)
          end do
       end do
    !  (3,2) block
       do i = imax-1, noc-1
          do ipp = imin, imax-2
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+2+ilow-1,ipp+1+ilow-1)
          end do
       end do
    !  (3,3) block
       do i = imax-1, noc-1
          do ipp = imax-1, noc-1
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+2+ilow-1,ipp+2+ilow-1)
          end do
       end do
    !  (1,4) block
       do i = 1, imin-1
          do ipp = noc, noc
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+ilow-1,i1+ilow-1)
          end do
       end do
    !  (2,4) block
       do i = imin, imax-2
          do ipp = noc, noc
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+1+ilow-1,i1+ilow-1)
          end do
       end do
    !  (3,4) block
       do i = imax-1, noc-1
          do ipp = noc, noc
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+2+ilow-1,i1+ilow-1)
          end do
       end do
    !  (4,1) block
       do i = noc, noc
          do ipp = 1, imin-1
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i2+ilow-1,ipp+ilow-1)
          end do
       end do
    !  (4,2) block
       do i = noc, noc
          do ipp = imin, imax-2
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i2+ilow-1,ipp+1+ilow-1)
          end do
       end do
    !  (4,3) block
       do i = noc, noc
          do ipp = imax-1, noc-1
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i2+ilow-1,ipp+2+ilow-1)
          end do
       end do
    !  (4,4) block
       do i = noc, noc
          do ipp = noc, noc
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i2+ilow-1,i1+ilow-1)
          end do
       end do
    !  Calculate alpha determinant
       temp1 = comp_det(ddet, noc)
       if (i1==i2) then
          return
       else if (i1/=i2) then
          temp1 = -1.0_dp*temp1
          return
       endif

    case (2)
    !  (1,1) block
       do i = 1, noc-1
          do ipp = 1, noc-1
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+ilow-1,ipp+ilow-1)
          end do
       end do
    !  (1,2) block
       ipp = noc
       do i = 1, noc-1
          iipp = (ipp-1)*noc+i
          ddet(iipp) = s_mo(i+ilow-1,ia2)
       end do
    !  (2,1) block
       i = noc
       do ipp = 1, noc-1
          iipp = (ipp-1)*noc+i
          ddet(iipp) = s_mo(ia1,ipp+ilow-1)
       end do
    !  (2,2) block
       i = noc
       ipp = noc
       iipp = (ipp-1)*noc+i
       ddet(iipp) = s_mo(ia1,ia2)

    !  Calculate 2 det
       temp1 = comp_det(ddet, noc)
       return

    case (3)
    !  (1,1) block
       do i = 1, i1-1
          do ipp = 1, i1-1
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i,ipp)
          end do
       end do
    !  (1,2) block
       do i = 1, i1-1
          do ipp = i1, noc-2
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i,ipp+1)
          end do
       end do
    !  (1,3) block
       do i = 1, i1-1
          ipp = noc-1
          iipp = (ipp-1)*noc+i
          ddet(iipp) = s_mo(i,i1)
       end do
    !  (1,4) block
       do i = 1, i1-1
          ipp = noc
          iipp = (ipp-1)*noc+i
          ddet(iipp) = s_mo(i,ia1)
       end do
    !  (2,1) block
       do i = i1, noc-2
          do ipp = 1, i1-1
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+1,ipp)
          end do
       end do
    !  (2,2) block
       do i = i1, noc-2
          do ipp = i1, noc-2
             iipp = (ipp-1)*noc+i
             ddet(iipp) = s_mo(i+1,ipp+1)
          end do
       end do
    !  (2,3) block
       do i = i1, noc-2
          ipp = noc-1
          iipp = (ipp-1)*noc+i
          ddet(iipp) = s_mo(i+1,i1)
       end do
    !  (2,4) block
       do i = i1, noc-2
          ipp = noc
          iipp = (ipp-1)*noc+i
          ddet(iipp) = s_mo(i+1,ia1)
       end do
    !  (3,1) block
       i = noc-1
       do ipp = 1, i1-1
          iipp = (ipp-1)*noc+i
          ddet(iipp) = s_mo(i+1,ipp)
       end do
    !  (3,2) block
       i = noc-1
       do ipp = i1, noc-2
          iipp = (ipp-1)*noc+i
          ddet(iipp) = s_mo(i+1,ipp+1)
       end do
    !  (3,3) block
       i = noc-1
       ipp = noc-1
       iipp = (ipp-1)*noc+i
       ddet(iipp) = s_mo(i+1,i1)
    !  (3,4) block
       i = noc-1
       ipp = noc
       iipp = (ipp-1)*noc+i
       ddet(iipp) = s_mo(i+1,ia1)
    !  (4,1) block
       i = noc
       do ipp = 1, i1-1
          iipp = (ipp-1)*noc+i
          ddet(iipp) = s_mo(i+1,ipp)
       end do
    !  (4,2) block
       i = noc
       do ipp = i1, noc-2
          iipp = (ipp-1)*noc+i
          ddet(iipp) = s_mo(i+1,ipp+1)
       end do
    !  (4,3) block
       i = noc
       ipp = noc-1
       iipp = (ipp-1)*noc+i
       ddet(iipp) = s_mo(i+1,i1)
    !  (4,4) block
       i = noc
       ipp = noc
       iipp = (ipp-1)*noc+i
       ddet(iipp) = s_mo(i+1,ia1)

    !  Calculate alpha determinant
       temp1 = comp_det(ddet, noc)
    end select

  end subroutine ov_exact