subroutine generate_init_config_this_proc(nwalkers, nwalkers_per_site_init, tOccupyDetermSpace)
! This routine will distribute nwalkers walkers uniformly across all possible determinants.
use CalcData, only: tSemiStochastic
use dSFMT_interface, only: genrand_real2_dSFMT
use FciMCData, only: HashIndex, TotWalkers, CurrentDets, HFDet
use FciMCData, only: TotParts, TotPartsOld, AllTotParts, AllTotPartsOld, ilutHF
use core_space_util, only: cs_replicas
use load_balance_calcnodes, only: DetermineDetNode
use hash, only: rm_unocc_dets_from_hash_table, hash_table_lookup
use hash, only: add_hash_table_entry
use hilbert_space_size, only: CreateRandomExcitLevDetUnbias, create_rand_heisenberg_det
use hilbert_space_size, only: create_rand_det_no_sym
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
logical, intent(in) :: tOccupyDetermSpace
integer :: proc, excit, nattempts, hash_val, det_ind, nI(nel)
integer :: ideterm, ndets
integer(n_int) :: ilut(0:NIfTot), int_sign(lenof_sign_kp)
real(dp) :: real_sign_1(lenof_sign_kp), real_sign_2(lenof_sign_kp)
real(dp) :: new_sign(lenof_sign_kp)
real(dp) :: r
logical :: tDetermAllOccupied, tDetFound
ilut = 0_n_int
ndets = 0
TotParts = 0.0_dp
ideterm = 0
tDetermAllOccupied = .false.
associate(rep => cs_replicas(1))
do
! If using the tOccupyDetermSpace option then we want to put walkers
! on states in the deterministic space first.
! If not, or if we have finished doing this, generate determinants
! randomly and uniformly.
if (tOccupyDetermSpace .and. (.not. tDetermAllOccupied)) then
ideterm = ideterm + 1
ilut = rep%core_space(:, ideterm + rep%determ_displs(iProcIndex))
call decode_bit_det(nI, ilut)
! If we have now occupied all deterministic states.
if (ideterm == rep%determ_sizes(iProcIndex)) tDetermAllOccupied = .true.
else
! Generate a random determinant (returned in 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
! If it doesn't belong to this processor, pick another.
call decode_bit_det(nI, ilut)
proc = DetermineDetNode(nel, nI, 0)
if (proc /= iProcIndex) cycle
end if
! Choose whether the walker should have a positive or negative amplitude, with
! 50% chance of each.
real_sign_1 = nwalkers_per_site_init
r = genrand_real2_dSFMT()
if (r < 0.5_dp) real_sign_1 = -1.0_dp * real_sign_1
int_sign = transfer(real_sign_1, int_sign)
tDetFound = .false.
! Search the hash table to see if this determinant is in CurrentDets
! already.
call hash_table_lookup(nI, ilut, nifd, HashIndex, CurrentDets, det_ind, hash_val, tDetFound)
if (tDetFound) then
! This determinant is already in CurrentDets. Just need to add
! the sign of this new walker on and update stats.
int_sign = CurrentDets(IlutBits%ind_pop:IlutBits%ind_pop + lenof_sign_kp - 1, det_ind)
real_sign_2 = transfer(int_sign, real_sign_2)
new_sign = real_sign_1 + real_sign_2
call encode_sign(CurrentDets(:, det_ind), new_sign)
TotParts = TotParts - abs(real_sign_2) + abs(new_sign)
else
! A new determiant needs to be added.
ndets = ndets + 1
det_ind = ndets
! Add the new determinant to the hash table.
call add_hash_table_entry(HashIndex, det_ind, hash_val)
! Copy determinant data across.
CurrentDets(0:nifd, det_ind) = ilut(0:nifd)
CurrentDets(IlutBits%ind_pop:IlutBits%ind_pop + lenof_sign_kp - 1, det_ind) = int_sign
CurrentDets(IlutBits%ind_flag, det_ind) = 0_n_int
TotParts = TotParts + abs(real_sign_1)
end if
! If we've reached the total requested number of walkers, finish.
if (TotParts(1) >= nwalkers) exit
end do
call MPISum(TotParts, AllTotParts)
TotPartsOld = TotParts
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(rep)
call add_core_states_currentdet_hash(core_run)
else
! Some determinants may have become occupied and then unoccupied in
! the course of the above. We need to remove the entries for these
! determinants from the hash table. For semi-stochastic calculations
! this is done in add_core_states_currentdet_hash.
call rm_unocc_dets_from_hash_table(HashIndex, CurrentDets, ndets)
end if
end associate
end subroutine generate_init_config_this_proc