ReadFromPopsfile Subroutine

public subroutine ReadFromPopsfile(EndPopsList, ReadBatch, CurrWalkers, CurrParts, CurrHF, Dets, DetsLen, pops_nnodes, pops_walkers, PopNIfSgn, PopNel, PopBalanceBlocks, gdata_read_handler, tCalcExtraInfo, filename_stem)

Arguments

Type IntentOptional Attributes Name
integer(kind=int64), intent(in) :: EndPopsList
integer, intent(in) :: ReadBatch
integer(kind=int64), intent(out) :: CurrWalkers
real(kind=dp) :: CurrParts(lenof_sign)
real(kind=dp), intent(out) :: CurrHF(lenof_sign)
integer(kind=n_int), intent(out) :: Dets(0:nIfTot,DetsLen)
integer, intent(in) :: DetsLen
integer, intent(in) :: pops_nnodes
integer(kind=int64), intent(in) :: pops_walkers(0:nProcessors-1)
integer :: PopNIfSgn
integer, intent(inout) :: PopNel
integer, intent(inout) :: PopBalanceBlocks
type(gdata_io_t), intent(in) :: gdata_read_handler
logical, intent(in) :: tCalcExtraInfo
character(len=*), intent(in), optional :: filename_stem

Contents

Source Code


