calcSingleOverlapMixed Subroutine

public subroutine calcSingleOverlapMixed(ilut, csf_i, excitInfo, excitations, nExcits, posSwitches, negSwitches)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(in) :: ilut(0:nifguga)
type(CSF_Info_t), intent(in) :: csf_i
type(ExcitationInformation_t), intent(inout) :: excitInfo
integer(kind=n_int), intent(out), allocatable :: excitations(:,:)
integer, intent(out) :: nExcits
real(kind=dp), intent(in) :: posSwitches(nSpatOrbs)
real(kind=dp), intent(in) :: negSwitches(nSpatOrbs)

Contents


Source Code

    subroutine calcSingleOverlapMixed(ilut, csf_i, excitInfo, excitations, nExcits, &
                                      posSwitches, negSwitches)
        ! control routine to calculate the single overlap excitation with
        ! mixed generators
        integer(n_int), intent(in) :: ilut(0:nifguga)
        type(CSF_Info_t), intent(in) :: csf_i
        type(ExcitationInformation_t), intent(inout) :: excitInfo
        integer(n_int), intent(out), allocatable :: excitations(:, :)
        integer, intent(out) :: nExcits
        real(dp), intent(in) :: posSwitches(nSpatOrbs), negSwitches(nSpatOrbs)
        character(*), parameter :: this_routine = "calcSingleOverlapMixed"

        type(WeightObj_t) :: weights
        integer :: iOrb, iEx, deltaB
        integer(n_int) :: t(0:nifguga)
        integer(n_int), allocatable :: tempExcits(:, :)
        real(dp) :: tempWeight

        ! todo: create weight objects, here are the normal single excitation
        ! weight operators.
        weights = init_singleWeight(csf_i, excitInfo%fullEnd)

        excitInfo%currentGen = excitInfo%firstGen
        ! create according start
        call createSingleStart(ilut, csf_i, excitInfo, posSwitches, negSwitches, &
                               weights, tempExcits, nExcits)

        ! set current gen
        ! do normal updates up to the overlap site
        do iOrb = excitInfo%fullStart + 1, excitInfo%secondStart - 1
            call singleUpdate(ilut, csf_i, iOrb, excitInfo, posSwitches, negSwitches, &
                              weights, tempExcits, nExcits)
        end do

        ! at single overlap site depending on which generator ends or start
        iOrb = excitInfo%secondStart

        if (excitInfo%firstGen == gen_type%L) then
            ! lowering ends
            ASSERT(isZero(ilut, iOrb))

            do iEx = 1, nExcits
                t = tempExcits(:, iEx)
                deltaB = getDeltaB(t)

                ! have to get double excitation matrix elements in here..
                call getDoubleMatrixElement(3, 0, deltaB, excitInfo%firstGen, &
                                            excitInfo%lastGen, csf_i%B_real(iOrb), 1.0_dp, tempWeight)

                call update_matrix_element(t, tempWeight, 1)

                ! change 0 -> 3
                set_orb(t, 2 * iOrb)
                set_orb(t, 2 * iOrb - 1)

                ! store
                tempExcits(:, iEx) = t

            end do

        else
            ! raising end
            ASSERT(isThree(ilut, iOrb))

            do iEx = 1, nExcits
                t = tempExcits(:, iEx)
                deltaB = getDeltaB(t)

                call getDoubleMatrixElement(0, 3, deltaB, excitInfo%firstGen, &
                                            excitInfo%lastGen, csf_i%B_real(iOrb), 1.0_dp, tempWeight)

                call update_matrix_element(t, tempWeight, 1)

                ! change 3 -> 0
                clr_orb(t, 2 * iOrb)
                clr_orb(t, 2 * iOrb - 1)

                ! store
                tempExcits(:, iEx) = t
            end do
        end if

        ! change excitInfo and continue singleUpdates
        excitInfo%currentGen = excitInfo%lastGen

        do iOrb = excitInfo%secondStart + 1, excitInfo%fullEnd - 1
            call singleUpdate(ilut, csf_i, iOrb, excitInfo, posSwitches, negSwitches, &
                              weights, tempExcits, nExcits)
        end do

        ! lets see if that calling works
        call singleEnd(ilut, csf_i, excitInfo, tempExcits, nExcits, excitations)

    end subroutine calcSingleOverlapMixed