CompressSpawnedList_simple Subroutine

private subroutine CompressSpawnedList_simple(ValidSpawned, iter_data)

Arguments

Type IntentOptional Attributes Name
integer :: ValidSpawned
type(fcimc_iter_data), intent(inout) :: iter_data

Contents


Source Code

    subroutine CompressSpawnedList_simple(ValidSpawned, iter_data)

        ! This sorts and compresses the spawned list to make it easier for the
        ! rest of the annihilation process. This is not essential, but should
        ! prove worthwhile.

        type(fcimc_iter_data), intent(inout) :: iter_data
        integer :: VecInd, ValidSpawned, DetsMerged, i, BeginningBlockDet, FirstInitIndex, CurrentBlockDet
        real(dp) :: SpawnedSign(lenof_sign), temp_sign(lenof_sign), temp_sign_2(lenof_sign)
        integer :: EndBlockDet, part_type, Parent_Array_Ind
        integer(n_int), pointer :: PointTemp(:, :)
        integer(n_int) :: cum_det(0:IlutBits%len_bcast), cum_det_cancel(0:IlutBits%len_bcast)
        logical :: any_allow, any_cancel
        character(len=*), parameter :: t_r = 'CompressSpawnedList_simple'
        type(timer), save :: Sort_time

        if (.not. bNodeRoot) return

        Sort_time%timer_name = 'Compress Sort interface'
        call set_timer(Sort_time, 20)

        call sort(SpawnedParts(0:IlutBits%len_bcast, 1:ValidSpawned), ilut_lt, ilut_gt)

        call halt_timer(Sort_time)

        if (tHistSpawn) HistMinInd2(1:NEl) = FCIDetIndex(1:NEl)

        !write(stdout,*) "SpawnedParts before:"
        !do i = 1, ValidSpawned
        !    call extract_sign (SpawnedParts(:, i), temp_sign)
        !    write(stdout,'(i7, 4x, i16, 4x, f18.7, 4x, l1)') i, SpawnedParts(0,i), temp_sign, &
        !        test_flag(SpawnedParts(:,i), get_initiator_flag(1))
        !end do

        ! First, we compress the list of spawned particles, so that they are
        ! only specified at most once in each processors list. During this, we
        ! transfer the particles from SpawnedParts to SpawnedParts2. If we are
        ! working with complex walkers, we essentially do the same thing twice,
        ! annihilating real and imaginary particles seperately.

        ! This is the index in the SpawnedParts2 array to copy the compressed
        ! walkers into.
        VecInd = 1
        ! BeginningBlockDet will indicate the index of the first entry for a
        ! given determinant in SpawnedParts.
        BeginningBlockDet = 1
        DetsMerged = 0
        Parent_Array_Ind = 1
        Spawned_Parts_Zero = 0

        do while (BeginningBlockDet <= ValidSpawned)

            ! Loop in blocks of the same determinant to the end of the list of
            ! walkers.

            FirstInitIndex = 0
            CurrentBlockDet = BeginningBlockDet + 1

            do while (CurrentBlockDet <= ValidSpawned)
                if (.not. (DetBitEQ(SpawnedParts(:, BeginningBlockDet), &
                                    SpawnedParts(:, CurrentBlockDet)))) exit
                ! Loop over walkers on the same determinant in SpawnedParts.
                CurrentBlockDet = CurrentBlockDet + 1
            end do

            ! EndBlockDet indicates that we have reached the end of the block
            ! of similar dets.
            EndBlockDet = CurrentBlockDet - 1

            if (EndBlockDet == BeginningBlockDet) then
                ! Optimisation: This block only consists of one entry. Simply
                !               copy it across rather than explicitly searching
                !               the list.

                ! If this one entry has no amplitude then don't add it to the
                ! compressed list, but just cycle.
                call extract_sign(SpawnedParts(:, BeginningBlockDet), temp_sign)
                if ((sum(abs(temp_sign)) < 1.e-12_dp)) then
                    DetsMerged = DetsMerged + 1
                    BeginningBlockDet = CurrentBlockDet
                    cycle
                end if

                ! Transfer all info to the other array.
                SpawnedParts2(:, VecInd) = SpawnedParts(:, BeginningBlockDet)

                VecInd = VecInd + 1
                ! Move onto the next block of determinants.
                BeginningBlockDet = CurrentBlockDet
                cycle ! Skip the rest of this block.
            end if

            ! Reset the cumulative determinant
            cum_det = 0_n_int
            cum_det(0:nifd) = SpawnedParts(0:nifd, BeginningBlockDet)

            ! This will only get used with the pure-initiator-space option.
            ! In this case, some spawnings to a site can be accepted, while
            ! others are rejected. The following is used to hold the rejected
            ! ones which can be used for constructing estimates later.
            cum_det_cancel = 0_n_int
            cum_det_cancel(0:nifd) = SpawnedParts(0:nifd, BeginningBlockDet)

            if (tPreCond .or. tReplicaEstimates) then
                cum_det(IlutBits%ind_hdiag) &
                    = SpawnedParts(IlutBits%ind_hdiag, BeginningBlockDet)
                cum_det_cancel(IlutBits%ind_hdiag) &
                    = SpawnedParts(IlutBits%ind_hdiag, BeginningBlockDet)
            end if

            do i = BeginningBlockDet, EndBlockDet
                ! if logged, accumulate the number of spawn events
                if (tLogNumSpawns) then
                    cum_det(IlutBits%ind_spawn) = cum_det(IlutBits%ind_spawn) + &
                                                  SpawnedParts(IlutBits%ind_spawn, i)
                end if
                ! Annihilate in this block seperately for walkers of different types.
                do part_type = 1, lenof_sign
                    if (tHistSpawn) then
                        call extract_sign(SpawnedParts(:, i), SpawnedSign)
                        call extract_sign(SpawnedParts(:, BeginningBlockDet), temp_sign)
                        call HistAnnihilEvent(SpawnedParts, SpawnedSign, temp_sign, part_type)
                    end if

                    ! If using a pure initiator space, then some spawnings to
                    ! a site can be accepted while others are rejected, unlike
                    ! the normal initiator approach. Here, check if this
                    ! particular spawning needs to be rejected.
                    if (test_flag(SpawnedParts(:, i), get_initiator_flag(part_type))) then
                        call FindResidualParticle_simple(cum_det, &
                                                         SpawnedParts(:, i), part_type, iter_data)
                    else
                        call FindResidualParticle_simple(cum_det_cancel, &
                                                         SpawnedParts(:, i), part_type, iter_data)
                    end if
                end do ! Over all spawns to the same determinant

            end do ! Loop over particle type.

            ! Copy details into the final array.
            call extract_sign(cum_det, temp_sign)
            call extract_sign(cum_det_cancel, temp_sign_2)

            any_allow = sum(abs(temp_sign)) > 1.e-12_dp
            any_cancel = sum(abs(temp_sign_2)) > 1.e-12_dp

            if (any_allow .and. any_cancel) then
                SpawnedParts2(0:NIfTot, VecInd) = cum_det(0:NIfTot)
                SpawnedParts2(0:NIfTot, VecInd + 1) = cum_det_cancel(0:NIfTot)
                if (tPreCond .or. tReplicaEstimates) then
                    SpawnedParts2(IlutBits%ind_hdiag, VecInd) &
                        = cum_det(IlutBits%ind_hdiag)
                    SpawnedParts2(IlutBits%ind_hdiag, VecInd + 1) &
                        = cum_det_cancel(IlutBits%ind_hdiag)
                end if

                VecInd = VecInd + 2
                DetsMerged = DetsMerged + EndBlockDet - BeginningBlockDet - 1
            else if (any_allow .and. (.not. any_cancel)) then
                SpawnedParts2(0:NIfTot, VecInd) = cum_det(0:NIfTot)
                if (tPreCond .or. tReplicaEstimates) then
                    SpawnedParts2(IlutBits%ind_hdiag, VecInd) &
                        = cum_det(IlutBits%ind_hdiag)
                end if

                VecInd = VecInd + 1
                DetsMerged = DetsMerged + EndBlockDet - BeginningBlockDet
            else if ((.not. any_allow) .and. any_cancel) then
                SpawnedParts2(0:NIfTot, VecInd) = cum_det_cancel(0:NIfTot)
                if (tPreCond .or. tReplicaEstimates) then
                    SpawnedParts2(IlutBits%ind_hdiag, VecInd) &
                        = cum_det_cancel(IlutBits%ind_hdiag)
                end if

                VecInd = VecInd + 1
                DetsMerged = DetsMerged + EndBlockDet - BeginningBlockDet
            else
                ! All particles from block have been annihilated.
                DetsMerged = DetsMerged + EndBlockDet - BeginningBlockDet + 1
            end if

            ! Move onto the next block of determinants.
            BeginningBlockDet = CurrentBlockDet

        end do

        ! This is the new number of unique spawned determinants on the processor.
        ValidSpawned = ValidSpawned - DetsMerged
        if (ValidSpawned /= (VecInd - 1)) then
            call stop_all(t_r, "Error in compression of spawned list")
        end if

        ! Want the compressed list in spawnedparts at the end of it - swap
        ! pointers around.
        PointTemp => SpawnedParts2
        SpawnedParts2 => SpawnedParts
        SpawnedParts => PointTemp

        !write(stdout,*) "SpawnedParts after:"
        !do i = 1, ValidSpawned
        !    if (SpawnedParts(0,i) == SpawnedParts(0,i+1)) write(stdout,*) "Here!"
        !    call extract_sign (SpawnedParts(:, i), temp_sign)
        !    write(stdout,'(i7, 4x, i16, 4x, f18.7, 4x, l1)') i, SpawnedParts(0,i), temp_sign, &
        !        test_flag(SpawnedParts(:,i), get_initiator_flag(1))
        !end do

    end subroutine CompressSpawnedList_simple