add_core_states_currentdet_hash Subroutine

public subroutine add_core_states_currentdet_hash(run)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: run

Contents


Source Code

    subroutine add_core_states_currentdet_hash(run)

        ! This routine adds the core states in SpawnedParts into CurrentDets. For all
        ! such states already in CurrentDets, we want to keep the amplitude (which
        ! may have come from a popsfile).

        ! This routine is for when the hashed walker main list. In this case,
        ! as all core states are always kept in the list, it is beneficial to keep
        ! them at the top always. So, in this routine, we move the non-core states
        ! in CurrentDets to the end and add the new core states in the gaps.

        ! WARNING: If there are any determinants in CurrentDets on input which are
        ! unoccupied then, for this function to work correctly, the determiant
        ! *must* have an entry in the hash table. Otherwise, these determinants
        ! will end up being repeated in CurrentDets is they are core determinants.
        ! This isn't ideal because when the FCIQMC calculation starts, such
        ! unoccupied determinants should *not* be in the hash table. During this
        ! routine, such determinants will be removed from the hash table and so
        ! on output, everything will be fine and ready for the FCIQMC calculation
        ! to start.

        use FciMCData, only: ll_node, HashIndex, nWalkerHashes
        use hash, only: clear_hash_table, FindWalkerHash
        use DetBitOps, only: tAccumEmptyDet

        integer, intent(in) :: run
        integer :: i, hash_val, PartInd, nwalkers, i_non_core
        integer :: nI(nel)
        real(dp) :: walker_sign(lenof_sign)
        type(ll_node), pointer :: temp_node
        logical :: tSuccess
        integer :: ierr
        character(*), parameter :: this_routine = 'add_core_states_currentdet'
        real(dp), allocatable :: gdata_buf(:, :)
        type(gdata_io_t) :: reorder_handler

        nwalkers = int(TotWalkers)

        associate(rep => cs_replicas(run))
            ! Test that SpawnedParts is going to be big enough
            if (rep%determ_sizes(iProcIndex) > MaxSpawned) then
#ifdef DEBUG_
                write(stderr, *) 'Spawned parts array will not be big enough for &
                    &Semi-Stochastic initialisation'
                write(stderr, *) 'Please increase MEMORYFACSPAWN'
#else
                write(stderr, *) 'Spawned parts array will not be big enough for &
                    &Semi-Stochastic initialisation on task ', iProcIndex
                write(stderr, *) 'Please increase MEMORYFACSPAWN'
#endif
                call stop_all(this_routine, "Insufficient memory assigned")
            end if

            call reorder_handler%init_gdata_io(tAutoAdaptiveShift, tScaleBlooms, tAccumPopsActive, &
              2 * inum_runs, 1, lenof_sign + 1)
            ! we need to reorder the adaptive shift data, too
            ! the maximally required buffer size is the current size of the
            ! determinant list plus the size of the semi-stochastic space (in case
            ! all core-dets are new)
            allocate(gdata_buf(reorder_handler%entry_size(), (nwalkers + rep%determ_sizes(iProcIndex))), &
                      stat=ierr)
            if (ierr /= 0) call stop_all(this_routine, &
                                         "Failed to allocate buffer for global det data")

            ! First find which CurrentDet states are in the core space.
            ! The warning above refers to this bit of code: If a core determinant is not in the
            ! hash table then they won't be found here and the deterministic flag won't be set!
            do i = 1, rep%determ_sizes(iProcIndex)

                tSuccess = .false.
                call decode_bit_det(nI, SpawnedParts(:, i))
                hash_val = FindWalkerHash(nI, nWalkerHashes)
                temp_node => HashIndex(hash_val)
                if (temp_node%ind /= 0) then
                    do while (associated(temp_node))
                        if (all(SpawnedParts(0:nifd, i) == CurrentDets(0:nifd, temp_node%ind))) then
                            tSuccess = .true.
                            PartInd = temp_node%ind
                            exit
                        end if
                        temp_node => temp_node%next
                    end do
                end if
                nullify (temp_node)

                ! Core state i is in CurrentDets.
                if (tSuccess) then
                    call set_flag(CurrentDets(:, PartInd), flag_deterministic(run))
                    ! Copy the amplitude of the state across to SpawnedParts.
                    call extract_sign(CurrentDets(:, PartInd), walker_sign)
                    call encode_sign(SpawnedParts(:, i), walker_sign)
                    ! Add up the already set flags to those to be set
                    SpawnedParts(IlutBits%ind_flag, i) = ior(CurrentDets(IlutBits%ind_flag, PartInd), &
                                                             SpawnedParts(IlutBits%ind_flag, i))
                    ! Cache the accumulated global det data
                    call reorder_handler%write_gdata(gdata_buf, 1, PartInd, i)
                else
                    ! This will be a new state added to CurrentDets.
                    nwalkers = nwalkers + 1
                    ! no auto-adaptive shift data available
                    gdata_buf(:, i) = 0.0_dp
                end if

            end do
            ! Next loop through CurrentDets and move all non-core states to after the last
            ! core state slot in SpawnedParts.
            i_non_core = rep%determ_sizes(iProcIndex)
            do i = 1, int(TotWalkers)
                if (.not. check_determ_flag(CurrentDets(:, i), run)) then
                    i_non_core = i_non_core + 1

                    ! Add a quick test in, to ensure that we don't overflow the
                    ! spawned parts array...
                    if (i_non_core > MaxSpawned) then
