assign_trial_states Subroutine

public subroutine assign_trial_states(replica_pairs, ilut_list, ilut_ht, trial_dets, trial_amps, trials_kept, energies, energies_kept)

Arguments

Type IntentOptional Attributes Name
logical, intent(in) :: replica_pairs
integer(kind=n_int), intent(in) :: ilut_list(0:,:)
type(ll_node), intent(inout), pointer :: ilut_ht(:)
integer(kind=n_int), intent(in) :: trial_dets(0:,:)
real(kind=dp), intent(in) :: trial_amps(:,:)
real(kind=dp), intent(out) :: trials_kept(:,:)
real(kind=dp), intent(in) :: energies(:)
real(kind=dp), intent(out) :: energies_kept(:)

Contents

Source Code


Source Code

    subroutine assign_trial_states(replica_pairs, ilut_list, ilut_ht, trial_dets, trial_amps, &
                                   trials_kept, energies, energies_kept)

        ! Calculate the overlaps between each trial state and FCIQMC replica
        ! pair. For each replica, keep the trial state which has the largest
        ! overlap (by magnitude).

        use FciMCData, only: ll_node
        use hash, only: hash_table_lookup

        logical, intent(in) :: replica_pairs
        integer(n_int), intent(in) :: ilut_list(0:, :)
        type(ll_node), pointer, intent(inout) :: ilut_ht(:)
        integer(n_int), intent(in) :: trial_dets(0:, :)
        HElement_t(dp), intent(in) :: trial_amps(:, :)
        HElement_t(dp), intent(out) :: trials_kept(:, :)
        real(dp), intent(in) :: energies(:)
        real(dp), intent(out) :: energies_kept(:)

        integer :: i, idet, itrial, ireplica, det_ind, hash_val
        integer :: nI(nel), best_trial(1)
        real(dp) :: fciqmc_amps_real(size(energies_kept)), all_fciqmc_amps(lenof_sign)
        real(dp) :: overlaps_real(size(energies_kept), size(trial_amps, 1))
        real(dp) :: all_overlaps_real(size(energies_kept), size(trial_amps, 1))
#ifdef CMPLX_
        real(dp) :: overlaps_imag(size(energies_kept), size(trial_amps, 1))
        real(dp) :: all_overlaps_imag(size(energies_kept), size(trial_amps, 1))
        real(dp) :: fciqmc_amps_imag(size(energies_kept))
#endif
        logical :: tDetFound

        overlaps_real = 0.0_dp
        all_overlaps_real = 0.0_dp
#ifdef CMPLX_
        overlaps_imag = 0.0_dp
        all_overlaps_imag = 0.0_dp
        unused_var(replica_pairs)
#endif

        ! Loop over all basis states (determinants) in the trial space.
        ! For each, add the overlap for each trial amplitude to a running
        ! total for replica-trial state combinations.

        do idet = 1, size(trial_amps, 2)
            ! Find if this determinant is occupied in any of the FCIQMC wave
            ! functions.
            call decode_bit_det(nI, trial_dets(0:NIfTot, idet))
            ! Search the hash table for this determinant.
            call hash_table_lookup(nI, trial_dets(:, idet), nifd, ilut_ht, ilut_list, det_ind, hash_val, tDetFound)
            if (tDetFound) then
                call extract_sign(ilut_list(:, det_ind), all_fciqmc_amps)
#ifdef CMPLX_
                do i = 1, inum_runs
                    fciqmc_amps_real(i) = all_fciqmc_amps(min_part_type(i))
                    fciqmc_amps_imag(i) = all_fciqmc_amps(max_part_type(i))
                end do
#else
                if (replica_pairs) then
#if defined(PROG_NUMRUNS_) || defined(DOUBLERUN_)
                    do i = 1, lenof_sign.div.2
                        ! When using pairs of replicas, average their amplitudes.
                        fciqmc_amps_real(i) = sum(all_fciqmc_amps(2 * i - 1:2 * i)) / 2.0_dp
                    end do
