create_crude_single Subroutine

public subroutine create_crude_single(ilut, csf_i, exc, branch_pgen, excitInfo_in)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(in) :: ilut(0:nifguga)
type(CSF_Info_t), intent(in) :: csf_i
integer(kind=n_int), intent(inout) :: exc(0:nifguga)
real(kind=dp), intent(out) :: branch_pgen
type(ExcitationInformation_t), intent(in), optional :: excitInfo_in

Contents

Source Code


Source Code

    subroutine create_crude_single(ilut, csf_i, exc, branch_pgen, excitInfo_in)
        integer(n_int), intent(in) :: ilut(0:nifguga)
        type(CSF_Info_t), intent(in) :: csf_i
        integer(n_int), intent(inout) :: exc(0:nifguga)
        real(dp), intent(out) :: branch_pgen
        type(ExcitationInformation_t), intent(in), optional :: excitInfo_in
        character(*), parameter :: this_routine = "create_crude_single"

        type(ExcitationInformation_t) :: excitInfo
        integer :: elec, orb
        real(dp) :: r

#ifdef DEBUG_
        ! also assert we are not calling it for a weight gen. accidently
        ASSERT(excitInfo%currentGen /= 0)
        ASSERT(isProperCSF_ilut(ilut))
        ! also check if calculated b vector really fits to ilut
        ASSERT(all(csf_i%B_real.isclose.calcB_vector_ilut(ilut(0:nifd))))
        if (excitInfo%currentGen == gen_type%R) then
            ASSERT(.not. isThree(ilut, excitInfo%fullStart))
        else if (excitInfo%currentGen == gen_type%L) then
            ASSERT(.not. isZero(ilut, excitInfo%fullStart))
        end if
#endif

        if (present(excitInfo_in)) then
            excitInfo = excitInfo_in
        else
            call stop_all(this_routine, "not yet implemented without excitInfo_in")
        end if

        elec = excitInfo%j
        orb = excitInfo%i

        exc = ilut

        select case (csf_i%stepvector(elec))

        case (0)
            call stop_all(this_routine, "empty orbital picked as electron!")

        case (1)

            clr_orb(exc, 2 * elec - 1)
            branch_pgen = 1.0_dp

            if (csf_i%stepvector(orb) == 0) then

                set_orb(exc, 2 * orb - 1)

            else if (csf_i%stepvector(orb) == 1) then

                branch_pgen = 0.0_dp

            else if (csf_i%stepvector(orb) == 2) then

                set_orb(exc, 2 * orb - 1)

            end if

        case (2)

            clr_orb(exc, 2 * elec)
            branch_pgen = 1.0_dp

            if (csf_i%stepvector(orb) == 0) then

                set_orb(exc, 2 * orb)

            else if (csf_i%stepvector(orb) == 1) then

                set_orb(exc, 2 * orb)

            else if (csf_i%stepvector(orb) == 2) then

                branch_pgen = 0.0_dp

            end if

        case (3)

            if (csf_i%stepvector(orb) == 0) then

                ! here i have to decide..
                branch_pgen = 0.5_dp

                r = genrand_real2_dSFMT()

                if (r < 0.5_dp) then

                    ! 1 -> 2
                    clr_orb(exc, 2 * elec)
                    set_orb(exc, 2 * orb)

                else
                    ! 2 -> 1
                    clr_orb(exc, 2 * elec - 1)
                    set_orb(exc, 2 * orb - 1)

                end if

            else
                branch_pgen = 1.0_dp

                if (csf_i%stepvector(orb) == 1) then

                    clr_orb(exc, 2 * elec)
                    set_orb(exc, 2 * orb)

                else if (csf_i%stepvector(orb) == 2) then

                    clr_orb(exc, 2 * elec - 1)
                    set_orb(exc, 2 * orb - 1)

                end if
            end if
        end select

        if (.not. isProperCSF_ilut(exc)) then
            ! i have to check if i created a proper CSF
            branch_pgen = 0.0_dp
        end if

    end subroutine create_crude_single