forced_mixed_start Subroutine

private subroutine forced_mixed_start(ilut, csf_i, excitInfo, t, probWeight)

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
integer(kind=n_int), intent(out) :: t(0:nifguga)
real(kind=dp), intent(out) :: probWeight

Contents

Source Code


Source Code

    subroutine forced_mixed_start(ilut, csf_i, excitInfo, t, probWeight)
        ! NOTE: mixed full-start matrix elements are stores in the same row
        ! as delta B = -1 ones -> so access getDoubleMatrixElement with
        ! db = -1 below!
        integer(n_int), intent(in) :: ilut(0:nifguga)
        type(CSF_Info_t), intent(in) :: csf_i
        type(ExcitationInformation_t), intent(in) :: excitInfo
        integer(n_int), intent(out) :: t(0:nifguga)
        real(dp), intent(out) :: probWeight
        character(*), parameter :: this_routine = "forced_mixed_start"

        real(dp) :: tempWeight, bVal, tempWeight_1
        integer :: st

        ASSERT(isProperCSF_ilut(ilut))
        ASSERT(.not. isZero(ilut, excitInfo%fullStart))
        ASSERT(.not. isThree(ilut, excitInfo%fullStart))

        ! cant be 0 or zero matrix element, but cant be 3 either or only
        ! deltaB=0 branch would have non-zero matrix element and that would
        ! to a single-excitation-like DE

        st = excitInfo%fullStart
        bVal = csf_i%B_real(st)

        t = ilut

        select case (csf_i%stepvector(st))
        case (1)
            if (csf_i%B_int(st) < 2) then
                ! actually not possible in this case, as only 0 branch valid..
                probWeight = 0.0_dp
                t = 0_n_int
                return
            end if

            ! otherwise switch 1 -> 2
            clr_one(t, st)
            set_two(t, st)

            call getDoubleMatrixElement(2, 1, -1, gen_type%L, gen_type%R, &
                                        bVal, 1.0_dp, x1_element=tempWeight_1)

            ! x0 matrix element:
            tempWeight = 0.0_dp

            call setDeltaB(2, t)

        case (2)
            ! choose -2 branch 2 -> 1
            clr_two(t, st)
            set_one(t, st)

            call getDoubleMatrixElement(1, 2, -1, gen_type%L, gen_type%R, &
                                        bVal, 1.0_dp, x1_element=tempWeight_1)

            tempWeight = 0.0_dp

            call setDeltaB(-2, t)
#ifdef DEBUG_
        case default
            call stop_all(this_routine, "wrong stepvalue!")
#endif
        end select

        call encode_matrix_element(t, tempWeight_1, 2)
        call encode_matrix_element(t, tempWeight, 1)

        if (near_zero(abs(tempWeight) + abs(tempWeight_1))) then
            probWeight = 0.0_dp
            t = 0_n_int
            return
        end if

        ! and since there is no choice: branch_pgen is 1
        probWeight = 1.0_dp

    end subroutine forced_mixed_start