get_qmc_trial_weights Subroutine

public subroutine get_qmc_trial_weights(trial_iluts, ntrial, qmc_iluts, qmc_ht, nexcit, paired_replicas, evecs_this_proc)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(in) :: trial_iluts(0:,:)
integer, intent(in) :: ntrial
integer(kind=n_int), intent(in) :: qmc_iluts(0:,:)
type(ll_node), intent(inout), pointer :: qmc_ht(:)
integer, intent(in) :: nexcit
logical, intent(in) :: paired_replicas
real(kind=dp), intent(out) :: evecs_this_proc(:,:)

Contents

Source Code


Source Code

    subroutine get_qmc_trial_weights(trial_iluts, ntrial, qmc_iluts, qmc_ht, nexcit, paired_replicas, evecs_this_proc)

        use bit_reps, only: decode_bit_det
        use hash, only: hash_table_lookup
        use Parallel_neci, only: MPISumAll
        use SystemData, only: nel

        integer(n_int), intent(in) :: trial_iluts(0:, :)
        integer, intent(in) :: ntrial
        integer(n_int), intent(in) :: qmc_iluts(0:, :)
        type(ll_node), pointer, intent(inout) :: qmc_ht(:)
        integer, intent(in) :: nexcit
        logical, intent(in) :: paired_replicas
        real(dp), intent(out) :: evecs_this_proc(:, :)

        integer :: i, j, ind, hash_val
        integer :: nI(nel)
        real(dp) :: qmc_sign(lenof_sign), trial_sign(nexcit)
        real(dp) :: norm(nexcit), tot_norm(nexcit)
        logical :: found
        character(*), parameter :: this_routine = 'get_qmc_trial_weights'

        if (paired_replicas) then
            ASSERT(nexcit == (lenof_sign.div.2))
        else
            ASSERT(nexcit == lenof_sign)
        end if

        norm = 0.0_dp

        ! Loop over all trial states.
        do i = 1, ntrial
            ! Check the QMC hash table to see if this state exists in the
            ! QMC list.
            call decode_bit_det(nI, trial_iluts(:, i))
            call hash_table_lookup(nI, trial_iluts(:, i), nifd, qmc_ht, &
                                   qmc_iluts, ind, hash_val, found)
            if (found) then
                call extract_sign(qmc_iluts(:, ind), qmc_sign)

                ! If using paired replicas then average the sign on each pair.
                if (paired_replicas) then
                    do j = 2, lenof_sign, 2
                        trial_sign(j / 2) = sum(qmc_sign(j - 1:j) / 2.0_dp)
                    end do
                else
                    trial_sign = qmc_sign
                end if

                norm = norm + trial_sign**2
                evecs_this_proc(:, i) = trial_sign
            else
                evecs_this_proc(:, i) = 0.0_dp
            end if
        end do

        call MPISumAll(norm, tot_norm)

        do i = 1, ntrial
            evecs_this_proc(:, i) = evecs_this_proc(:, i) / sqrt(tot_norm)
        end do

    end subroutine get_qmc_trial_weights