AnnihilateDiagParts Subroutine

public subroutine AnnihilateDiagParts(ValidSpawned, TotWalkersNew, iter_data)

Arguments

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

Contents

Source Code


Source Code

    subroutine AnnihilateDiagParts(ValidSpawned, TotWalkersNew, iter_data)
        ! this is the new "annihilation" routine which mimics the actual
        ! diagonal death-step, between the reloaded CurrentDets y(n) and the
        ! diagonal parts in k2, DiagParts
        integer, intent(inout) :: ValidSpawned, TotWalkersNew
        type(fcimc_iter_data), intent(inout) :: iter_data
        character(*), parameter :: this_routine = "AnnihilateDiagParts"

        integer :: PartInd, i, j
        real(dp), dimension(lenof_sign) :: CurrentSign, SpawnedSign, SignTemp
        real(dp), dimension(lenof_sign) :: SignProd
        real(dp) :: HDiag
        HElement_t(dp) :: HOffDiag
        integer :: DetHash, nJ(nel)
        logical :: tSuccess, tDetermState
        integer :: run, err, pos

        ! rewrite the original Annihilation routine to fit the new
        ! requirements here

        ! this routine updated the NoDied variables etc.. for the 2nd RK step
        ! so i dont think i need to change much here
        ! Only node roots to do this.
        if (.not. bNodeRoot) return

        do i = 1, ValidSpawned

            call decode_bit_det(nJ, DiagParts(:, i))

            ! Search the hash table HashIndex for the determinant defined by
            ! nJ and DiagParts(:,i). If it is found, tSuccess will be
            ! returned .true. and PartInd will hold the position of the
            ! determinant in CurrentDets. Else, tSuccess will be returned
            ! .false. (and PartInd shouldn't be accessed).
            ! Also, the hash value, DetHash, is returned by this routine.
            ! tSuccess will determine whether the particle has been found or not.
            call hash_table_lookup(nJ, DiagParts(:, i), nifd, HashIndex, &
                                   CurrentDets, PartInd, DetHash, tSuccess)

            tDetermState = .false.
            CurrentSign = 0.0_dp

!            write(stdout,*) 'i,DiagParts(:,i)',i,DiagParts(:,i)

            if (tSuccess) then

                ! Our DiagParts determinant is found in CurrentDets.

                call extract_sign(CurrentDets(:, PartInd), CurrentSign)
                call extract_sign(DiagParts(:, i), SpawnedSign)

                SignProd = CurrentSign * SpawnedSign

                tDetermState = test_flag_multi(CurrentDets(:, PartInd), flag_deterministic)

                if (sum(abs(CurrentSign)) >= 1.e-12_dp .or. tDetermState) then
                    ! Transfer new sign across.
                    call encode_sign(CurrentDets(:, PartInd), SpawnedSign + CurrentSign)
                    ! I dont see why DiagParts has to be nullified. In the verlet scheme
                    ! warmup, we might want to use DiagParts afterwards
                    call encode_sign(DiagParts(:, i), null_part)
                    do j = 1, lenof_sign
                        run = part_type_to_run(j)
                        if (is_run_unnocc(CurrentSign, run)) then
                            ! This determinant is actually *unoccupied* for the
                            ! walker type/set we're considering. We need to
                            ! decide whether to abort it or not.
                            ! rt_iter_adapt : also count those spawned onto an initiator
                            ! this is not necessary in the normal version as walkers are
                            ! counted there on spawn
                            ! also, they are born for any sign
                            iter_data%nborn(j) = iter_data%nborn(j) + &
                                                 abs(SpawnedSign(j))
                            NoBorn(run) = NoBorn(run) + abs(SpawnedSign(j))

                        else if (SignProd(j) < 0) then
                            ! in the real-time for the final combination
                            ! y(n) + k2 i have to check if the "spawned"
                            ! particle is actually a diagonal death/born
                            ! walker
                            ! remember the - sign when filling up the DiagParts
                            ! list -> so opposite sign here means a death!
                            iter_data%ndied(j) = iter_data%ndied(j) + &
                                                 min(abs(CurrentSign(j)), abs(SpawnedSign(j)))

                            NoDied(run) = NoDied(run) + min(abs(CurrentSign(j)), &
                                                            abs(SpawnedSign(j)))
                            ! and if the Spawned sign magnitude is even higher
                            ! then the currentsign -> born anti-particles

                            iter_data%nborn(j) = iter_data%nborn(j) + &
                                                 max(abs(SpawnedSign(j)) - abs(CurrentSign(j)), 0.0_dp)

                            NoBorn(run) = NoBorn(run) + max(abs(SpawnedSign(j)) - &
                                                            abs(CurrentSign(j)), 0.0_dp)
                        else
                            ! if it has the same sign i have to keep track of
                            ! the born particles, or as in the original
                            ! walker_death, reduce the number of died parts
                            iter_data%ndied(j) = iter_data%ndied(j) - &
                                                 abs(SpawnedSign(j))

                        end if

                    end do ! Over all components of the sign.

                    if (.not. tDetermState) then
                        call extract_sign(CurrentDets(:, PartInd), SignTemp)
                        if (IsUnoccDet(SignTemp)) then
                            ! All walkers in this main list have been annihilated
                            ! away. Remove it from the hash index array so that
                            ! no others find it (it is impossible to have another
                            ! spawned walker yet to find this determinant).
                            call remove_hash_table_entry(HashIndex, nJ, PartInd)
                            ! Add to "freeslot" list so it can be filled in.
                            iEndFreeSlot = iEndFreeSlot + 1
                            FreeSlot(iEndFreeSlot) = PartInd
                        end if
                    end if

                    ! RT_M_Merge: Disabled rdm functionality
                    ! We must use the instantaneous value for the off-diagonal
                    ! contribution. However, we can't just use CurrentSign from
                    ! the previous iteration, as this has been subject to death
                    ! but not the new walkers. We must add on SpawnedSign, so
                    ! we're effectively taking the instantaneous value from the
                    ! next iter. This is fine as it's from the other population,
                    ! and the Di and Dj signs are already strictly uncorrelated.

                end if

            end if

            if (((.not. tSuccess) .or. (tSuccess .and. sum(abs(CurrentSign)) < 1.e-12_dp .and. (.not. tDetermState)))) then

                ! Running the full, non-initiator scheme.
                ! Determinant in newly spawned list is not found in
                ! CurrentDets. If coeff <1, apply removal criterion.
                call extract_sign(DiagParts(:, i), SignTemp)

                if (.not. IsUnoccDet(SignTemp)) then
                    ! Walkers have not been aborted and so we should copy the
                    ! determinant straight over to the main list. We do not
                    ! need to recompute the hash, since this should be the
                    ! same one as was generated at the beginning of the loop.
                    ! also here treat those new walkers as born particles

                    iter_data%nborn = iter_data%nborn + abs(SignTemp)
                    do run = 1, inum_runs
                        NoBorn(run) = NoBorn(run) + sum(abs(SignTemp( &
                                                            min_part_type(run):max_part_type(run))))
                    end do
                    HDiag = get_diagonal_matel(nJ, DiagParts(:, i))
                    HOffDiag = get_off_diagonal_matel(nJ, DiagParts(:, i))
                    call AddNewHashDet(TotWalkersNew, DiagParts(:, i), DetHash, nJ, HDiag, HOffDiag, pos, err)
                    if (err /= 0) exit

                end if
            end if
        end do
        ! Update remaining number of holes in list for walkers stats.
        if (iStartFreeSlot > iEndFreeSlot) then
            ! All slots filled
            HolesInList = 0
        else
            HolesInList = iEndFreeSlot - (iStartFreeSlot - 1)
        end if

    end subroutine AnnihilateDiagParts