generate_fval_pop_hist Subroutine

public subroutine generate_fval_pop_hist(hist, histPop, histAccRate, popPoints, accRatePoints, allHist)

Create the data written out in the histogram of shift factor over population. The generated data can be passed to print_2d_hist. This is a synchronizing routine. @param[out] hist on return, histogram data of this proc only @param[out] histPop on return, population axis of the histogram @param[out] histAccRate on return, shift factor axis of the histogram @param[in] popPoints number of bins on the population axis @param[in] accRatePoints number of bins on the shift factor axis @param[out] allHist on return, histogram data over all procs

Arguments

Type IntentOptional Attributes Name
integer, intent(out) :: hist(popPoints,accRatePoints)
real(kind=dp), intent(out) :: histPop(popPoints)
real(kind=dp), intent(out) :: histAccRate(accRatePoints)
integer, intent(in) :: popPoints
integer, intent(in) :: accRatePoints
integer, intent(out) :: allHist(:,:)

Contents


Source Code

    subroutine generate_fval_pop_hist(hist, histPop, histAccRate, &
        popPoints, accRatePoints, allHist)
        use global_det_data, only: get_acc_spawns, get_tot_spawns
        ! count the acceptance ratio per energy and create
        ! a histogram hist(:,:) with axes histEnergy(:), histAccRate(:)
        ! entries of hist(:,:) : numer of occurances
        ! entries of histPop(:) : populations of histogram entries (with tolerance, first dimension)
        ! entries of histAccRate(:) : acc. rates of histogram entries (second dimension)
        implicit none
        ! number of energy/acc rate windows in the histogram
        integer, intent(in) :: popPoints, accRatePoints
        integer, intent(out) :: hist(popPoints,accRatePoints)
        integer, intent(out) :: allHist(:,:)
        real(dp), intent(out) :: histPop(popPoints), histAccRate(accRatePoints)
        integer :: i, run
        integer :: popInd, arInd
        real(dp) :: maxPop, locMaxPop, totSpawn, pop
        real(dp), dimension(lenof_sign) :: sgn

        ! The bins of the population has width of 1, except the last bin which
        ! constains all populations larger than popPoints-1
        hist = 0
        locMaxPop = 0.0
        do i = 1, int(TotWalkers)
            call extract_sign(CurrentDets(:,i), sgn)
            do run = 1, inum_runs
                pop = mag_of_run(sgn,run)
                if(pop > locMaxPop) locMaxPop = pop
                totSpawn = get_tot_spawns(i,run)
                if(abs(totSpawn) > eps) then
                    popInd = getHistIndex(pop, 0.0_dp, popPoints, 1.0_dp)
                    arInd = getHistIndex(get_acc_spawns(i,run)/totSpawn, &
                        0.0_dp, accRatePoints, 1.0_dp/real(accRatePoints,dp))
                    hist(popInd, arInd) = hist(popInd, arInd) + 1
                end if
            end do
        end do

        call MPIAllReduce(locMaxPop, MPI_MAX, maxPop)
        do i = 1, popPoints+1
            histPop(i) = real(i-1)
        end do
        if(maxPop>histPop(popPoints+1))  histPop(popPoints+1)= maxPop

        do i = 1, accRatePoints+1
            histAccRate(i) = (i-1)/real(accRatePoints,dp)
        end do

        ! communicate the histogram
        call MPISum(hist, allHist)

    end subroutine generate_fval_pop_hist