generate_nGAS_single Subroutine

private subroutine generate_nGAS_single(this, nI, ilutI, nJ, ilutJ, ex, par, pgen)

Type Bound

GAS_disc_ExcGenerator_t

Arguments

Type IntentOptional Attributes Name
class(GAS_disc_ExcGenerator_t), intent(in) :: this
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(out) :: ex(2,2)
logical, intent(out) :: par
real(kind=dp), intent(out) :: pgen

Contents

Source Code


Source Code

    subroutine generate_nGAS_single(this, nI, ilutI, nJ, ilutJ, ex, par, pgen)
        class(GAS_disc_ExcGenerator_t), intent(in) :: this
        integer, intent(in) :: nI(nel)
        integer(n_int), intent(in) :: ilutI(0:NIfTot)
        integer, intent(out) :: nJ(nel), ex(2, 2)
        integer(n_int), intent(out) :: ilutJ(0:NIfTot)
        logical, intent(out) :: par
        real(dp), intent(out) :: pgen

        integer :: elec, tgt, src, spin_idx
        real(dp) :: r

        ! we assume that each active space has a possible single excitation
        ! (this is very mild because active spaces without such are trivial)
        ! -> pick any random electron
        r = genrand_real2_dSFMT()
        elec = int(r * nel) + 1
        src = nI(elec)
        ! adjust pgen
        pgen = 1.0_dp / nel
        ! from the same active space, get a hole
        spin_idx = get_spin(src)
        tgt = this%pick_weighted_hole(nI, Excite_1_t(src), spin_idx, this%GAS_table(src), pgen)

        if (tgt == 0) then
            nJ(1) = 0
            ilutJ = 0_n_int
            return
        end if

        call make_single(nI, nJ, elec, tgt, ex, par)
        ilutJ = ilutI
        clr_orb(ilutJ, src)
        set_orb(ilutJ, tgt)
    end subroutine generate_nGAS_single