subroutine AddNewHashDet(TotWalkersNew, iLutCurr, DetHash, nJ, HDiag, HOffDiag, DetPosition, err)
! Add a new determinant to the main list. This involves updating the
! list length, copying it across, updating its flag, adding its diagonal
! helement (if neccessary). We also need to update the hash table to
! point at it correctly.
integer, intent(inout) :: TotWalkersNew
integer(n_int), intent(inout) :: iLutCurr(0:NIfTot)
integer, intent(in) :: DetHash, nJ(nel)
integer, intent(out) :: DetPosition
real(dp), intent(in) :: HDiag
HElement_t(dp), intent(in) :: HOffDiag
integer, intent(out) :: err
HElement_t(dp) :: trial_amps(ntrial_excits)
logical :: tTrial, tCon
real(dp), dimension(lenof_sign) :: SignCurr
err = 0
if (iStartFreeSlot <= iEndFreeSlot) then
! We can slot it into a free slot in the main list, rather than increase its length
DetPosition = FreeSlot(iStartFreeSlot)
CurrentDets(:, DetPosition) = iLutCurr(:)
iStartFreeSlot = iStartFreeSlot + 1
else
! We must increase the length of the main list to slot the new walker in
TotWalkersNew = TotWalkersNew + 1
DetPosition = TotWalkersNew
if (TotWalkersNew >= MaxWalkersPart) then
! return with an error
err = 1
return
end if
CurrentDets(:, DetPosition) = iLutCurr(:)
! if the list is almost full, activate the walker decay
if (t_prone_walkers .and. TotWalkersNew > 0.95_dp * real(MaxWalkersPart, dp)) then
t_activate_decay = .true.
write(stderr, *) "Warning: Starting to randomly kill singly-spawned walkers"
end if
end if
CurrentDets(:, DetPosition) = iLutCurr(:)
! For the RDM code we need to set all of the elements of CurrentH to 0,
! except the first one, holding the diagonal Hamiltonian element.
global_determinant_data(:, DetPosition) = 0.0_dp
call set_det_diagH(DetPosition, real(HDiag, dp) - Hii)
call set_det_offdiagH(DetPosition, HOffDiag)
! we reset the death timer, so this determinant can linger again if
! it died before
! we add the determinant to the cache
call store_decoding(DetPosition, nJ)
if (tSeniorInitiators) then
call extract_sign(ilutCurr, SignCurr)
call set_all_spawn_pops(DetPosition, SignCurr)
call reset_all_tau_ints(DetPosition)
call reset_all_shift_ints(DetPosition)
end if
if (tAutoAdaptiveShift) then
call reset_all_tot_spawns(DetPosition)
call reset_all_acc_spawns(DetPosition)
end if
! If using a trial wavefunction, search to see if this state is in
! either the trial or connected space. If so, *_search_trial returns
! the corresponding amplitude, which is stored.
!
! n.b. if this routine is called from load balancing whilst loading a popsfile, the
! trial wavefunction code will not be enabled yet (despite being enabled).
! Therefore, skip this functionality
if (tTrialWavefunction .and. allocated(current_trial_amps)) then
! Search to see if this is a trial or connected state, and
! retreive the corresponding amplitude (zero if neither a trial or
! connected state).
if (tTrialHash) then
call hash_search_trial(CurrentDets(:, DetPosition), nJ, trial_amps, tTrial, tCon)
else
call bin_search_trial(CurrentDets(:, DetPosition), trial_amps, tTrial, tCon)
end if
! Set the appropraite flag (if any). Unset flags which aren't
! appropriate, just in case.
if (tTrial) then
call set_flag(CurrentDets(:, DetPosition), flag_trial, .true.)
call set_flag(CurrentDets(:, DetPosition), flag_connected, .false.)
else if (tCon) then
call set_flag(CurrentDets(:, DetPosition), flag_trial, .false.)
call set_flag(CurrentDets(:, DetPosition), flag_connected, .true.)
else
call set_flag(CurrentDets(:, DetPosition), flag_trial, .false.)
call set_flag(CurrentDets(:, DetPosition), flag_connected, .false.)
end if
! Set the amplitude (which may be zero).
current_trial_amps(:, DetPosition) = trial_amps
end if
! If we are storing spawning rates for continuous time propagation, do
! it here
if (tContTimeFCIMC .and. tContTimeFull) then
call set_spawn_rate(DetPosition, spawn_rate_full(nJ, ilutCurr))
end if
!In case we are filling a hole, clear the removed flag
call set_flag(CurrentDets(:, DetPosition), flag_removed, .false.)
if (associated(lookup_supergroup_indexer)) then
call set_supergroup_idx(DetPosition, lookup_supergroup_indexer%idx_nI(nJ))
end if
! Add the new determinant to the hash table.
call add_hash_table_entry(HashIndex, DetPosition, DetHash)
end subroutine AddNewHashDet