generate_init_config_basic Subroutine

private subroutine generate_init_config_basic(nwalkers, nwalkers_per_site_init)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nwalkers
real(kind=dp) :: nwalkers_per_site_init

Contents


Source Code

    subroutine generate_init_config_basic(nwalkers, nwalkers_per_site_init)

        ! This routine will distribute nwalkers walkers uniformly across all possible determinants.

        use AnnihilationMod, only: SendProcNewParts, CompressSpawnedList
        use CalcData, only: tTruncInitiator, tSemiStochastic
        use dSFMT_interface, only: genrand_real2_dSFMT
        use fcimc_helper, only: create_particle
        use FciMCData, only: nWalkerHashes, fcimc_iter_data, HashIndex, SpawnedParts, SpawnedParts2
        use FciMCData, only: TotWalkers, CurrentDets, InitialSpawnedSlots, TotParts, TotPartsOld
        use FciMCData, only: AllTotParts, AllTotPartsOld, ValidSpawnedList, ilutHF, HFDet
        use hash, only: fill_in_hash_table
        use hilbert_space_size, only: CreateRandomExcitLevDetUnbias, create_rand_heisenberg_det
        use hilbert_space_size, only: create_rand_det_no_sym
        use replica_data, only: allocate_iter_data, clean_iter_data
        use semi_stoch_procs, only: copy_core_dets_to_spawnedparts, add_core_states_currentdet_hash
        use SystemData, only: nel, tHeisenberg, tAllSymSectors

        integer, intent(in) :: nwalkers
        real(dp) :: nwalkers_per_site_init
        integer :: i, ireplica, excit, nattempts, DetHash
        integer :: nspawns, ndets, err
        integer(n_int) :: ilut(0:NIfTot)
        integer :: nI(nel)
        real(dp) :: r, walker_amp, walker_sign(lenof_sign_kp)
        logical :: tInitiatorTemp
        type(fcimc_iter_data) :: unused_data
        integer(n_int), pointer :: PointTemp(:, :)
        HElement_t(dp) :: hdiag_spawn

        call allocate_iter_data(unused_data)

        ! Turn off the initiator method for the annihilation steps to be used here.
        tInitiatorTemp = tTruncInitiator
        tTruncInitiator = .false.

        ! Set the spawning slots to their starting positions.
        ValidSpawnedList = InitialSpawnedSlots

        ilut = 0_n_int
        nspawns = ceiling(real(nwalkers, dp) / nwalkers_per_site_init)

        do i = 1, nspawns
            ! Generate a determinant (output to ilut).
            if (tAllSymSectors) then
                call create_rand_det_no_sym(ilut)
            else
                if (tHeisenberg) then
                    call create_rand_heisenberg_det(ilut)
                else
                    call CreateRandomExcitLevDetUnbias(nel, HFDet, ilutHF, ilut, excit, nattempts)
                end if
            end if
            call decode_bit_det(nI, ilut)

            ! Choose whether the walker should have a positive or negative amplitude, with
            ! 50% chance of each.
            walker_amp = nwalkers_per_site_init
            r = genrand_real2_dSFMT()
            if (r < 0.5_dp) walker_amp = -1.0_dp * walker_amp

            do ireplica = 1, inum_runs
                walker_sign = 0.0_dp
                walker_sign(ireplica) = walker_amp
                call create_particle(nI, ilut, walker_sign, ireplica, hdiag_spawn, err)
            end do
        end do

        ! Perform annihilation steps:
        ! Send the walkers to their correct processors. The resulting walkers will be stored in
        ! SpawnedParts2.
        call SendProcNewParts(ndets, tSingleProc=.false.)
        ! CompressSpawnedList works on SpawnedParts, not SpawnedParts2, so swap the pointers around.
        PointTemp => SpawnedParts2
        SpawnedParts2 => SpawnedParts
        SpawnedParts => PointTemp
        call CompressSpawnedList(ndets, unused_data)

        ! Finally, add the determinants in the spawned walker list to the main walker list.
        ! Copy the determinants themselves to CurrentDets.
        TotParts = 0.0_dp
        do i = 1, ndets
            CurrentDets(:, i) = SpawnedParts(0:NIfTot, i)
            walker_sign = transfer(CurrentDets(IlutBits%ind_pop:IlutBits%ind_pop + lenof_sign_kp - 1, i), walker_sign)
            TotParts = TotParts + abs(walker_sign)
        end do
        TotPartsOld = TotParts

        ! Add the entries into the hash table.
        call fill_in_hash_table(HashIndex, nWalkerHashes, CurrentDets, ndets, .true.)

        call MPISum(TotParts, AllTotParts)
        AllTotPartsOld = AllTotParts
        TotWalkers = int(ndets, int64)

        if (tSemiStochastic) then
            ! Always need the core determinants to be at the top of CurrentDets, even when unoccupied.
            ! These routines will do this.
            call copy_core_dets_to_spawnedparts(cs_replicas(core_run))
            call add_core_states_currentdet_hash(core_run)
        end if

        ValidSpawnedList = InitialSpawnedSlots

        ! Turn the initiator method back on, if it was turned off at the start of this routine.
        tTruncInitiator = tInitiatorTemp

        call clean_iter_data(unused_data)

    end subroutine generate_init_config_basic