walker_death_realtime Subroutine

public subroutine walker_death_realtime(iter_data, DetCurr, iLutCurr, Kii, RealwSign, DetPosition, walkExcitLevel)

Arguments

Type IntentOptional Attributes Name
type(fcimc_iter_data), intent(inout) :: iter_data
integer, intent(in) :: DetCurr(nel)
integer(kind=n_int), intent(in) :: iLutCurr(0:niftot)
real(kind=dp), intent(in) :: Kii
real(kind=dp), intent(in), dimension(lenof_sign) :: RealwSign
integer, intent(in) :: DetPosition
integer, intent(in) :: walkExcitLevel

Contents

Source Code


Source Code

    subroutine walker_death_realtime(iter_data, DetCurr, iLutCurr, Kii, RealwSign, &
                                     DetPosition, walkExcitLevel)
        ! need new walker_death routine for the real-time fciqmc, since i
        ! have complex diagonal influence: due to real-time formulation and
        ! the use of a imaginary energy damping factot -ie
        ! (n_i + ic_i)' = -i(H_ii -E0 - ie)(n_i + ic_i)
        ! n_i' = -e n_i + (H_ii - E0) c_i
        ! c_i' = -e c_i - (H_ii - E0) n_i
        ! so the complex and imaginary walkers mix in the diagonal death
        ! step also! but with opposite sign between Re <-> Im spawn
        ! so the 'death' as annihilation only happens after 2 steps..

        ! i may also need a second version of this, where the diagonal step
        ! is no directly executed on the current list, but builds it into
        ! the spawned list array, since i need the original y(n) list to
        ! combine with k2: y(n+1) = y(n) + k2
        ! i guess that could be done..
        ! this routine is now only used in the first loop of the 2nd order
        ! RK! in the 2nd loop i store the diagonal events in a specific
        ! list, which later gets merged with the original y(n) list with
        ! the Annihilation routine

        integer, intent(in) :: DetCurr(nel)
        real(dp), dimension(lenof_sign), intent(in) :: RealwSign
        integer(kind=n_int), intent(in) :: iLutCurr(0:niftot)
        real(dp), intent(in) :: Kii
        integer, intent(in) :: DetPosition
        type(fcimc_iter_data), intent(inout) :: iter_data
        real(dp), dimension(lenof_sign) :: ndie
        real(dp), dimension(lenof_sign) :: CopySign
        integer, intent(in) :: walkExcitLevel
        integer :: i, j
        real(dp) :: r

        character(len=*), parameter :: t_r = "walker_death_realtime"

        integer(n_int), pointer :: ilut(:)

        ilut => CurrentDets(:, DetPosition)

        ndie = attempt_die_realtime(Kii, realwSign, walkExcitLevel)

        ! this routine only gets called in the first runge-kutta step ->
        ! so only update the stats for the first here!
        do i = 1, lenof_sign
            ! check if the parent and ndie have the same sign
            if (sign(1.0_dp, RealwSign(i)) .isclose.sign(1.0_dp, ndie(i))) then
                ! then the entries in ndie kill the parent, but only maximally
                ! the already occupying walkers can get killed
                iter_data%ndied(i) = iter_data%ndied(i) + &
                                     abs(min(abs(RealwSign(i)), abs(ndie(i))))

                ! if ndie is bigger than the original occupation i am actually
                ! spawning 'anti-particles' which i have to count as born
                ! and reduce the ndied number.. or not?
                ! hm the old code is actually not counting births, due to
                ! the shift.. interesting.. but just subtracts that from
                ! the ndied quantity...
                iter_data%nborn(i) = iter_data%nborn(i) + &
                                     max(abs(ndie(i)) - abs(RealwSign(i)), 0.0_dp)

            else
                ! if they have opposite sign, as in the original algorithm
                ! reduce the number of ndied by that amount
                iter_data%ndied(i) = iter_data%ndied(i) - abs(ndie(i))

            end if
        end do

        CopySign = RealwSign - nDie

        if (any(.not. near_zero(CopySign))) then
            ! For the hashed walker main list, the particles don't move.
            ! Therefore just adjust the weight.
            call encode_sign(CurrentDets(:, DetPosition), CopySign)
            ! rotated_time_setup: there is only one initiator flag for
            ! both complex and real walkers -> initiator flag is kept
        else

            if (tTruncInitiator) then
                ! All particles on this determinant have gone. If the determinant was an initiator, update the stats
                do j = 1, inum_runs
                    if (test_flag(iLutCurr, get_initiator_flag_by_run(j))) then
                        NoAddedInitiators(j) = NoAddedInitiators(j) - 1
                    end if
                end do
            end if

            ! Remove the determinant from the indexing list
            call remove_hash_table_entry(HashIndex, DetCurr, DetPosition)
            ! Add to the "freeslot" list
            iEndFreeSlot = iEndFreeSlot + 1
            FreeSlot(iEndFreeSlot) = DetPosition
            ! Encode a null det to be picked up
            call encode_sign(CurrentDets(:, DetPosition), null_part)
        end if

    end subroutine walker_death_realtime