rezero_iter_stats_each_iter Subroutine

public subroutine rezero_iter_stats_each_iter(iter_data, rdm_defs)

Arguments

Type IntentOptional Attributes Name
type(fcimc_iter_data), intent(inout) :: iter_data
type(rdm_definitions_t), intent(in) :: rdm_defs

Contents


Source Code

    subroutine rezero_iter_stats_each_iter(iter_data, rdm_defs)

        use global_det_data, only: len_av_sgn_tot
        use rdm_data, only: rdm_definitions_t

        type(fcimc_iter_data), intent(inout) :: iter_data
        type(rdm_definitions_t), intent(in) :: rdm_defs

        real(dp) :: prev_AvNoatHF(len_av_sgn_tot), AllInstNoatHF(lenof_sign)
        integer :: irdm, av_ind_1, av_ind_2, part_type

        NoInitDets = 0_int64
        NoNonInitDets = 0_int64
        NoInitWalk = 0.0_dp
        NoNonInitWalk = 0.0_dp
        InitRemoved = 0_int64

        ! replica-initiator info

        NoAborted = 0.0_dp
        NoRemoved = 0.0_dp
        NoatHF = 0.0_dp
        NoatDoubs = 0.0_dp
        if (tLogEXLEVELStats) EXLEVEL_WNorm = 0.0_dp

        iter_data%nborn = 0.0_dp
        iter_data%ndied = 0.0_dp
        iter_data%nannihil = 0.0_dp
        iter_data%naborted = 0.0_dp
        iter_data%nremoved = 0.0_dp

        call InitHistMin()

        if (tFillingStochRDMonFly) then
            associate(ind => rdm_defs%sim_labels)

                call MPISumAll(InstNoatHF, AllInstNoAtHF)
                InstNoAtHF = AllInstNoAtHF

                if (tFullHFAv) then
                    Prev_AvNoatHF = AvNoatHF

                    do irdm = 1, rdm_defs%nrdms
                        if (abs(IterRDM_HF(irdm)) > 1.0e-12_dp) then
                            AvNoatHF(irdm) = ((real((Iter + PreviousCycles - IterRDM_HF(irdm)), dp) * Prev_AvNoatHF(irdm)) &
                                              + InstNoatHF(irdm)) / real((Iter + PreviousCycles - IterRDM_HF(irdm)) + 1, dp)
                        end if
                    end do

                else
                    if (((Iter + PreviousCycles - IterRDMStart) > 0) .and. &
                        & (mod(((Iter - 1) + PreviousCycles - IterRDMStart + 1), RDMEnergyIter) == 0)) then
                        ! The previous iteration was one where we added in diagonal
                        ! elements To keep things unbiased, we need to set up a new
                        ! averaging block now.
                        do irdm = 1, rdm_defs%nrdms
                            av_ind_1 = irdm * 2 - 1
                            av_ind_2 = irdm * 2

                            AvNoAtHF(av_ind_1) = InstNoAtHF(ind(1, irdm))
                            IterRDM_HF(av_ind_1) = Iter + PreviousCycles
                            AvNoAtHF(av_ind_2) = InstNoAtHF(ind(2, irdm))
                            IterRDM_HF(av_ind_2) = Iter + PreviousCycles
                        end do
                    else
                        do irdm = 1, rdm_defs%nrdms
                            ! The indicies of the first and second replicas in this
                            ! particular pair, in the *average* sign arrays.
                            av_ind_1 = irdm * 2 - 1
                            av_ind_2 = irdm * 2

                            if ((abs(InstNoAtHF(ind(1, irdm))) < 1.0e-12_dp .and. abs(IterRDM_HF(av_ind_1)) > 1.0e-12_dp) .or. &
                                (abs(InstNoAtHF(ind(2, irdm))) < 1.0e-12_dp .and. abs(IterRDM_HF(av_ind_2)) > 1.0e-12_dp)) then
                                ! At least one of the populations has just become
                                ! zero. Start a new averaging block.
                                IterRDM_HF(av_ind_1) = Iter + PreviousCycles
                                IterRDM_HF(av_ind_2) = Iter + PreviousCycles
                                AvNoatHF(av_ind_1) = InstNoAtHF(ind(1, irdm))
                                AvNoatHF(av_ind_2) = InstNoAtHF(ind(2, irdm))
                                if (abs(InstNoAtHF(ind(1, irdm))) < 1.0e-12_dp) IterRDM_HF(av_ind_1) = 0
                                if (abs(InstNoAtHF(ind(2, irdm))) < 1.0e-12_dp) IterRDM_HF(av_ind_2) = 0

                            else if ((abs(InstNoAtHF(ind(1, irdm))) > 1.0e-12_dp .and. abs(IterRDM_HF(av_ind_1)) < 1.0e-12_dp) .or. &
                                     (abs(InstNoAtHF(ind(2, irdm))) > 1.0e-12_dp .and. abs(IterRDM_HF(av_ind_2)) < 1.0e-12_dp)) then
                                ! At least one of the populations has just
                                ! become occupied. Start a new block here.
                                IterRDM_HF(av_ind_1) = Iter + PreviousCycles
                                IterRDM_HF(av_ind_2) = Iter + PreviousCycles
                                AvNoAtHF(av_ind_1) = InstNoAtHF(ind(1, irdm))
                                AvNoAtHF(av_ind_2) = InstNoAtHF(ind(2, irdm))
                                if (abs(InstNoAtHF(ind(1, irdm))) < 1.0e-12_dp) IterRDM_HF(av_ind_1) = 0
                                if (abs(InstNoAtHF(ind(2, irdm))) < 1.0e-12_dp) IterRDM_HF(av_ind_2) = 0
                            else
                                Prev_AvNoatHF = AvNoatHF

                                if (abs(IterRDM_HF(av_ind_1)) > 1.0e-12_dp) then
                                    AvNoatHF(av_ind_1) = ((real((Iter + PreviousCycles - IterRDM_HF(av_ind_1)), dp) &
                                                           * Prev_AvNoatHF(av_ind_1)) + InstNoatHF(ind(1, irdm))) &
                                                         / real((Iter + PreviousCycles - IterRDM_HF(av_ind_1)) + 1, dp)
                                end if
                                if (abs(IterRDM_HF(av_ind_2)) > 1.0e-12_dp) then
                                    AvNoatHF(av_ind_2) = ((real((Iter + PreviousCycles - IterRDM_HF(av_ind_2)), dp) &
                                                           * Prev_AvNoatHF(av_ind_2)) + InstNoatHF(ind(2, irdm))) &
                                                         / real((Iter + PreviousCycles - IterRDM_HF(av_ind_2)) + 1, dp)
                                end if

                            end if
                        end do

                    end if
                end if

            end associate
        end if

        HFInd = 0

        min_trial_ind = 1
        min_conn_ind = 1

        trial_num_inst = 0.0_dp
        trial_denom_inst = 0.0_dp

    end subroutine rezero_iter_stats_each_iter