apply_hamiltonian Subroutine

public subroutine apply_hamiltonian(population, popsize, tGetFreeSlots, tGetInitflags, tSumE)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(in) :: population(0:,1:)
integer, intent(in) :: popsize
logical, intent(in) :: tGetFreeSlots
logical, intent(in) :: tGetInitflags
logical, intent(in) :: tSumE

Contents

Source Code


Source Code

    subroutine apply_hamiltonian(population, popsize, tGetFreeSlots, tGetInitFlags, tSumE)
        ! this subroutine performs spawning from population to spawnVec by applying
        ! delta t H once. No annihilation is performed and no other steps performed
        ! keep it minimalistic and stick to the SRP principle

        ! by default writes into SpawnedParts
        integer(n_int), intent(in) :: population(0:, 1:)
        integer, intent(in) :: popsize
        ! store the free slots in population
        logical, intent(in) :: tGetInitflags, tGetFreeSlots, tSumE
        integer :: idet, nI(nel), determ_index, unused_flags, ex_level
        real(dp) :: parent_sign(lenof_sign), hdiag
        logical :: tEmptyDet, tCoreDet
        HElement_t(dp) :: hoffdiag

        ! where to put this?
        ! attempt_create => attempt_create_normal

        ValidSpawnedList = InitialSpawnedSlots
        ! we do not know what is in spawn_ht, but we clearly want it to be empty now
        call clear_hash_table(spawn_ht)

        determ_index = 1

        do idet = 1, popsize
            ! apply spawn and death for each walker
            fcimc_excit_gen_store%tFilled = .false.
            unused_flags = 0

            call extract_bit_rep(population(:, idet), nI, parent_sign, unused_flags, idet, &
                                 fcimc_excit_gen_store)

            tEmptyDet = IsUnoccDet(parent_sign)
            ! we collect free slots only in the first application
            if (tEmptyDet) then
                if (tGetFreeSlots) then
                    iEndFreeSlot = iEndFreeSlot + 1
                    FreeSlot(iEndFreeSlot) = idet
                end if
                cycle
            end if

            ! also, initiator flags are only reset in the first iteration
            ! when using the initiator criterium
            ! else, we just let the flags be
            hdiag = 0.0_dp
            if (tGetInitFlags) call CalcParentFlag(idet, unused_flags)
            hdiag = get_diagonal_matel(nI, population(:, idet)) - Hii
            ! if desired, sum in the energy contribution
            if (tSumE) then
                ex_level = FindBitExcitLevel(ilutRef, population(:, idet), max_calc_ex_level)
                hoffdiag = get_off_diagonal_matel(nI, population(:, idet))
                call SumEContrib(nI, ex_level, parent_sign, population(:, idet), &
                                 hdiag, hoffdiag, 1.0_dp, tPairedReplicas, idet)
            end if

            tCoreDet = check_determ_flag(population(:, idet))

            if (tCoreDet) then
                associate(rep => cs_replicas(core_run))
                    rep%indices_of_determ_states(determ_index) = idet
                    rep%partial_determ_vecs(:, determ_index) = parent_sign
                end associate
                determ_index = determ_index + 1
                if (IsUnoccDet(parent_sign)) cycle
            end if

            ! the initiator flags are set upon the first iteration

            call perform_spawn(population(:, idet), nI, parent_sign, hdiag, tCoreDet)
        end do

        if (tSemiStochastic) call real_time_determ_projection()
    end subroutine apply_hamiltonian