calcFullStopR2L_stochastic Subroutine

public subroutine calcFullStopR2L_stochastic(ilut, csf_i, excitInfo, t, pgen, posSwitches, negSwitches, opt_weight)

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) :: t(0:nifguga)
real(kind=dp), intent(out) :: pgen
real(kind=dp), intent(in) :: posSwitches(nSpatOrbs)
real(kind=dp), intent(in) :: negSwitches(nSpatOrbs)
type(WeightObj_t), intent(in), optional :: opt_weight

Contents


Source Code

    subroutine calcFullStopR2L_stochastic(ilut, csf_i, excitInfo, t, pgen, &
                                          posSwitches, negSwitches, opt_weight)
        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) :: t(0:nifguga)
        real(dp), intent(out) :: pgen
        real(dp), intent(in) :: posSwitches(nSpatOrbs), negSwitches(nSpatOrbs)
        type(WeightObj_t), intent(in), optional :: opt_weight

        type(WeightObj_t) :: weights
        integer :: st, se, en, i, j, k, l, elecInd, holeInd
        real(dp) :: branch_pgen, &
                    temp_pgen, p_orig, rdm_mat, orb_pgen
        HElement_t(dp) :: integral

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

        ! init weights
        if (present(opt_weight)) then
            weights = opt_weight
        else
            if (t_approx_exchange .or. (t_approx_exchange_noninits .and. (.not. is_init_guga))) then
                weights = init_forced_end_semistart_weight(csf_i, se, en, negSwitches(se), &
                                                           posSwitches(se), csf_i%B_real(se))
            else
                weights = init_semiStartWeight(csf_i, se, en, negSwitches(se), &
                                               posSwitches(se), csf_i%B_real(se))
            end if
        end if

        ! create start
        call createStochasticStart_single(ilut, csf_i, excitInfo, weights, posSwitches, &
                                          negSwitches, t, branch_pgen)

        ! in case of early exit pgen should be set to 0
        pgen = 0.0_dp
        ! check validity

        check_abort_excit(branch_pgen, t)

        do i = st + 1, se - 1
            call singleStochasticUpdate(ilut, csf_i, i, excitInfo, weights, posSwitches, &
                                        negSwitches, t, temp_pgen)
            branch_pgen = branch_pgen * temp_pgen
            ! check validity

            check_abort_excit(branch_pgen, t)
        end do

        ! do the specific semi-start
        weights = weights%ptr

        call calcLoweringSemiStartStochastic(ilut, csf_i, excitInfo, weights, negSwitches, &
                                             posSwitches, t, branch_pgen)

        ! check validity
        check_abort_excit(branch_pgen, t)

        ! do the specific double update to ensure a switch
        ! although switch can also happen at end only...
        ! actually that would be, in the full-stop case, temporary measure...
        ! but would unjust favor certain types of excitations..
        do i = se + 1, en - 1
            call doubleUpdateStochastic(ilut, csf_i, i, excitInfo, &
                                        weights, negSwitches, posSwitches, t, branch_pgen)
            ! also here there has to be a switch at some point so check x1
            if (near_zero(extract_matrix_element(t, 2)) .or. near_zero(branch_pgen)) then
                t = 0_n_int
                return
            end if
        end do

        call mixedFullStopStochastic(ilut, csf_i, excitInfo, t)

        ! check if matrix element is non-zero and if a switch happened
        if (.not. near_zero(extract_matrix_element(t, 1))) then
            t = 0_n_int
            return
        end if
        if (t_approx_exchange .or. (t_approx_exchange_noninits .and. (.not. is_init_guga))) then
            ! make it crude for now, that we only check if the delta B value
            ! is non-zero at the end, otherwise abort this spawn..
            if (getDeltaB(t) == 0) then
                t = 0_n_int
                return
            end if
        end if

        if (near_zero(extract_matrix_element(t, 2))) then
            t = 0_n_int
            return
        end if

        ! if we do RDMs also store the x0 and x1 coupling coeffs
        ! and I need to do it before the routines below since excitInfo
        ! gets changed there
        if (tFillingStochRDMOnFly) then
            ! i need to unbias against the total pgen later on in the
            ! RDM sampling otherwise the rdm-bias factor is not correct!
            ! encode the necessary information in the rdm-matele!
            i = excitInfo%i
            j = excitInfo%j
            k = excitInfo%k
            l = excitInfo%l
            elecInd = se
            holeInd = st
            rdm_mat = extract_matrix_element(t, 2)
            call calc_orbital_pgen_contrib_end(&
                    csf_i, [2 * elecInd, 2 * en], holeInd, orb_pgen)
            p_orig = orb_pgen * branch_pgen / real(ElecPairs, dp)
            if (csf_i%stepvector(elecInd) == 3) p_orig = p_orig * 2.0_dp

        end if

        ! the x1-element is still encoded in the second entry..
        ! move it to the first elemen
        call encode_matrix_element(t, extract_matrix_element(t, 2), 1)

        global_excitInfo = excitInfo

        if (t_approx_exchange .or. (t_approx_exchange_noninits .and. (.not. is_init_guga))) then

            call calc_mixed_end_contr_approx(t, csf_i, excitInfo, integral)
            pgen = branch_pgen

        else
            call calc_mixed_end_contr_integral(ilut, csf_i, t, excitInfo, &
                integral)
            call calc_mixed_end_contr_pgen(ilut, csf_i, t, excitInfo, &
                branch_pgen, pgen)

!             call calc_mixed_end_contr_sym(ilut, csf_i, t, excitInfo, branch_pgen, pgen, integral)
        end if

        if (tFillingStochRDMOnFly) then
            if (.not. near_zero(p_orig)) then
                call encode_stochastic_rdm_info(GugaBits, t, rdm_ind= &
                        contract_2_rdm_ind(i, j, k, l, excit_lvl=2, &
                                           excit_typ=excitInfo%typ), x0=0.0_dp, &
                                                x1=rdm_mat * pgen / p_orig)
            end if
        end if

        call encode_matrix_element(t, 0.0_dp, 2)
        call update_matrix_element(t, (get_umat_el(en, st, se, en) + &
                                       get_umat_el(st, en, en, se)) / 2.0_dp + integral, 1)

    end subroutine calcFullStopR2L_stochastic