!if(tRDMonFly.and.(.not.tExplicitAllRDM)) then ! ! And now for CurrentH ! j = int(WalkersonNodes(i)) * (1+2*lenof_sign) ! call MPIRecv (AllCurrentH(:, 1:WalkersonNodes(i)), j, & ! NodeRoots(i), Tag2, error) !endif
!if(tRDMonFly.and.(.not.tExplicitAllRDM)) then ! j = int(nDets) * (1+2lenof_sign) ! call MPISend (CurrentH(1:1+2lenof_sign, 1:nDets), j, root, Tag2, error) !endif
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=n_int), | intent(inout) | :: | Dets(0:nIfTot,1:nDets) | |||
integer(kind=int64), | intent(in) | :: | nDets | |||
character(len=255), | intent(in), | optional | :: | filename_stem |
SUBROUTINE WriteToPopsfileParOneArr(Dets, nDets, filename_stem)
use constants, only: size_n_int, n_int
use CalcData, only: iPopsFileNoWrite
use MemoryManager, only: TagIntType
integer(int64), intent(in) :: nDets !The number of occupied entries in Dets
integer(kind=n_int), intent(inout) :: Dets(0:nIfTot, 1:nDets)
character(255), intent(in), optional :: filename_stem
INTEGER :: error
integer(int64) :: WalkersonNodes(0:nNodes - 1), writeoutdet
integer(int64) :: node_write_attempts(0:nNodes - 1)
INTEGER :: Tag, Tag2, fTag
INTEGER :: i, j
INTEGER(KIND=n_int), ALLOCATABLE :: Parts(:, :)
real(dp), allocatable :: gdata(:, :)
INTEGER(TagIntType) :: PartsTag = 0
integer :: nMaxDets
integer :: iunit, iunit_2, Initiator_Count, nwrite
integer(int64) :: write_count, write_count_sum
CHARACTER(len=*), PARAMETER :: this_routine = 'WriteToPopsfileParOneArr'
character(255) :: popsfile
real(dp) :: TempSign(lenof_sign)
character(1024) :: out_tmp
character(255) :: identifier
character(12) :: num_tmp
type(timer), save :: write_timer
! manages the output of global determinant data
type(gdata_io_t) :: gdata_write_handler
! size of the global determinant data per determinant
integer :: gdata_size
! Tme the overall popsfile read in
write_timer%timer_name = 'POPS-write'
call set_timer(write_timer)
if (present(filename_stem)) then
identifier = filename_stem
else
identifier = 'POPSFILE'
endif
! If we are using the new popsfile format, then use it!
if (tHDF5PopsWrite) then
call write_popsfile_hdf5()
! And stop timing
call halt_timer(write_timer)
return
end if
CALL MPIBarrier(error) !sync
! write(stdout,*) "Get Here",nDets
! CALL neci_flush(stdout)
!First, make sure we have up-to-date information - again collect AllTotWalkers
! ,AllSumNoatHF and AllSumENum...
!Calculate the energy by summing all on HF and doubles - convert number at HF
! to a real since no int*8 MPI data type
if (.not. t_real_time_fciqmc) then
! Not required for real-time popsfiles: those values have no meaning there
CALL MPISum(SumNoatHF, 1, AllSumNoatHF)
CALL MPISum(SumENum, 1, AllSumENum)
endif
Tag = 125
Tag2 = 126
fTag = 127
!We have to make the distinction here between the number of entries to expect,
!and the number of determinants we are writing out. Since the list is not
!necessarily contiguous any more, we have to calculate Alltotwalkers seperately.
Writeoutdet = 0
do i = 1, int(nDets)
call extract_sign(Dets(:, i), TempSign)
if (sum(abs(tempsign)) > binarypops_min_weight) then
!Count this det in AllTotWalkers
Writeoutdet = Writeoutdet + 1
endif
enddo
writeoutdet = int(writeoutdet / iPopsPartEvery)
call mpisum(writeoutdet, 1, AllTotWalkers)
! We also need to tell the root processor how many particles to expect
! from each node - these are gathered into WalkersonNodes
if (bNodeRoot) then
call MPIAllGather(ndets, node_write_attempts, error, Roots)
call MPIAllGather(writeoutdet, WalkersonNodes, error, Roots)
end if
! We only want to be using space/speed saving devices if we are doing
! them to the best of our ability.
if (tSplitPops .and. .not. tBinPops) then
call stop_all(this_routine, "Split pops files only works with &
&binary pops files")
end if
! set the module variables for global det data i/o
call gdata_write_handler%init_gdata_io(tAutoAdaptiveShift, tScaleBlooms, &
tAccumPopsActive, fvals_size, max_ratio_size, apvals_size)
if (iProcIndex == root) then
! Construct an output string to give feedback about what is
! being done.
write (stdout, *)
write (stdout, '("*********************************")')
write (out_tmp, '("Writing a ", i2, "-bit")') bits_n_int
if (iPopsPartEvery /= 1) out_tmp = trim(out_tmp)//" reduced"
out_tmp = trim(out_tmp)//" POPSFILE"
if (tBinPops) out_tmp = trim(out_tmp)//"BIN"
out_tmp = trim(out_tmp)//"..."
write (stdout, '(a)') trim(adjustl(out_tmp))
! If we aren't outputting all of the walkers, indicate how many
! we are going to output.
!if (iPopsPartEvery /= 1) then
write (num_tmp, '(i12)') AllTotWalkers
write (stdout, '("Writing a total of ",a," determinants.")') &
trim(adjustl(num_tmp))
!end if
! If we are only outputting determinants with a certain weight
if (tBinPops .and. abs(binarypops_min_weight) > 1.0e-12_dp) then
write (num_tmp, '(f12.3)') binarypops_min_weight
write (stdout, '("Only outputting determinants holding >= ",a,&
&" walkers")') trim(adjustl(num_tmp))
end if
write (stdout, '("*********************************")')
write (stdout, *)
! With a normal popsfile, the header is written at the start.
! Thus we need to do that now.
if (.not. tBinPops) then
call get_unique_filename(trim(identifier), tIncrementPops, .true., &
iPopsFileNoWrite, popsfile)
print *, "WRITING", popsfile, tIncrementPops
iunit = get_free_unit()
! We set recl=50000, which allows the line length written to be
! up to 50000 characters long. This allows popsfiles to be
! written in jobs with up to around 2940 MPI processes (the
! default value caused crashes when using over 2000 processes).
open (iunit, file=popsfile, status='replace', recl=50000)
call write_popsfile_header(iunit, AllTotWalkers, &
WalkersonNodes)
end if
end if
write_count = 0
write_count_sum = 0
nMaxDets = int(maxval(node_write_attempts))
gdata_size = gdata_write_handler%entry_size()
allocate (gdata(gdata_size, nMaxDets), stat=error)
call gdata_write_handler%write_gdata(gdata, int(ndets))
if ((tSplitPops .and. bNodeRoot) .or. iProcIndex == root) then
! For a binary popsfile, the actual data is stored more
! compactly.
if (tBinPops) then
if (tSplitPops) then
write (num_tmp, '(i12)') iProcIndex
out_tmp = 'POPSFILEBIN-'//trim(adjustl(num_tmp))
else
out_tmp = trim(identifier)//'BIN'
end if
! We need to get another unit, as not all nodes will go
! through the header-output section.
iunit = get_free_unit()
call get_unique_filename(trim(out_tmp), tIncrementPops, &
.true., iPopsFileNoWrite, popsfile)
open (iunit, file=popsfile, status='replace', &
form='unformatted')
end if
! Are we outputting initiators as well?
if (tPrintInitiators) then
iunit_2 = get_free_unit()
if (tSplitPops) then
write (num_tmp, '(i12)') iProcIndex
out_tmp = 'INITIATORS-'//trim(adjustl(num_tmp))
else
out_tmp = 'INITIATORS'
end if
open (iunit_2, file=trim(out_tmp), status='UNKNOWN')
Initiator_Count = 0
end if
!if(tRDMonFly.and.(.not.tExplicitAllRDM)) then
!iunit_3 = get_free_unit()
!
!if (tSplitPops) then
! write(num_tmp, '(i12)') iProcIndex
! out_tmp = 'RDM_AV_POP-' // trim(adjustl(num_tmp))
! else
! out_tmp = 'RDM_AV_POP'
! end if
! open(iunit_3, file=trim(out_tmp), status='UNKNOWN', form='unformatted')
!endif
! Write out the dets from this node (the head node, unless we
! are in split-pops mode.
do j = 1, int(ndets)
! Count the number of written particles
if (write_pops_det(iunit, iunit_2, Dets(:, j), j, gdata)) then
write_count = write_count + 1
end if
end do
if (.not. tSplitPops) then
! Allocate a temporary array to store the data being received
! from the other nodes.
! BEWARE - this is a potential crash point with a big
! calculation.
! TODO: If this is the end of the run (say on a big machine),
! should we have deallocated spawnedparts by here to
! ensure we have room, or does the deallocated space from
! dealing with freezing give us plenty of room?
allocate (Parts(0:NIfTot, nMaxDets), stat=error)
call LogMemAlloc('Parts', int(nMaxDets, int32) * (NIfTot + 1), &
size_n_int, this_routine, PartsTag, error)
! Loop over the other nodes in the system sequentially, receive
! their walker lists, and output to the popsfiles.
do i = 1, nNodes - 1
! Calculate the amount of data to receive, and get it.
j = int(node_write_attempts(i)) * (NIfTot + 1)
call MPIRecv(Parts(:, 1:node_write_attempts(i)), j, &
NodeRoots(i), Tag, error)
j = int(node_write_attempts(i)) * gdata_size
call MPIRecv(gdata(1:gdata_size, 1:node_write_attempts(i)), j, NodeRoots(i), &
fTag, error)
! SDS: Catherine seems to have disabled writing these out
! so no need to communicate them.
!!!if(tRDMonFly.and.(.not.tExplicitAllRDM)) then
!!! ! And now for CurrentH
!!! j = int(WalkersonNodes(i)) * (1+2*lenof_sign)
!!! call MPIRecv (AllCurrentH(:, 1:WalkersonNodes(i)), j, &
!!! NodeRoots(i), Tag2, error)
!!!endif
! Then write it out in the same way as above.
nwrite = int(node_write_attempts(i))
do j = 1, nwrite
if (write_pops_det(iunit, iunit_2, Parts(:, j), j, gdata)) then
!if(tRDMonFly.and.(.not.tExplicitAllRDM)) then
! write(iunit_3) AllCurrentH(1:1+2*lenof_sign,j)
!endif
write_count = write_count + 1
end if
end do
end do
end if
! Close the output files
call neci_flush(iunit)
close (iunit)
if (tPrintInitiators) close (iunit_2)
else if (bNodeRoot) then
! Send the data on this processor to the root node for outputting
! in a combined popsfile
ASSERT(.not. tSplitPops)
j = int(nDets) * (NIfTot + 1)
call MPISend(Dets(0:NIfTot, 1:nDets), j, root, Tag, error)
! gdata have already been written
j = int(nDets) * gdata_size
call MPISend(gdata(1:gdata_size, 1:nDets), j, root, fTag, error)
!!!if(tRDMonFly.and.(.not.tExplicitAllRDM)) then
!!! j = int(nDets) * (1+2*lenof_sign)
!!! call MPISend (CurrentH(1:1+2*lenof_sign, 1:nDets), j, root, Tag2, error)
!!!endif
end if
! With binary popsfiles, the header is written once the popsfiles
! have been created, so that we can store the number of particles
! actually written, rather than the total within the system.
! Note that for non-binary popsfiles, all particles are written.
if (tBinPops) then
! Get a count of the number of particles written
call MPISum(write_count, write_count_sum)
if (abs(binarypops_min_weight) < 1.0e-12_dp .and. &
write_count_sum /= AllTotwalkers) then
write (stdout, *) 'WARNING: Number of particles written (', &
write_count_sum, ') does not equal AllTotWalkers (', &
AllTotWalkers, ')'
end if
! Now we need to actually write the header
if (iProcIndex == root) then
call get_unique_filename(trim(identifier)//'HEAD', tIncrementPops, &
.true., iPopsFileNoWrite, popsfile)
iunit = get_free_unit()
! We set recl=50000, which allows the line length written to be
! up to 50000 characters long. This allows popsfiles to be
! written in jobs with up to around 2940 MPI processes (the
! default value caused crashes when using over 2000 processes).
open (iunit, file=popsfile, status='replace', recl=50000)
call write_popsfile_header(iunit, write_count_sum, &
WalkersonNodes)
close (iunit)
end if
end if
! And stop timing
call halt_timer(write_timer)
if (allocated(Parts)) then
deallocate (Parts)
call LogMemDealloc('Popsfile', PartsTag)
end if
if (allocated(gdata)) deallocate (gdata)
! Reset some globals
AllSumNoatHF = 0
AllSumENum = 0
AllTotWalkers = 0
end subroutine WriteToPopsfileParOneArr