calcAllExcitations_double Subroutine

public subroutine calcAllExcitations_double(ilut, csf_i, i, j, k, l, excitations, nExcits, t_full)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(in) :: ilut(0:nifguga)
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 :: excitations(:,:)
integer, intent(out) :: nExcits
logical, intent(in), optional :: t_full

Contents


Source Code

    subroutine calcAllExcitations_double(ilut, csf_i, i, j, k, l, excitations, nExcits,&
            t_full)
        ! function to calculate all possible double excitation for a CSF
        ! given in (ilut) format and indices (i,j,k,l).
        ! used to calculate the action of the Hamiltonian H|D> to calculate
        ! the projected energy
        integer(n_int), intent(in) :: ilut(0:nifguga)
        type(CSF_Info_t), intent(in) :: csf_i
        integer, intent(in) :: i, j, k, l
        integer(n_int), intent(out), allocatable :: excitations(:, :)
        integer, intent(out) :: nExcits
        logical, intent(in), optional :: t_full
        character(*), parameter :: this_routine = "calcAllExcitations_double"

        real(dp) :: posSwitches(nSpatOrbs), negSwitches(nSpatOrbs)
        HElement_t(dp) :: umat
        integer :: ierr, n, exlevel
        type(ExcitationInformation_t) :: excitInfo
        logical :: compFlag, t_full_
        def_default(t_full_, t_full, .true.)

        ! if called with k = l = 0 -> call single version of function
        if (k == 0 .and. l == 0) then
            call calcAllExcitations(ilut, csf_i, i, j, excitations, nExcits, t_full_)
            return
        end if

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

        ! default:
        nExcits = 0

        ! first check two-particle integral
        if (t_full_) then
            umat = get_umat_el(i, k, j, l)
        else
            umat = h_cast(1.0_dp)
        end if

        if (near_zero(umat)) then
            allocate(excitations(0, 0), stat=ierr)
            return
        end if

        ! 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)

        if (.not. compFlag) then
            allocate(excitations(0, 0), stat=ierr)
            return
        end if

        ! with the more involved excitation identification, i should write
        ! really specific double excitation creators

        ! maybe I should adapt this here also for other type of lattice
        ! models with restricted 2-body interaction..
        if (t_mixed_hubbard) then
            select case (excitInfo%typ)
            case (excit_type%single, &
                  excit_type%raising, &
                  excit_type%lowering, &
                  excit_type%single_overlap_lowering, &
                  excit_type%single_overlap_raising)

                allocate(excitations(0, 0), stat=ierr)
                return
            end select
        end if

        if (t_heisenberg_model) then
            if (excitInfo%typ /= excit_type%fullstart_stop_mixed) then
                allocate(excitations(0,0), stat = ierr)
                return
            end if
        end if

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

        case (excit_type%raising) ! weight + lowering 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, excitations, &
                                                 nExcits, posSwitches, negSwitches)

            exlevel = 1

        case (excit_type%lowering) ! weight + raising gen
            call calcDoubleExcitation_withWeight(ilut, csf_i, excitInfo, excitations, &
                                                 nExcits, posSwitches, negSwitches)

            exlevel = 1

        case (excit_type%non_overlap) ! non overlap
            call calcNonOverlapDouble(ilut, csf_i, excitInfo, excitations, nExcits, &
                                      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, excitations, nExcits, &
                                           posSwitches, negSwitches)

            exlevel = 1

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

            exlevel = 1

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

            exlevel = 2

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

            exlevel = 2

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

            exlevel = 2

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

            exlevel = 2

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

            exlevel = 2

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

            exlevel = 2

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

            exlevel = 2

        case (excit_type%double_R_to_L) ! raising into lowering double
            call calcDoubleR2L(ilut, csf_i, excitInfo, excitations, nExcits, &
                               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, excitations, nExcits, &
                                      posSwitches, negSwitches)

            exlevel = 2

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

            exlevel = 2

        case (excit_type%fullstop_L_to_R) ! full stop lowering into raising
            call calcFullStopL2R(ilut, csf_i, excitInfo, excitations, nExcits, &
                                 posSwitches, negSwitches)

            ! 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
            exlevel = 2

        case (excit_type%fullstop_R_to_L) ! full stop raising into lowering
            call calcFullStopR2L(ilut, csf_i, excitInfo, excitations, nExcits, &
                                 posSwitches, negSwitches)

            ! same as for 16
            exlevel = 2

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

            exlevel = 2

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

            exlevel = 2

        case (excit_type%fullStart_L_to_R) ! full start lowering into raising
            call calcFullStartL2R(ilut, csf_i, excitInfo, excitations, nExcits, &
                                  posSwitches, negSwitches)

            ! same as for 16
            exlevel = 2

        case (excit_type%fullstart_R_to_L) ! full start raising into lowering
            call calcFullStartR2L(ilut, csf_i, excitInfo, excitations, nExcits, &
                                  posSwitches, negSwitches)

            ! same as for 16
            exlevel = 2

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

            exlevel = 2

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

            ! same as for 16
            exlevel = 2

        end select

        ! for checking purposes encode umat/2 here to compare with sandeeps
        ! results for now..
        do n = 1, nExcits

            call encode_matrix_element(excitations(:, n), 0.0_dp, 2)
            if (t_full_) then
                call update_matrix_element(excitations(:, n), umat / 2.0_dp, 1)
            else
                call encode_matrix_element(excitations(:, n), 1.0_dp, 1)
            end if
            ! and also use the deltaB value for finished excitations to
            ! indicate the level of excitation IC for the remaining NECI code
            call setDeltaB(exlevel, excitations(:, n))

        end do

    end subroutine calcAllExcitations_double