calcAllExcitations_excitInfo_single Subroutine

private subroutine calcAllExcitations_excitInfo_single(ilut, csf_i, excitInfo, posSwitches, negSwitches, tmatFlag, excitations, nExcits)

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(in) :: excitInfo
real(kind=dp), intent(in) :: posSwitches(nSpatOrbs)
real(kind=dp), intent(in) :: negSwitches(nSpatOrbs)
logical, intent(in) :: tmatFlag
integer(kind=n_int), intent(out), allocatable :: excitations(:,:)
integer, intent(out) :: nExcits

Contents


Source Code

    subroutine calcAllExcitations_excitInfo_single(ilut, csf_i, excitInfo, posSwitches, &
               negSwitches, tmatFlag, excitations, nExcits)
        ! excitation calculation if excitInfo is already calculated
        integer(n_int), intent(in) :: ilut(0:nifguga)
        type(CSF_Info_t), intent(in) :: csf_i
        type(ExcitationInformation_t), intent(in) :: excitInfo
        logical, intent(in) :: tmatFlag
        real(dp), intent(in) :: posSwitches(nSpatOrbs), negSwitches(nSpatOrbs)
        integer(n_int), intent(out), allocatable :: excitations(:, :)
        integer, intent(out) :: nExcits

        HElement_t(dp) :: tmat
        integer :: iOrb
        type(WeightObj_t) :: weights
        integer(n_int), allocatable :: tempExcits(:, :)

        ! to do combine (i,j) version and this version to one

        ! if tmatFlag is false do not check or include the tmat element
        ! for double excitations which essentially only involve sinlge excits

        ! i think this is always called with the flag as false... hm,
        ! maybe redundant and remove.
        ! when is this called actually? because this influences if i can use
        ! the csf_i%stepvector info and stuff..
        ! i am pretty sure that in all the calls to this routine
        ! the necessary quantities are known!
        if (tmatFlag) then
            tmat = getTmatEl(2 * excitInfo%i, 2 * excitInfo%j)
            if (near_zero(tmat)) then
                nExcits = 0
                allocate(excitations(0, 0))
                return
            end if
        else
            tmat = 1.0_dp
        end if

        ! also do not need to check if excitation is compatible, since this
        ! has already been done
        weights = init_singleWeight(csf_i, excitInfo%fullEnd)

        ! have to give probabilistic weight object as input, to deal
        call createSingleStart(ilut, csf_i, excitInfo, posSwitches, &
                               negSwitches, weights, tempExcits, nExcits)

        ! to not call getTmatEl again in createSingleStart loop over
        ! the atmost two excitations here and multiply with tmat
        do iOrb = excitInfo%fullStart + 1, excitInfo%fullEnd - 1
            call singleUpdate(ilut, csf_i, iOrb, excitInfo, posSwitches, &
                              negSwitches, weights, tempExcits, nExcits)
        end do

        call singleEnd(ilut, csf_i, excitInfo, tempExcits, &
                       nExcits, excitations)

    end subroutine calcAllExcitations_excitInfo_single