calcFullStopR2L Subroutine

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

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)
logical, intent(in), optional :: t_no_singles_opt

Contents

Source Code


Source Code

    subroutine calcFullStopR2L(ilut, csf_i, excitInfo, excitations, nExcits, &
                               posSwitches, negSwitches, t_no_singles_opt)
        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)
        logical, intent(in), optional :: t_no_singles_opt

        ! not sure if single or double weight is necessary here...
        type(WeightObj_t) :: weights
        integer(n_int), allocatable :: tempExcits(:, :)
        integer(n_int) :: t(0:nifguga)
        integer :: iOrb, iEx, cnt, st, se, en, ierr
        real(dp) :: minusWeight, plusWeight, zeroWeight
        logical :: t_no_singles

        st = excitInfo%fullStart
        se = excitInfo%secondStart
        en = excitInfo%fullEnd

        if (present(t_no_singles_opt)) then
            t_no_singles = t_no_singles_opt
        else
            t_no_singles = .false.
        end if

        if (t_no_singles .and. csf_i%stepvector(en) == 3) then
            allocate(excitations(0, 0), stat=ierr)
            nExcits = 0
            return
        end if

        ! init weights
        weights = init_semiStartWeight(csf_i, se, en, negSwitches(se), &
                                       posSwitches(se), csf_i%B_real(se))

        excitInfo%currentGen = excitInfo%firstGen

        ! create start
        call createSingleStart(ilut, csf_i, excitInfo, posSwitches, negSwitches, &
                               weights, tempExcits, nExcits)

        ! loop until semi-start
        do iOrb = excitInfo%fullStart + 1, excitInfo%secondStart - 1
            call singleUpdate(ilut, csf_i, iOrb, excitInfo, posSwitches, negSwitches, &
                              weights, tempExcits, nExcits)
        end do

        ! can i write a general raisingSemiStart function, to reuse in
        ! other types of excitations?
        ! have to init specific prob. weights
        weights = weights%ptr
        ! depending if there is a 3 at the fullend there can be no switching
        ! possible
        if (csf_i%stepvector(en) == 3) then
            zeroWeight = weights%proc%zero(0.0_dp, 0.0_dp, csf_i%B_real(se), weights%dat)
            ! is plus and minus weight zero then? i think so
            minusWeight = 0.0_dp
            plusWeight = 0.0_dp

            if (t_mixed_hubbard) then
                nExcits = 0
                allocate(excitations(0, 0), stat=ierr)
                return
            end if
        else
            minusWeight = weights%proc%minus(negSwitches(se), csf_i%B_real(se), weights%dat)
            plusWeight = weights%proc%plus(posSwitches(se), csf_i%B_real(se), weights%dat)
            zeroWeight = weights%proc%zero(negSwitches(se), posSwitches(se), &
                                           csf_i%B_real(se), weights%dat)
        end if
        ! do i have to give them posSwitches and negSwitches or could I
        ! just put in actual weight values?
        call calcLoweringSemiStart(ilut, csf_i, excitInfo, tempExcits, nExcits, &
                                   plusWeight, minusWeight, zeroWeight)

        if (csf_i%stepvector(en) == 3) then
            ! in mixed generator case there is no sign coming from the
            ! overlap region, so only finish up the excitation, by just
            ! giving it the correct matrix element

            ! and have to check if non-zero and store in final
            ! excitations list

            cnt = 1
            do iEx = 1, nExcits
                t = tempExcits(:, iEx)

                if (near_zero(extract_matrix_element(t, 1))) cycle

                ! also no change in stepvector in this case
                call update_matrix_element(t, Root2, 1)
                call encode_matrix_element(t, 0.0_dp, 2)

                tempExcits(:, cnt) = t
                cnt = cnt + 1
            end do

            nExcits = cnt - 1

            allocate(excitations(0:nifguga, nExcits), stat=ierr)
            excitations = tempExcits(:, 1:nExcits)

            deallocate(tempExcits)

        else
            ! then i have to do a proper double excitation until the end
            ! with the new weights also
            do iOrb = excitInfo%secondStart + 1, excitInfo%fullEnd - 1
                call doubleUpdate(ilut, csf_i, iOrb, excitInfo, weights, tempExcits, nExcits, &
                                  negSwitches, posSwitches)
            end do

            ! and then do a mixed fullstop. also write a general function for that
            call mixedFullStop(ilut, csf_i, excitInfo, tempExcits, nExcits, excitations, &
                               t_no_singles)

        end if

    end subroutine calcFullStopR2L