attempt_die_realtime Function

public function attempt_die_realtime(Kii, RealwSign, walkExcitLevel) result(ndie)

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(in) :: Kii
real(kind=dp), intent(in), dimension(lenof_sign) :: RealwSign
integer, intent(in) :: walkExcitLevel

Return Value real(kind=dp), dimension(lenof_sign)


Contents

Source Code


Source Code

    function attempt_die_realtime(Kii, RealwSign, walkExcitLevel) &
        result(ndie)
        ! also write a function, which calculates the new "signs"(weights) of
        ! the real and complex walkers for the diagonal death/cloning step
        ! since i need that for both 1st and 2nd loop of RK, but at different
        ! points
        implicit none
        real(dp), dimension(lenof_sign), intent(in) :: RealwSign
        real(dp), intent(in) :: Kii
        real(dp), dimension(lenof_sign) :: ndie
        integer, intent(in) :: walkExcitLevel
        integer :: run
        real(dp) :: fac(lenof_sign), rat, r

        character(*), parameter :: this_routine = "attempt_die_realtime"

        ! do it dependent on if there is damping, since otherwise i can save
        ! some effort..

        ! do i need a minus sign here?? just from the convention
        ! in the rest of the neci code? yes!
        ! rmneci_setup: there is no reason to use an imaginary shift
        ! except if when dealing with rotated times)
        ! TODO: Energies should be taken with respect to the N-particle ground state energy

        ! important: the matrix element Kii does not contain the reference energy,
        ! therefore it has to be added manually
        do run = 1, inum_runs
            fac(min_part_type(run)) = tau_real * (Kii + Hii - gs_energy(run))
            fac(max_part_type(run)) = 0.0_dp
        end do

        if (abs(real_time_info%damping) < EPS .and. .not. t_rotated_time) then
            if (any(fac > 1.0_dp)) then
                if (any(fac > 2.0_dp)) then
                    if ((tau_search_method /= possible_tau_search_methods%OFF) .or. t_scale_tau_to_death) then
                        ! If we are early in the calculation, and are using tau
                        ! searching, then this is not a big deal. Just let the
                        ! searching deal with it
                        write(stdout, '("** WARNING ** Death probability > 2: Algorithm unstable.")')
                        write(stdout, '("** WARNING ** Truncating spawn to ensure stability")')
                        do run = 1, lenof_sign
                            fac(run) = min(2.0_dp, fac(run))
                        end do
                    else
                        call stop_all(this_routine, "Death probability > 2: Algorithm unstable. Reduce timestep.")
                    end if
                else
                    do run = 1, inum_runs
                        write(stdout, '(1X,f13.7)', advance='no') fac(run)
                    end do
                    write(stdout, '()')
                end if
            end if
            do run = 1, inum_runs
                if (((tRealCoeffByExcitLevel .and. (WalkExcitLevel <= RealCoeffExcitThresh)) &
                     .or. tAllRealCoeff)) then
                    ndie(min_part_type(run)) = -fac(min_part_type(run)) &
                                               * realwSign(max_part_type(run))
                    ! and - from Re -> Im
                    ! does this give the correct sign compared to the parent sign?
                    ! additional -1 is added in postprocessing to convert ndie -> nborn
                    ndie(max_part_type(run)) = fac(min_part_type(run)) &
                                               * realwSign(min_part_type(run))

                else
                    ! if not exact i have to round stochastically
                    ! Im -> Re
                    rat = -fac(min_part_type(run)) * RealwSign(max_part_type(run))

                    ndie(min_part_type(run)) = real(int(rat), dp)
                    rat = rat - ndie(min_part_type(run))

                    r = genrand_real2_dSFMT()
                    if (abs(rat) > r) ndie(min_part_type(run)) = &
                        ndie(min_part_type(run)) + real(nint(sign(1.0_dp, rat)), dp)

                    ! Re -> Im
                    rat = fac(min_part_type(run)) * RealwSign(min_part_type(run))
                    ndie(max_part_type(run)) = real(int(rat), dp)
                    rat = rat - ndie(max_part_type(run))
                    r = genrand_real2_dSFMT()
                    if (abs(rat) > r) ndie(max_part_type(run)) = &
                        ndie(max_part_type(run)) + real(nint(sign(1.0_dp, rat)), dp)
                end if
            end do

            ! here i only have influence from the other type of walkers, which
            ! is already stored in the ndie() array
        else
            ! there is an imaginary energy damping factor.. -> rotated time
            ! so i have mixed influences for the diagonal step
            ! n_i' = -e n_i + H_ii c_i
            ! c_i' = -e c_i - H_ii n_i

            ! temporarily use the 2 entries of fac to store both influences
            ! not sure about the sign of this fac factor but the 2nd entry
            ! as it is implemented right now has to have the opposite sign
            ! of the first on above
            ! the diagnonal matrix elements are always real, so there is
            ! no contribution from tau_real via Kii (and no influence of
            ! tau_imag on fac(min_part_type(run))
            do run = 1, inum_runs
                ! tau_real and tau_imag have opposite signs -> shift changed sign
                ! when moved from real to imaginary part
                fac(max_part_type(run)) = tau_real * (real_time_info%damping) &
                                          + tau_imag * (Kii + Hii - gs_energy(run) - DiagSft(run))
            end do

            ! and also about the fac restrictions.. for now but it here anyway..
            if (any(fac > 1.0_dp)) then
                if (any(fac > 2.0_dp)) then
                    if ((tau_search_method /= possible_tau_search_methods%OFF) .or. t_scale_tau_to_death) then
                        ! If we are early in the calculation, and are using tau
                        ! searching, then this is not a big deal. Just let the
                        ! searching deal with it
                        write(stdout, '("** WARNING ** Death probability > 2: Algorithm unstable.")')
                        write(stdout, '("** WARNING ** Truncating spawn to ensure stability")')
                        do run = 1, lenof_sign
                            fac(run) = min(2.0_dp, fac(run))
                        end do
                    else
                        call stop_all(this_routine, "Death probability > 2: Algorithm unstable. Reduce timestep.")
                    end if
                else
                    write(stdout, *) "Warning, diagonal spawn probability > 1. Reduce timestep."
                    do run = 1, inum_runs
                        write(stdout, '(1X,f13.7)', advance='no') fac(run)
                    end do
                    write(stdout, '()')
                end if
            end if
            do run = 1, inum_runs
                if (((tRealCoeffByExcitLevel .and. (WalkExcitLevel <= RealCoeffExcitThresh)) &
                     .or. tAllRealCoeff) .and. .not. tVerletSweep) then
                    ! this gives a huge overhead in the verlet scheme, as all diagonal
                    ! events are classified as valid -> #spawns ~ #dets
                    ! i exact.. just get the weights, this should also still work.
                    ! but have to exchange the weights to come from the other
                    ! type of particles..
                    ! the number of deaths has +sign from Im -> Re
                    ! can i just add the other contribution here?
                    ! and also include the sign of the parent occupation here
                    ! already.
                    ndie(min_part_type(run)) = -fac(min_part_type(run)) &
                                               * realwSign(max_part_type(run)) - &
                                               fac(max_part_type(run)) * realwSign(min_part_type(run))
                    ! and - from Re -> Im
                    ! does this give the correct sign compared to the parent sign?
                    ndie(max_part_type(run)) = fac(min_part_type(run)) * &
                                               (realwSign(min_part_type(run))) - &
                                               fac(max_part_type(run)) * (RealwSign(max_part_type(run)))

                else
                    ! if not exact i have to round stochastically
                    ! is this ok here to just add the second contribution? todo
                    !  -> Re
                    rat = -fac(min_part_type(run)) * RealwSign(max_part_type(run)) &
                          - fac(max_part_type(run)) * RealwSign(min_part_type(run))

                    ndie(min_part_type(run)) = real(int(rat), dp)
                    rat = rat - ndie(min_part_type(run))

                    r = genrand_real2_dSFMT()
                    if (abs(rat) > r) ndie(min_part_type(run)) = ndie(min_part_type(run)) &
                                                                 + real(nint(sign(1.0_dp, rat)), dp)

                    !  -> Im
                    rat = fac(min_part_type(run)) * RealwSign(min_part_type(run)) &
                          - fac(max_part_type(run)) * RealwSign(max_part_type(run))

                    ndie(max_part_type(run)) = real(int(rat), dp)
                    rat = rat - ndie(max_part_type(run))

                    r = genrand_real2_dSFMT()
                    if (abs(rat) > r) ndie(max_part_type(run)) = &
                        ndie(max_part_type(run)) + real(nint(sign(1.0_dp, rat)), dp)

                end if
            end do
        end if

    end function attempt_die_realtime