gen_single_back_spawn Subroutine

private subroutine gen_single_back_spawn(nI, ilutI, part_type, nJ, ilutJ, ex, tPar, pgen)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nI(nel)
integer(kind=n_int), intent(in) :: ilutI(0:niftot)
integer, intent(in) :: part_type
integer, intent(out) :: nJ(nel)
integer(kind=n_int), intent(out) :: ilutJ(0:niftot)
integer, intent(out) :: ex(2,maxExcit)
logical, intent(out) :: tPar
real(kind=dp), intent(out) :: pgen

Contents

Source Code


Source Code

    subroutine gen_single_back_spawn(nI, ilutI, part_type, nJ, ilutJ, ex, tPar, pgen)
        ! specialised single excitation routine for the back-spawn method
        ! for the moment i still have to decide, which back-spawn method is
        ! in use..
        integer, intent(in) :: nI(nel)
        integer(n_int), intent(in) :: ilutI(0:niftot)
        integer, intent(in) :: part_type
        integer, intent(out) :: nJ(nel), ex(2, maxExcit)
        integer(n_int), intent(out) :: ilutJ(0:niftot)
        logical, intent(out) :: tPar
        real(dp), intent(out) :: pgen

        integer :: elec, src, cc_index, loc, tgt
        real(dp) :: pgen_elec

        ! depending on the method we pick electrons accordingly
        if (t_back_spawn_flex) then
            elec = 1 + floor(genrand_real2_dSFMT() * nel)

            loc = check_electron_location([nI(elec), 0], 1, part_type)

            pgen_elec = 1.0_dp / real(nel, dp)
        else
            call pick_virtual_electron_single(nI, part_type, elec, pgen_elec)
            loc = -1

        end if

        src = nI(elec)

        cc_index = ClassCountInd(get_spin(src), SpinOrbSymLabel(src), &
                                 G1(src)%Ml)

        ! i have to make the logic easier here at some point..
        ! we just need to do some extensive testing and decide on one
        ! back-spawn method and abandon the rest..
        ! i hope this logic is correct: otherwise check on the bottom
        if ((t_back_spawn_occ_virt) .or. (t_back_spawn_flex .and. ( &
                                          (loc == 2 .and. occ_virt_level /= -1) .or. occ_virt_level == 2))) then

            call pick_occupied_orbital_single(ilutI, src, cc_index, part_type, pgen, tgt)

        else

            tgt = select_orb_sing(nI, ilutI, src, cc_index, pgen)

        end if
!
!         if (t_back_spawn_occ_virt) then
!             call pick_occupied_orbital_single(nI, src, cc_index, pgen, tgt)
!
!         else if (t_back_spawn_flex) then
!
!             if (loc == 2 .and. occ_virt_level /= -1) then
!                 call pick_occupied_orbital_single(nI, src, cc_index, pgen, tgt)
!
!             else
!                 if (occ_virt_level == 2) then
!                     call pick_occupied_orbital_single(nI, src, cc_index, pgen, tgt)
!                 else
!                     tgt = select_orb_sing(nI, ilutI, src, cc_index, pgen)
!                 end if
!             end if
!         else
!             tgt = select_orb_sing(nI, ilutI, src, cc_index, pgen)
!         end if

        if (tgt == 0) then
            nJ(1) = 0
            pgen = 0.0_dp
            return
        end if

        call make_single(nI, nJ, elec, tgt, ex, tPar)

        ilutJ = make_ilutJ(ilutI, ex, 1)

        pgen = pgen * pgen_elec

    end subroutine gen_single_back_spawn