#ifdef DEBUG_
                        write(stderr, *) 'Spawned parts array too small for &
                            &semi-stochastic initialisation'
                        write(stderr, *) 'Please increase MEMORYFACSPAWN'
#else
                        write(stderr, *) 'Spawned parts array too small for &
                            &semi-stochastic initialisation on task ', iProcIndex
                        write(stderr, *) 'Please increase MEMORYFACSPAWN'
#endif
                        call stop_all(this_routine, 'Insufficient memory assigned')
                    end if

                    SpawnedParts(0:NIfTot, i_non_core) = CurrentDets(0:NIfTot, i)
                    call reorder_handler%write_gdata(gdata_buf, 1, i, i_non_core)
                end if
            end do
            ! Now copy all the core states in SpawnedParts into CurrentDets.
            ! Note that the amplitude in CurrentDets was copied across, so this is fine.
            do i = 1, nwalkers
                CurrentDets(0:NifTot, i) = SpawnedParts(0:NIfTot, i)
            end do
            ! Re-assign the reordered global det data cached in gdata_buf
            call reorder_handler%read_gdata(gdata_buf, nwalkers)

            ! After reordering the dets, we have to reset all supergroup
            ! idx in the global_det_data
            if (associated(lookup_supergroup_indexer)) then
                do i = 1, nwalkers
                    call decode_bit_det(nI, CurrentDets(:, i))
                    call set_supergroup_idx(&
                        i, lookup_supergroup_indexer%idx_nI(nI))
                end do
            end if

            call clear_hash_table(HashIndex)

            ! Finally, add the indices back into the hash index array.
            do i = 1, nwalkers
                call extract_sign(CurrentDets(:, i), walker_sign)
                ! Don't add the determinant to the hash table if its unoccupied and not
                ! in the core space and not accumulated.
                if (IsUnoccDet(walker_sign) .and. (.not. check_determ_flag(CurrentDets(:, i))) &
                    .and. .not. tAccumEmptyDet(CurrentDets(:, i))) cycle
                call decode_bit_det(nI, CurrentDets(:, i))
                hash_val = FindWalkerHash(nI, nWalkerHashes)
                temp_node => HashIndex(hash_val)
                ! If the first element in the list has not been used.
                if (temp_node%ind == 0) then
                    temp_node%ind = i
                else
                    do while (associated(temp_node%next))
                        temp_node => temp_node%next
                    end do
                    allocate(temp_node%next)
                    nullify (temp_node%next%next)
                    temp_node%next%ind = i
                end if
                nullify (temp_node)

                ! These core states will always stay in the same position.
                if (i <= rep%determ_sizes(iProcIndex)) rep%indices_of_determ_states(i) = i
            end do
        end associate
        TotWalkers = int(nwalkers, int64)
    end subroutine add_core_states_currentdet_hash