Source Code

    subroutine ReadFromPopsfile(EndPopsList, ReadBatch, CurrWalkers, &
                                CurrParts, CurrHF, Dets, DetsLen, &
                                pops_nnodes, pops_walkers, PopNifSgn, &
                                PopNel, PopBalanceBlocks, gdata_read_handler, &
                                tCalcExtraInfo, filename_stem)

        integer(int64), intent(in) :: EndPopsList  !Number of entries in the POPSFILE.
        integer, intent(in) :: ReadBatch       !Size of the batch of determinants to read in in one go.
        integer(int64), intent(out) :: CurrWalkers    !Number of determinants which end up on a given processor.
        real(dp), intent(out) :: CurrHF(lenof_sign)
        integer, intent(inout) :: PopNel
        ! If true, then calculate the diagonal Hamiltonian elements of the
        ! read-in popsfile, and, if using the linear scaling algorithm,
        ! fill in the walker hash table.
        logical, intent(in) :: tCalcExtraInfo
        integer, intent(in) :: pops_nnodes
        integer, intent(inout) :: PopBalanceBlocks
        integer(int64), intent(in) :: pops_walkers(0:nProcessors - 1)
        type(gdata_io_t), intent(in) :: gdata_read_handler
        real(dp) :: CurrParts(lenof_sign)
        integer :: iunit
        integer(int64) :: i
        integer :: err
        integer(int64) :: AllCurrWalkers
        logical :: formpops, binpops
        real(dp), dimension(lenof_sign) :: SignTemp
        integer :: PopsVersion
        character(len=*), parameter :: this_routine = 'ReadFromPopsfile'
        HElement_t(dp) :: HElemTemp, HElemTemp2
        character(255) :: popsfile
        !variables from header file
        logical :: tPopHPHF, tPopLz, tPop64Bit
        integer :: iPopLenof_sign, iPopIter, PopNIfD, PopNIfSgn, PopNIfFlag
        integer :: PopNIfTot, PopBlockingIter, read_nnodes, Popinum_runs
        integer :: PopRandomHash(2056)
        integer(int64) :: iPopAllTotWalkers
        real(dp) :: PopDiagSft(inum_runs), read_tau, read_psingles
        real(dp) :: read_pparallel, read_ptriples
        real(dp), dimension(lenof_sign) :: PopSumNoatHF
        integer(int64) :: read_walkers_on_nodes(0:nProcessors - 1)
        integer, intent(in) :: DetsLen
        INTEGER(kind=n_int), intent(out) :: Dets(0:nIfTot, DetsLen)
        character(*), intent(in), optional :: filename_stem
        character(12) :: tmp_num
        character(255) :: tmp_char, identifier
        HElement_t(dp) :: PopAllSumENum(inum_runs)
        integer, allocatable :: TempnI(:)
        logical :: trimmed_parts
        real(dp) :: pops_norm_temp
        type(timer), save :: read_timer, process_timer

        ! Tme the overall popsfile read in
        read_timer%timer_name = 'POPS-read'
        process_timer%timer_name = 'POPS-process'
        if (present(filename_stem)) then
            identifier = filename_stem
        else
            identifier = 'POPSFILE'
        endif

        call set_timer(read_timer)
        ! These values might not be present in some popsfile versions, causing
        ! potential undefined behaviour
        read_psingles = 0.0_dp
        read_ptriples = 0.0_dp
        read_pparallel = 0.0_dp
        read_tau = 0.0_dp
        PopDiagSft = 0.0_dp

        ! For the following operations, the mapping of determinants to processors
        ! has to be initialized
        if (.not. t_initialized_roi) call stop_all(this_routine, &
                                                   "Random orbital mapping indices un-initialized at popsfile read")

        if (tHDF5PopsRead) then

            CurrWalkers = read_popsfile_hdf5(dets)

        else

            call open_pops_head(iunit, formpops, binpops, identifier)
            ! Determine version number.
            PopsVersion = FindPopsfileVersion(iunit)
            IF (FormPops) THEN
                if (PopsVersion == 3) then
                    call ReadPopsHeadv3(iunit, tPop64Bit, tPopHPHF, tPopLz, iPopLenof_Sign, PopNel, &
                                        iPopAllTotWalkers, PopDiagSft, PopSumNoatHF, PopAllSumENum, iPopIter, &
                                        PopNIfD, PopNIfSgn, PopNIfFlag, PopNIfTot)
                else
                    call ReadPopsHeadv4(iunit, tPop64Bit, tPopHPHF, tPopLz, iPopLenof_Sign, PopNel, &
                                        iPopAllTotWalkers, PopDiagSft, PopSumNoatHF, PopAllSumENum, iPopIter, &
                                        PopNIfD, PopNIfSgn, Popinum_runs, PopNIfFlag, PopNIfTot, read_tau, &
                                        PopBlockingIter, PopRandomHash, read_psingles, read_pparallel, &
                                        read_ptriples, read_nnodes, read_walkers_on_nodes, PopBalanceBlocks)
                endif

                if (EndPopsList /= iPopAllTotWalkers) then
                    call stop_all(this_routine, "Error in assessing number of entries in POPSFILE")
                endif

            else if (BinPops) then
                ! If we are reading a binary popsfile, then we want to close the
                ! header file and open the file containing the determinants.
                ! We don't need to bother reading in the header, it has already
                ! been done (as it has been for non-binary popsfiles, too).
                if (iProcIndex == root .or. (bNodeRoot .and. tSplitPops)) then
                    ! Close the header file.
                    close (iunit)

                    ! Use the correct pops file name.
                    tmp_char = trim(identifier)//'BIN'
                    if (tSplitPops) then
                        write (tmp_num, '(i12)') iProcIndex
                        tmp_char = trim(tmp_char)//'-'//trim(adjustl(tmp_num))
                    end if
                    call get_unique_filename(trim(tmp_char), tIncrementPops, &
                                             .false., iPopsFileNoRead, popsfile)
                    open (iunit, file=popsfile, status='old', form='unformatted')
                endif

                ! We need to consider the same parameters for particle
                ! distribution
                read_nnodes = pops_nnodes
                read_walkers_on_nodes = pops_walkers
            end if

            allocate (TempnI(PopNel))

            call mpibarrier(err)

            IF (iProcIndex == Root) THEN
                IF (abs(iWeightPopRead) > 1.0e-12_dp) THEN
                    write (stdout, "(A,I15,A,es17.10,A)") "Although ", EndPopsList, &
                        " configurations will be read in, only determinants &
                        &with a weight of over ", iWeightPopRead, " will be stored."
                else
                    write (stdout, "(A,I15,A)") "Reading in a total of ", EndPopsList, " configurations from POPSFILE."
                ENDIF
                !if(abs(ScaleWalkers - 1) > 1.0e-12_dp) then
                !call warning_neci(this_routine,"ScaleWalkers parameter found, but not implemented in POPSFILE v3 - ignoring.")
                !endif
                call neci_flush(stdout)
            ENDIF

            ! If read in particles are removed due to being unoccupied, or
            ! below the threshold, this is recorded here.
            trimmed_parts = .false.

            ! Which of the POPSFILE reading routines are we going to make use of?
            if (tSplitPops) then
                CurrWalkers = read_pops_splitpops(iunit, PopNel, TempnI, BinPops, &
                                                  Dets, DetsLen, PopNIfSgn, &
                                                  gdata_read_handler, trimmed_parts)
            else
                CurrWalkers = read_pops_general(iunit, PopNel, TempnI, BinPops, Dets, &
                                                DetsLen, ReadBatch, EndPopsList, &
                                                PopNIfSgn, gdata_read_handler, trimmed_parts)
                ! The walkers will be redistributed in a default manner. Ignore
                ! the balancing information in the POPSFILE
                PopBalanceBlocks = -1
            end if

            ! Add the contributions to the norm of the popsfile wave function from
            ! all processes.
            pops_norm_temp = pops_norm
            call MPISumAll(pops_norm_temp, pops_norm)

            ! Close the popsfiles.
            if (iProcIndex == Root .or. (tSplitPops .and. bNodeRoot)) then
                close (iunit)
            endif

            ! Test we have still got all determinants
            write (stdout, *) "initial number of walker read-in: CurrWalkers: ", CurrWalkers
            call MPISum(CurrWalkers, 1, AllCurrWalkers)
            if (iProcIndex == Root) then
                if (AllCurrWalkers /= EndPopsList .and. .not. trimmed_parts) then
                    write (stdout, *) "AllCurrWalkers: ", AllCurrWalkers
                    write (stdout, *) "EndPopsList: ", EndPopsList

                    if (tSplitPops) then
                        write (stdout, *)
                        write (stdout, *) '*******************************************'
                        write (stdout, *) 'Currently using: ', nProcessors, ' processors'
                        write (stdout, *)
                        write (stdout, *) 'Using pre-split popsfiles (SPLIT-POPS, &
                                   &POPSFILEBIN-*).'
                        write (stdout, *) 'Ensure that the number of processors and &
                                   &number of POPSFILEs match.'
                        write (stdout, *) '*******************************************'
                    end if

                    call stop_All(this_routine, "Not all walkers accounted for &
                                                &when reading in")
                endif
            endif

        end if

        ! Clear all deterministic and trial flags so that they can be changed later.
        do i = 1, CurrWalkers
            call clr_flag_multi(Dets(:, i), flag_deterministic)
            call clr_flag(Dets(:, i), flag_determ_parent)
            call clr_flag(Dets(:, i), flag_trial)
            call clr_flag(Dets(:, i), flag_connected)
            call clr_flag(Dets(:, i), flag_removed)

            ! store the determinant
            if (tStoredDets) then
                call decode_bit_det(TempnI, Dets(:, i))
                call store_decoding(int(i), TempnI)
            end if
        end do

        call halt_timer(read_timer)
        call set_timer(process_timer)

        if (tCalcExtraInfo) then

            call clear_hash_table(HashIndex)
            call fill_in_hash_table(HashIndex, nWalkerHashes, Dets, int(CurrWalkers), .true.)

            ! Run through all determinants on each node, and calculate the total number of walkers, and noathf
            CurrHF = 0.0_dp
            CurrParts = 0.0_dp
            do i = 1, CurrWalkers
                call extract_sign(Dets(:, i), SignTemp)

                CurrParts = CurrParts + abs(SignTemp)
                if (DetBitEQ(Dets(:, i), iLutRef(:, 1), nifd)) then
                    if (abs(CurrHF(1)) > 1.0e-12_dp) then
                        call stop_all(this_routine, "HF already found, but shouldn't have")
                    endif
                    CurrHF = CurrHF + SignTemp
                    if (.not. tSemiStochastic) then
                        call set_det_diagH(int(i), 0.0_dp)
                        call set_det_offdiagH(int(i), h_cast(0_dp))
                    end if
                else
                    if (.not. tSemiStochastic) then
                        ! Calculate diagonal matrix element
                        call decode_bit_det(TempnI, Dets(:, i))
                        HElemTemp = get_diagonal_matel(TempnI, Dets(:,i))
                        HElemTemp2 = get_off_diagonal_matel(TempnI, Dets(:,i))
                        call set_det_diagH(int(i), real(HElemTemp, dp) - Hii)
                        call set_det_offdiagH(int(i), HElemTemp2)
                    endif
                endif
            enddo

        end if

        ! If we are doing load balancing, do an initial balance now that we
        ! have read the particles in
        ! n.b. This must be done after the hash tables have been initialised
        ! properly, or everything will break!
        if (tLoadBalanceBlocks .and. tCalcExtraInfo) &
            call pops_init_balance_blocks(PopBalanceBlocks)

        call mpibarrier(err)

        if (allocated(TempnI)) deallocate (TempnI)

        call halt_timer(process_timer)

    end subroutine ReadFromPopsfile