FindResidualParticle_simple Subroutine

private subroutine FindResidualParticle_simple(cum_det, new_det, part_type, iter_data)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(inout) :: cum_det(0:nIfTot)
integer(kind=n_int), intent(in) :: new_det(0:niftot+nifd+2)
integer, intent(in) :: part_type
type(fcimc_iter_data), intent(inout) :: iter_data

Contents


Source Code

    subroutine FindResidualParticle_simple(cum_det, new_det, part_type, iter_data)

        ! This routine is called whilst compressing the spawned list during
        ! annihilation. It considers the sign and flags from two particles
        ! on the same determinant, and calculates the correct sign/flags
        ! for the compressed particle.
        !
        ! --> The information is stored within the first particle in a block
        ! --> Should be called for real/imaginary particles seperately

        integer(n_int), intent(inout) :: cum_det(0:nIfTot)
        integer(n_int), intent(in) :: new_det(0:niftot + nifd + 2)
        integer, intent(in) :: part_type
        !integer, intent(inout) :: Parent_Array_Ind
        type(fcimc_iter_data), intent(inout) :: iter_data

        real(dp) :: new_sgn, cum_sgn, updated_sign
        real(dp) :: sgn_prod
        integer :: run

        new_sgn = extract_part_sign(new_det, part_type)
        cum_sgn = extract_part_sign(cum_det, part_type)

        ! If the cumulative and new signs for this replica are both non-zero
        ! then there have been at least two spawning events to this site, so
        ! set the initiator flag.
        if (tTruncInitiator) then
            if (test_flag(new_det, get_initiator_flag(part_type))) &
                call set_flag(cum_det, get_initiator_flag(part_type))
        end if

        sgn_prod = cum_sgn * new_sgn

        ! Update annihilation statistics.
        if (.not. tPrecond) then
            if (sgn_prod < 0.0_dp) then
                run = part_type_to_run(part_type)
                Annihilated(run) = Annihilated(run) + 2 * min(abs(cum_sgn), abs(new_sgn))
                iter_data%nannihil(part_type) = iter_data%nannihil(part_type) &
                                                + 2 * min(abs(cum_sgn), abs(new_sgn))
            end if
        end if

        ! Update the cumulative sign count.
        updated_sign = cum_sgn + new_sgn
        call encode_part_sign(cum_det, updated_sign, part_type)

    end subroutine FindResidualParticle_simple