calc_all_excits_guga_rdm_doubles Subroutine

public subroutine calc_all_excits_guga_rdm_doubles(ilut, csf_i, i, j, k, l, excits, n_excits)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(in) :: ilut(0:GugaBits%len_tot)
type(CSF_Info_t), intent(in) :: csf_i
integer, intent(in) :: i
integer, intent(in) :: j
integer, intent(in) :: k
integer, intent(in) :: l
integer(kind=n_int), intent(out), allocatable :: excits(:,:)
integer, intent(out) :: n_excits

Contents


Source Code

    subroutine calc_all_excits_guga_rdm_doubles(ilut, csf_i, i, j, k, l, excits, n_excits)
        integer(n_int), intent(in) :: ilut(0:GugaBits%len_tot)
        type(CSF_Info_t), intent(in) :: csf_i
        integer, intent(in) :: i, j, k, l
        integer(n_int), intent(out), allocatable :: excits(:, :)
        integer, intent(out) :: n_excits
        character(*), parameter :: this_routine = "calc_all_excits_guga_rdm_doubles"

        real(dp) :: posSwitches(nSpatOrbs), negSwitches(nSpatOrbs)
        integer :: ierr, n, exlevel
        type(ExcitationInformation_t) :: excitInfo
        logical :: compFlag

        ASSERT(isProperCSF_ilut(ilut))
        ASSERT(i > 0 .and. i <= nSpatOrbs)
        ASSERT(j > 0 .and. j <= nSpatOrbs)

        ! if called with k = l = 0 -> call single version of function
        if (k == 0 .and. l == 0) then
            call calc_all_excits_guga_rdm_singles(ilut, csf_i, i, j, excits, n_excits)
            return
        end if

        ASSERT(k > 0 .and. k <= nSpatOrbs)
        ASSERT(l > 0 .and. l <= nSpatOrbs)

        ! default:
        n_excits = 0

        ! otherwise get excitation information
        excitInfo = excitationIdentifier(i, j, k, l)

        ! screw it. for now write a function which checks if indices and ilut
        ! are compatible, and not initiate the excitation right away
        ! but check cases again
        ! in checkCompatibility the number of switches is already
        ! calulated to check if the probabilistic weights fit... maybe but
        ! that out and reuse.. to not waste any effort.
        call checkCompatibility(csf_i, excitInfo, compFlag, posSwitches, negSwitches)

        ! for mixed type full starts and/or full stops i have to consider
        ! the possible diagonal/single excitations here!
        ! although I think the full-stops are already correctly taken into
        ! account..
        ! and also the the full-starts are maybe correct already..
        ! so it was just the full-start into full-stop mixed!
        if (.not. compFlag) then
            allocate(excits(0, 0), stat=ierr)
            return
        end if

        select case (excitInfo%typ)

        case (excit_type%raising, &
              excit_type%lowering, &
              excit_type%single_overlap_lowering, &
              excit_type%single_overlap_raising)

            ! in the case of mimicking stochasitic
            ! excitation generation, we should abort here for
            ! these type of excitations!
            allocate(excits(0, 0), stat=ierr)
            n_excits = 0
            return

        end select

        select case (excitInfo%typ)
        case (excit_type%single)
            ! shouldnt be here.. onyl single excits and full weight gens
            allocate(excits(0, 0), stat=ierr)
            return

        case (excit_type%raising) ! weight + raising gen.
            ! can be treated almost like a single excitation
            ! essentially the same, except if d(w) == 3 in the excitaton regime

            call calcDoubleExcitation_withWeight(ilut, csf_i, excitInfo, excits, &
                                                 n_excits, posSwitches, negSwitches)

            exlevel = 1

        case (excit_type%lowering) ! weight + lowering gen
            call calcDoubleExcitation_withWeight(ilut, csf_i, excitInfo, excits, &
                                                 n_excits, posSwitches, negSwitches)

            exlevel = 1

        case (excit_type%non_overlap) ! non overlap
            call calcNonOverlapDouble(ilut, csf_i, excitInfo, excits, n_excits, &
                                      posSwitches, negSwitches)

            exlevel = 2

        case (excit_type%single_overlap_lowering) ! single overlap two lowering
            ! how can i efficiently adress that?
            ! can i write that efficiently in one function or do i need more?
            ! probably need more... i already determined
            call calcSingleOverlapLowering(ilut, csf_i, excitInfo, excits, n_excits, &
                                           posSwitches, negSwitches)

            exlevel = 1

        case (excit_type%single_overlap_raising) ! single overlap raising
            call calcSingleOverlapRaising(ilut, csf_i, excitInfo, excits, n_excits, &
                                          posSwitches, negSwitches)

            exlevel = 1

        case (excit_type%single_overlap_L_to_R) ! single overlap lowering into raising
            call calcSingleOverlapMixed(ilut, csf_i, excitInfo, excits, n_excits, &
                                        posSwitches, negSwitches)

            exlevel = 2

        case (excit_type%single_overlap_R_to_L) ! single overlap raising into lowering
            call calcSingleOverlapMixed(ilut, csf_i, excitInfo, excits, n_excits, &
                                        posSwitches, negSwitches)

            exlevel = 2

        case (excit_type%double_lowering) ! normal double overlap two lowering
            call calcDoubleLowering(ilut, csf_i, excitInfo, excits, n_excits, &
                                    posSwitches, negSwitches)

            exlevel = 2

        case (excit_type%double_raising) ! normal double overlap two raising
            call calcDoubleRaising(ilut, csf_i, excitInfo, excits, n_excits, &
                                   posSwitches, negSwitches)

            exlevel = 2

        case (excit_type%double_L_to_R_to_L) ! lowering into raising into lowering
            call calcDoubleRaising(ilut, csf_i, excitInfo, excits, n_excits, &
                                   posSwitches, negSwitches)

            exlevel = 2

        case (excit_type%double_R_to_L_to_R) ! raising into lowering into raising
            call calcDoubleLowering(ilut, csf_i, excitInfo, excits, n_excits, &
                                    posSwitches, negSwitches)

            exlevel = 2

        case (excit_type%double_L_to_R) ! lowering into raising double
            call calcDoubleL2R(ilut, csf_i, excitInfo, excits, n_excits, &
                               posSwitches, negSwitches)

            exlevel = 2

        case (excit_type%double_R_to_L) ! raising into lowering double
            call calcDoubleR2L(ilut, csf_i, excitInfo, excits, n_excits, &
                               posSwitches, negSwitches)

            exlevel = 2

        case (excit_type%fullstop_lowering) ! full stop 2 lowering
            ! can i write a function for both alike generator combinations
            ! i think i can
            call calcFullstopLowering(ilut, csf_i, excitInfo, excits, n_excits, &
                                      posSwitches, negSwitches)

            exlevel = 2

        case (excit_type%fullstop_raising) ! full stop 2 raising
            call calcFullstopRaising(ilut, csf_i, excitInfo, excits, n_excits, &
                                     posSwitches, negSwitches)

            exlevel = 2

        case (excit_type%fullstop_L_to_R) ! full stop lowering into raising
            call calcFullStopL2R(ilut, csf_i, excitInfo, excits, n_excits, &
                                 posSwitches, negSwitches, t_no_singles_opt=.true.)

            ! in this case there is also the possibility for one single-like
            ! excitation if there is no change in the double overlap region!
            ! todo! how to fix that? is that so important? its only max. 1

            ! depending on the full-stop step-value:
            ! if it is 3, all excitations are single like, and should be
            ! disregarded if we mimic stochastic sampling.

            ! if the end step value is 1 or 2, there is on Delta B = 0
            ! branch with ofc. multiple possible singles
            ! associated with it.
            ! i think i should use a optional flag to indicate that I want
            ! to mimick stochastic exctiations generation
            exlevel = 2

        case (excit_type%fullstop_R_to_L) ! full stop raising into lowering
            call calcFullStopR2L(ilut, csf_i, excitInfo, excits, n_excits, &
                                 posSwitches, negSwitches, t_no_singles_opt=.true.)

            ! same as for 16
            exlevel = 2

        case (excit_type%fullstart_lowering) ! full start 2 lowering
            call calcFullStartLowering(ilut, csf_i, excitInfo, excits, n_excits, &
                                       posSwitches, negSwitches)

            exlevel = 2

        case (excit_type%fullstart_raising) ! full start 2 raising
            call calcFulLStartRaising(ilut, csf_i, excitInfo, excits, n_excits, &
                                      posSwitches, negSwitches)

            exlevel = 2

        case (excit_type%fullstart_L_to_R) ! full start lowering into raising
            call calcFullStartL2R(ilut, csf_i, excitInfo, excits, n_excits, &
                                  posSwitches, negSwitches, t_no_singles_opt=.true.)

            ! same as for 16
            exlevel = 2

        case (excit_type%fullstart_R_to_L) ! full start raising into lowering
            call calcFullStartR2L(ilut, csf_i, excitInfo, excits, n_excits, &
                                  posSwitches, negSwitches, t_no_singles_opt=.true.)

            ! same as for 16
            exlevel = 2

        case (excit_type%fullstart_stop_alike) ! full start into full stop alike
            call calcFullStartFullStopAlike(ilut, csf_i, excitInfo, excits)
            n_excits = 1

            exlevel = 2

        case (excit_type%fullstart_stop_mixed) ! full start into full stop mixed
            call calcFullStartFullStopMixed(ilut, csf_i, excitInfo, excits, n_excits, &
                                            posSwitches, negSwitches)

            ! same as for 16
            exlevel = 2

        end select

        ! indicate the level of excitation IC for the remaining NECI code
        if (n_excits > 0) then
            do n = 1, n_excits
                ! for consistency also encode x0 and x1 in the ind_rdm_x0
                ! and ind_rdm_x1 entries
                ! damn.. i already added x0 and x1 in the routines above
                ! i think...
                call encode_rdm_ind(excits(:, n), contract_2_rdm_ind(i, j, k, l))
            end do
        end if

    end subroutine calc_all_excits_guga_rdm_doubles