gen_excit_mixed_k_space_hub_transcorr Subroutine

public subroutine gen_excit_mixed_k_space_hub_transcorr(nI, ilutI, nJ, ilutJ, exFlag, ic, ex, tParity, pgen, hel, store, run)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nI(nel)
integer(kind=n_int), intent(in) :: ilutI(0:NIfTot)
integer, intent(out) :: nJ(nel)
integer(kind=n_int), intent(out) :: ilutJ(0:NIfTot)
integer, intent(in) :: exFlag
integer, intent(out) :: ic
integer, intent(out) :: ex(2,3)
logical, intent(out) :: tParity
real(kind=dp), intent(out) :: pgen
real(kind=dp), intent(out) :: hel
type(excit_gen_store_type), intent(inout), target :: store
integer, intent(in), optional :: run

Contents


Source Code

    subroutine gen_excit_mixed_k_space_hub_transcorr(nI, ilutI, nJ, ilutJ, &
                                                     exFlag, ic, ex, tParity, pgen, hel, store, run)
        ! excitation generator, which makes doubles uniform and triples in a
        ! weighted fashion to reduce cost, but at the same time increase
        ! the time-step and the sampling
        integer, intent(in) :: nI(nel), exFlag
        integer(n_int), intent(in) :: ilutI(0:NIfTot)
        integer, intent(out) :: nJ(nel), ic, ex(2, 3)
        integer(n_int), intent(out) :: ilutJ(0:NIfTot)
        real(dp), intent(out) :: pgen
        logical, intent(out) :: tParity
        HElement_t(dp), intent(out) :: hel
        type(excit_gen_store_type), intent(inout), target :: store
        integer, intent(in), optional :: run
#ifdef DEBUG_
        character(*), parameter :: this_routine = "gen_excit_mixed_k_space_hub_transcorr"
#endif

        unused_var(exFlag)
        unused_var(store)
        unused_var(run)

        hel = h_cast(0.0_dp)
        ilutJ = 0_n_int
        ic = 0
        nJ(1) = 0

        if (genrand_real2_dsfmt() < pDoubles) then

            ! double excitation
            ic = 2

            if (genrand_real2_dsfmt() < pParallel) then

                call gen_uniform_double_para(nI, ilutI, nJ, ilutJ, ex, tParity, pgen)

                pgen = pgen * pDoubles * pParallel

            else

                call gen_uniform_double_anti(nI, ilutI, nJ, ilutJ, ex, tParity, pgen)

                pgen = pgen * pDoubles * (1.0_dp - pParallel)

            end if

        else
            ! triple excitation
            ic = 3
            call gen_triple_hubbard(nI, ilutI, nJ, ilutJ, ex, tParity, pgen)
            pgen = pgen * (1.0_dp - pDoubles)

        end if

#ifdef DEBUG_
        if (nJ(1) /= 0) then
            if (abs(pgen - calc_pgen_mixed_k_space_hub_transcorr(nI, ilutI, ex, ic)) > EPS) then
                print *, "nI: ", nI
                print *, "nJ: ", nJ
                print *, "ic: ", ic
                print *, "calc. pgen: ", calc_pgen_mixed_k_space_hub_transcorr(nI, ilutI, ex, ic)
                print *, "prd. pgen: ", pgen
                call stop_all(this_routine, "pgens wrong!")
            end if
        end if
#endif

    end subroutine gen_excit_mixed_k_space_hub_transcorr