fill_rdm_diag_wrapper Subroutine

public subroutine fill_rdm_diag_wrapper(rdm_defs, spawn, one_rdms, ilut_list, ndets, tNonInit, tLagrCorr)

Arguments

Type IntentOptional Attributes Name
type(rdm_definitions_t), intent(in) :: rdm_defs
type(rdm_spawn_t), intent(inout) :: spawn
type(one_rdm_t), intent(inout) :: one_rdms(:)
integer(kind=n_int), intent(in) :: ilut_list(:,:)
integer, intent(in) :: ndets
logical, intent(in), optional :: tNonInit
logical, intent(in), optional :: tLagrCorr

Contents

Source Code


Source Code

    subroutine fill_rdm_diag_wrapper(rdm_defs, spawn, one_rdms, ilut_list, ndets, tNonInit, &
                                     tLagrCorr)

        ! Loop over all states in ilut_list and see if any signs have just
        ! become unoccupied or become reoccupied. In which case, we have
        ! started a new averaging block, so we need to add in the
        ! contributions from the last block to the corresponding RDMs.

        use bit_rep_data, only: extract_sign
        use bit_reps, only: all_runs_are_initiator
        use CalcData, only: tPairedReplicas
        use global_det_data, only: get_iter_occ_tot, get_av_sgn_tot
        use global_det_data, only: len_av_sgn_tot, len_iter_occ_tot
        use rdm_data, only: one_rdm_t, rdm_definitions_t

        type(rdm_definitions_t), intent(in) :: rdm_defs
        type(rdm_spawn_t), intent(inout) :: spawn
        type(one_rdm_t), intent(inout) :: one_rdms(:)
        integer(n_int), intent(in) :: ilut_list(:, :)
        integer, intent(in) :: ndets
        logical, intent(in), optional :: tNonInit, tLagrCorr

        integer :: idet, irdm, av_ind_1, av_ind_2
        real(dp) :: curr_sign(lenof_sign), adapted_sign(len_av_sgn_tot)
        real(dp) :: av_sign(len_av_sgn_tot), iter_occ(len_iter_occ_tot)
        logical :: tAllContribs, tLC

        def_default(tAllContribs, tNonInit, .true.)
        def_default(tLC, tLagrCorr, .true.)

        associate(ind => rdm_defs%sim_labels)

            do idet = 1, ndets

                call extract_sign(ilut_list(:, idet), curr_sign)
                ! All average sign from all RDMs.
                av_sign = get_av_sgn_tot(idet)
                ! The iteration on which each replica became occupied.
                iter_occ = get_iter_occ_tot(idet)

                adapted_sign = 0.0_dp

                do irdm = 1, rdm_defs%nrdms
                    ! The indicies of the first and second replicas in this
                    ! particular pair, in the *average* sign arrays (and
                    ! therefore also for the iter_occ array).
                    av_ind_1 = irdm * 2 - 1
                    av_ind_2 = irdm * 2

                    if ((abs(curr_sign(ind(1, irdm))) < 1.0e-10_dp .and. abs(iter_occ(av_ind_1)) > 1.0e-10_dp) .or. &
                        (abs(curr_sign(ind(2, irdm))) < 1.0e-10_dp .and. abs(iter_occ(av_ind_2)) > 1.0e-10_dp) .or. &
                        (abs(curr_sign(ind(1, irdm))) > 1.0e-10_dp .and. abs(iter_occ(av_ind_1)) < 1.0e-10_dp) .or. &
                        (abs(curr_sign(ind(2, irdm))) > 1.0e-10_dp .and. abs(iter_occ(av_ind_2)) < 1.0e-10_dp)) then

                        ! In this case we want to include this diagonal element,
                        ! so transfer the sign.
                        adapted_sign(av_ind_1:av_ind_2) = av_sign(av_ind_1:av_ind_2)
                    end if
                end do

                ! At least one of the signs has just gone to zero or just become
                ! reoccupied, so we need to add in diagonal elements and connections to HF
                if (any(abs(adapted_sign) > 1.e-12_dp)) then
                    if (tAllContribs .or. all_runs_are_initiator(ilut_list(:, idet))) &
                        call det_removed_fill_diag_rdm(spawn, one_rdms, ilut_list(:, idet), adapted_sign, iter_occ, tLC)
                end if

            end do

        end associate

    end subroutine fill_rdm_diag_wrapper