perform_real_time_iteration Subroutine

public subroutine perform_real_time_iteration(err)

Arguments

Type IntentOptional Attributes Name
integer, intent(out) :: err

Contents


Source Code

    subroutine perform_real_time_iteration(err)
        ! routine which performs one real-time fciqmc iteration
        implicit none
        integer, intent(out) :: err
        character(*), parameter :: this_routine = "perform_real_time_iteration"
        integer :: TotWalkersNew, run, MaxIndex
        real(dp) :: tmp_sign(lenof_sign), tau_real_tmp, tau_imag_tmp
        logical :: both, rkone, rktwo
        ! 0)
        ! do all the necessary preperation(resetting pointers etc.)
        ! concerning the statistics: i could use the "normal" iter_data
        ! for the first spawn, except change it for the new death-step, as
        ! there the particles also change from Re <-> Im

        call init_real_time_iteration(iter_data_fciqmc, second_spawn_iter_data)
        ! 1)
        ! do a "normal" spawning step and combination to y(n) + k1/2
        ! into CurrentDets:

        tau_real_tmp = tau_real
        tau_imag_tmp = tau_imag
        ! The factor corresponding to the quadratic damping is added to the initial half factor (defaults to 0)
        tau_real = (real_time_info%quad_damp_fac + 0.5d0)*tau_real
        tau_imag = (real_time_info%quad_damp_fac + 0.5d0)*tau_imag


        call first_real_time_spawn(err)
        if (catch_error()) return

        tau_real = tau_real_tmp
        tau_imag = tau_imag_tmp

        if (iProcIndex == root .and. .false.) then
            print *, "ValidSpawnedList", ValidSpawnedList
            print *, "TotParts and totDets after first spawn: ", TotParts, TotWalkers
            print *, "=========================="
        end if

! #ifdef DEBUG_
! TODO: This has to switched on again and fixed!
!         call check_update_growth(iter_data_fciqmc, "Error in first RK step")
! #endif

        ! for now update the iter data here, although in the final
        ! implementation i only should do that after the 2nd RK step
        ! or keep track of two different iter_datas for first and second
        ! spawning..

        call update_iter_data(iter_data_fciqmc)

        ! 2)
        ! reset the spawned list and do a second spawning step to create
        ! the spawend list k2
        ! but DO NOT yet recombine with stored walker list
        ! if i want to keep track of the two distinct spawns in 2 different
        ! iter_datas i probably have to reset some values before the
        ! second spawn.. to then keep the new values in the 2nd list

        call reset_spawned_list()

        NoBorn = 0
        Annihilated = 0
        NoDied = 0

        ! create a second spawned list from y(n) + k1/2
        ! have to think to exclude death_step here and store this
        ! information into the spawned k2 list..
        ! quick solution would be to loop again over reloaded y(n)
        ! and do a death step for wach walker
        call second_real_time_spawn(err)
        if (catch_error()) return

        ! 3)
        ! reload stored temp_det_list y(n) into CurrentDets
        ! have to figure out how to effectively save the previous hash_table
        ! or maybe just use two with different types of update functions..

        call reload_current_dets()

        ! 4)
        ! for the death_step for now: loop once again over the walker list
        ! and do a death_step for each walker..
        ! meh.. that seems really inefficient, better do it more cleverly
        ! in the creation of the k2 spawned list + annihilation!

!         do idet = 1, int(TotWalkers)
!
!             ! do death_step only.. maybe..
!         end do

        ! combine y(n+1) = y(n) + k2 into CurrentDets to finish time-step
        ! this should be done with a single Annihilation step between
        ! y(n) = CurrentDets and k2 = SpawnedWalkers

        ! UPDATE! have changed the 2nd diagonal event, so these particles get
        ! stored in a seperate DiagParts array -> so i have to do two
        ! annihilation events, first with the diagonal list and then with
        ! the actual spawned particles, to best mimick the old algorithm and
        ! also to correctly keep the stats of the events!

        TotWalkersNew = int(TotWalkers)

        ! also have to set the SumWalkersCyc before the "proper" annihilaiton
        do run = 1, inum_runs
            SumWalkersCyc(run) = SumWalkersCyc(run) + &
                                 sum(TotParts(min_part_type(run):max_part_type(run)))
        end do

        call DirectAnnihilation_diag(TotWalkersNew, second_spawn_iter_data)
        TotWalkersNew = int(TotWalkersNew)

        ! and then do the "normal" annihilation with the SpawnedParts array!
        ! Annihilation is done after loop over walkers
        call communicate_and_merge_spawns(MaxIndex, second_spawn_iter_data, .false.)
        call DirectAnnihilation(TotWalkersNew, MaxIndex, second_spawn_iter_data, err)
        if (catch_error()) return

#ifdef DEBUG_
        call check_update_growth(second_spawn_iter_data, "Error in second RK step")
#endif

        TotWalkers = int(TotWalkersNew)

        ! also do the update on the second_spawn_iter_data to combine both of
        ! them outside this function

        call update_iter_data(second_spawn_iter_data)

    contains
        function catch_error() result(throw)
            logical :: throw
            integer :: all_err

            call MPISumAll(err, all_err)
            err = all_err
            throw = all_err /= 0
        end function catch_error

    end subroutine perform_real_time_iteration