generate_fval_energy_hist Subroutine

public subroutine generate_fval_energy_hist(hist, histEnergy, histAccRate, enPoints, accRatePoints, allHist)

Create the data written out in the histogram of shift factor over energy. 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] histEnergy on return, energy axis of the histogram @param[out] histAccRate on return, shift factor axis of the histogram @param[in] enPoints number of bins on the energy 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(enPoints,accRatePoints)
real(kind=dp), intent(out) :: histEnergy(enPoints)
real(kind=dp), intent(out) :: histAccRate(accRatePoints)
integer, intent(in) :: enPoints
integer, intent(in) :: accRatePoints
integer, intent(out) :: allHist(:,:)

Contents


Source Code

    subroutine generate_fval_energy_hist(hist, histEnergy, histAccRate, &
        enPoints, accRatePoints, allHist)
        use global_det_data, only: det_diagH, 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 histEnergy(:) : energies 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) :: enPoints, accRatePoints
        integer, intent(out) :: hist(enPoints,accRatePoints)
        integer, intent(out) :: allHist(:,:)
        real(dp), intent(out) :: histEnergy(enPoints), histAccRate(accRatePoints)
        integer :: i, run
        integer :: enInd, arInd
        real(dp) :: minEn, maxEn, enWindow, arWindow, locMinEn, locMaxEn, totSpawn

        ! get the energy window size (the acc. rate window size is just 1.0/(accRatePoints+1))
        locMinEn = det_diagH(1)
        locMaxEn = det_diagH(1)
        do i = 2, int(TotWalkers)
            if(det_diagH(i) > locMaxEn) locMaxEn = det_diagH(i)
            if(det_diagH(i) < locMinEn) locMinEn = det_diagH(i)
        end do
        ! communicate the energy window
        call MPIAllReduce(locMinEn, MPI_MIN, minEn)
        call MPIAllReduce(locMaxEn, MPI_MAX, maxEn)
        enWindow = (maxEn - minEn) / real(enPoints,dp)

        if(enWindow > eps) then
            ! set up the histogram axes
            arWindow = 1.0_dp / real(accRatePoints,dp)
            do i = 1, accRatePoints+1
                histAccRate(i) = (i-1)*arWindow
            end do

            do i = 1, enPoints+1
                histEnergy(i) = minEn + (i-1)*enWindow
            end do

            ! then, fill the histogram itself
            hist = 0
            do i = 1, int(TotWalkers)
                do run = 1, inum_runs
                    totSpawn = get_tot_spawns(i,run)
                    if(abs(totSpawn) > eps) then
                        enInd = getHistIndex(det_diagH(i), minEn, enPoints, enWindow)
                        arInd = getHistIndex(get_acc_spawns(i,run)/totSpawn, &
                            0.0_dp, accRatePoints, arWindow)
                        hist(enInd, arInd) = hist(enInd, arInd) + 1
                    end if
                end do
            end do

            ! communicate the histogram
            call MPISum(hist, allHist)
        else
            write(stdout,*) "WARNING: Empty energy histogram of acceptance rates"
        end if
    end subroutine generate_fval_energy_hist