#endif
                else
                    fciqmc_amps_real = all_fciqmc_amps
                end if
#endif
                ! Add in the outer product between fciqmc_amps and the trial
                ! state amplitudes.
                do itrial = 1, size(trial_amps, 1)
#ifdef CMPLX_
                    ! (a+ib)(c+id) = ac-bd +i(ad+bc)
                    overlaps_real(:, itrial) = overlaps_real(:, itrial) &
                                    + real(trial_amps(itrial, idet)) * fciqmc_amps_real - aimag(trial_amps(itrial, idet)) * fciqmc_amps_imag
                    overlaps_imag(:, itrial) = overlaps_imag(:, itrial) &
                                    + real(trial_amps(itrial, idet)) * fciqmc_amps_imag + aimag(trial_amps(itrial, idet)) * fciqmc_amps_real
#else
                    overlaps_real(:, itrial) = overlaps_real(:, itrial) + trial_amps(itrial, idet) * fciqmc_amps_real
#endif
                end do
            end if
        end do

        call MPISumAll(overlaps_real, all_overlaps_real)
#ifdef CMPLX_
        call MPISumAll(overlaps_imag, all_overlaps_imag)
#endif

        ! Now, find the best trial state for each FCIQMC replica:
        if (t_choose_trial_state) then

#ifdef CMPLX_
            do ireplica = 1, inum_runs
                trials_kept(ireplica, :) = trial_amps(trial_excit_choice(ireplica), :)
                energies_kept(ireplica) = energies(trial_excit_choice(ireplica))
            end do
#else
            if (replica_pairs) then
                do ireplica = 1, lenof_sign.div.2
                    trials_kept(ireplica, :) = trial_amps(trial_excit_choice(ireplica), :)
                    energies_kept(ireplica) = energies(trial_excit_choice(ireplica))

                    root_print "trial state: ", trial_excit_choice(ireplica), &
                        " chosen for replica: ", ireplica, &
                        " chosen by input, with energy: ", energies(trial_excit_choice(ireplica))
                end do
            else
                do ireplica = 1, lenof_sign
                    trials_kept(ireplica, :) = trial_amps(trial_excit_choice(ireplica), :)
                    energies_kept(ireplica) = energies(trial_excit_choice(ireplica))

                    root_print "trial state: ", trial_excit_choice(ireplica), &
                        " chosen for replica: ", ireplica, &
                        " chosen by input, with energy: ", energies(trial_excit_choice(ireplica))

                end do
            end if
#endif
        else
#ifdef CMPLX_
            do ireplica = 1, inum_runs
                best_trial = maxloc(abs(all_overlaps_real(ireplica, :)**2 + all_overlaps_imag(ireplica, :)**2))
                trials_kept(ireplica, :) = trial_amps(best_trial(1), :)
                energies_kept(ireplica) = energies(best_trial(1))
            end do
#else
            if (replica_pairs) then
                do ireplica = 1, lenof_sign.div.2
                    best_trial = maxloc(abs(all_overlaps_real(ireplica, :)))
                    trials_kept(ireplica, :) = trial_amps(best_trial(1), :)
                    energies_kept(ireplica) = energies(best_trial(1))
#ifdef DEBUG_
                    root_print "trial state: ", best_trial, " kept for replica ", ireplica, &
                        " based on overlap, with energy: ", energies(best_trial(1))
#endif
                end do
            else
                do ireplica = 1, lenof_sign
                    best_trial = maxloc(abs(all_overlaps_real(ireplica, :)))
                    trials_kept(ireplica, :) = trial_amps(best_trial(1), :)
                    energies_kept(ireplica) = energies(best_trial(1))

                    root_print "trial state: ", best_trial, " kept for replica ", ireplica, &
                        " based on overlap, with energy: ", energies(best_trial(1))

                end do
            end if
#endif
        end if

    end subroutine assign_trial_states