Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_VisMF.H
Go to the documentation of this file.
1
2#ifndef AMREX_VISMF_H_
3#define AMREX_VISMF_H_
4#include <AMReX_Config.H>
5
6#include <AMReX_AsyncOut.H>
7#include <AMReX_BLProfiler.H>
8#include <AMReX_FabArray.H>
9#include <AMReX_FArrayBox.H>
10#include <AMReX_FabConv.H>
11#include <AMReX_NFiles.H>
13#include <AMReX_VisMFBuffer.H>
14
15#include <fstream>
16#include <iterator>
17#include <iostream>
18#include <sstream>
19#include <deque>
20#include <map>
21#include <numeric>
22#include <string>
23#include <type_traits>
24
25namespace amrex {
26
27class IArrayBox;
28
34class VisMF
35 : public VisMFBuffer
36{
37public:
50 explicit VisMF (std::string fafab_name);
51 ~VisMF () = default;
52 VisMF (const VisMF&) = delete;
53 VisMF (VisMF&&) = delete;
54 VisMF& operator= (const VisMF&) = delete;
55 VisMF& operator= (VisMF&&) = delete;
57 struct FabOnDisk
58 {
60 FabOnDisk () = default;
62 FabOnDisk (std::string name, Long offset);
64 std::string m_name;
66 };
68 struct Header
69 {
85 Header (const FabArray<FArrayBox>& mf, VisMF::How how, Version version = Version_v1,
86 bool calcMinMax = true, MPI_Comm = ParallelDescriptor::Communicator());
87
88 ~Header () = default;
89 Header (Header&& rhs) noexcept = default;
90 Header (Header const&) = delete;
91 Header& operator= (Header const&) = delete;
92 Header& operator= (Header &&) = delete;
93
96 int procToWrite = ParallelDescriptor::IOProcessorNumber(),
98 //
99 // The data.
100 //
107 //
108 // These are not defined if they are not in the header
109 //
115 };
116
119 {
120 int rankToRead{-1};
121 int faIndex{-1};
124
125 FabReadLink() = default;
126 FabReadLink(int ranktoread, int faindex, Long fileoffset, const Box &b);
127 };
128
144
150 static std::ifstream *OpenStream(const std::string &fileName);
151 static void CloseStream(const std::string &fileName, bool forceClose = false);
152 static void DeleteStream(const std::string &fileName);
153 static void CloseAllStreams();
154 static bool NoFabHeader(const VisMF::Header &hdr);
155
157 [[nodiscard]] int nComp () const;
159 [[nodiscard]] int nGrow () const;
160 [[nodiscard]] IntVect nGrowVect () const;
162 [[nodiscard]] int size () const;
164 [[nodiscard]] const BoxArray& boxArray () const;
166 [[nodiscard]] Real min (int fabIndex, int nComp) const;
168 [[nodiscard]] Real min (int nComp) const;
170 [[nodiscard]] Real max (int fabIndex, int nComp) const;
172 [[nodiscard]] Real max (int nComp) const;
173
179 [[nodiscard]] const FArrayBox& GetFab (int fabIndex, int compIndex) const;
181 void clear (int fabIndex, int compIndex);
183 void clear (int fabIndex);
185 void clear ();
193 static Long Write (const FabArray<FArrayBox> &mf,
194 const std::string& name,
195 VisMF::How how = NFiles,
196 bool set_ghost = false);
197
198 static void AsyncWrite (const FabArray<FArrayBox>& mf, const std::string& mf_name,
199 bool valid_cells_only = false);
200 static void AsyncWrite (FabArray<FArrayBox>&& mf, const std::string& mf_name,
201 bool valid_cells_only = false);
202
210 static Long WriteOnlyHeader (const FabArray<FArrayBox> & mf,
211 const std::string & mf_name,
212 VisMF::How how = NFiles);
214 static void RemoveFiles(const std::string &name, bool verbose = false);
215
225 static void Read (FabArray<FArrayBox> &mf,
226 const std::string &name,
227 const char *faHeader = nullptr,
228 int coordinatorProc = ParallelDescriptor::IOProcessorNumber(),
229 int allow_empty_mf = 0);
230
232 static bool Exist (const std::string &name);
233
235 static void ReadFAHeader (const std::string &fafabName,
236 Vector<char> &header);
237
239 static bool Check (const std::string &name);
241 static Long FileOffset (std::ostream& os);
243 FArrayBox* readFAB (int idx, const std::string& mf_name);
245 FArrayBox* readFAB (int idx, int icomp);
246
247 static int GetNOutFiles ();
248 static void SetNOutFiles (int noutfiles, MPI_Comm comm = ParallelDescriptor::Communicator());
249
250 static int GetMFFileInStreams () { return nMFFileInStreams; }
251 static void SetMFFileInStreams (int nstreams, MPI_Comm comm = ParallelDescriptor::Communicator());
252
253 static int GetVerbose () { return verbose; }
254 static void SetVerbose (int v) { verbose = v; }
255
256 static VisMF::Header::Version GetHeaderVersion () { return currentVersion; }
258 { currentVersion = version; }
259
260 static bool GetGroupSets () { return groupSets; }
261 static void SetGroupSets (bool groupsets) { groupSets = groupsets; }
262
263 static bool GetSetBuf () { return setBuf; }
264 static void SetSetBuf (bool setbuf) { setBuf = setbuf; }
265
266 static bool GetUseSingleRead () { return useSingleRead; }
267 static void SetUseSingleRead (bool usesingleread) { useSingleRead = usesingleread; }
268
269 static bool GetUseSingleWrite () { return useSingleWrite; }
270 static void SetUseSingleWrite (bool usesinglewrite) { useSingleWrite = usesinglewrite; }
271
272 static bool GetCheckFilePositions () { return checkFilePositions; }
273 static void SetCheckFilePositions (bool cfp) { checkFilePositions = cfp; }
274
275 static bool GetUsePersistentIFStreams () { return usePersistentIFStreams; }
276 static void SetUsePersistentIFStreams (bool usepifs) { usePersistentIFStreams = usepifs; }
277
278 static bool GetUseSynchronousReads () { return useSynchronousReads; }
279 static void SetUseSynchronousReads (bool usepsr) { useSynchronousReads = usepsr; }
280
281 static bool GetUseDynamicSetSelection () { return useDynamicSetSelection; }
282 static void SetUseDynamicSetSelection (bool usedss) { useDynamicSetSelection = usedss; }
283
284 static bool GetNoFlushAfterWrite () { return noFlushAfterWrite; }
285 static void SetNoFlushAfterWrite (bool nfaw) { noFlushAfterWrite = nfaw; }
286
287 static bool GetBarrierAfterLevel () { return barrierAfterLevel; }
288 static void SetBarrierAfterLevel (bool bal) { barrierAfterLevel = bal; }
289
290 static std::string DirName (const std::string& filename);
291 static std::string BaseName (const std::string& filename);
292
293 static void Initialize ();
294 static void Finalize ();
295
296private:
297 static FabOnDisk Write (const FArrayBox& fab,
298 const std::string& filename,
299 std::ostream& os,
300 Long& bytes);
301
302 static Long WriteHeaderDoit (const std::string &mf_name,
303 VisMF::Header const &hdr);
304
305 static Long WriteHeader (const std::string &mf_name,
306 VisMF::Header &hdr,
307 int procToWrite = ParallelDescriptor::IOProcessorNumber(),
309
311 static void FindOffsets (const FabArray<FArrayBox> &mf,
312 const std::string &filePrefix,
313 VisMF::Header &hdr,
314 VisMF::Header::Version whichVersion,
315 NFilesIter &nfi,
324 static FArrayBox *readFAB (int idx,
325 const std::string &mf_name,
326 const Header &hdr,
327 int whichComp = -1);
329 static void readFAB (FabArray<FArrayBox> &mf,
330 int idx,
331 const std::string &mf_name,
332 const Header& hdr);
333
334 static void AsyncWriteDoit (const FabArray<FArrayBox>& mf, const std::string& mf_name,
335 bool is_rvalue, bool valid_cells_only);
336
338 std::string m_fafabname;
340 Header m_hdr;
342 mutable Vector< Vector<FArrayBox*> > m_pa;
348 static AMREX_EXPORT std::map<std::string, VisMF::PersistentIFStream> persistentIFStreams;
350 static AMREX_EXPORT int nOutFiles;
351 static AMREX_EXPORT int nMFFileInStreams;
352
353 static AMREX_EXPORT int verbose;
354 static AMREX_EXPORT VisMF::Header::Version currentVersion;
355 static AMREX_EXPORT bool groupSets;
356 static AMREX_EXPORT bool setBuf;
357 static AMREX_EXPORT bool useSingleRead;
358 static AMREX_EXPORT bool useSingleWrite;
359 static AMREX_EXPORT bool checkFilePositions;
360 static AMREX_EXPORT bool usePersistentIFStreams;
361 static AMREX_EXPORT bool useSynchronousReads;
362 static AMREX_EXPORT bool useDynamicSetSelection;
363 static AMREX_EXPORT bool allowSparseWrites;
364 static AMREX_EXPORT bool noFlushAfterWrite;
365 static AMREX_EXPORT bool barrierAfterLevel;
366};
367
369std::ostream& operator<< (std::ostream& os, const VisMF::FabOnDisk& fod);
371std::istream& operator>> (std::istream& is, VisMF::FabOnDisk& fod);
373std::ostream& operator<< (std::ostream& os, const Vector<VisMF::FabOnDisk>& fa);
375std::istream& operator>> (std::istream& is, Vector<VisMF::FabOnDisk>& fa);
377std::ostream& operator<< (std::ostream& os, const VisMF::Header& hd);
379std::istream& operator>> (std::istream& is, VisMF::Header& hd);
380
390template <typename FAB>
391// This function does work for MultiFab, but let's disable it to avoid confusion.
392std::enable_if_t<std::is_same_v<FAB,IArrayBox>>
393Write (const FabArray<FAB>& fa, const std::string& name)
394{
395 BL_PROFILE("Write(FabArray)");
396 AMREX_ASSERT(!name.empty() && name.back() != '/');
397
398 auto data_descriptor = FAB::getDataDescriptor();
399 int data_bytes = data_descriptor->numBytes();
400
401 bool useSparseFPP = false;
402 const Vector<int> &pmap = fa.DistributionMap().ProcessorMap();
403 std::set<int> procsWithData;
404 Vector<int> procsWithDataVector;
405 for(int i : pmap) {
406 procsWithData.insert(i);
407 }
408 const int nOutFiles = VisMF::GetNOutFiles();
409 if (std::ssize(procsWithData) < nOutFiles) {
410 useSparseFPP = true;
411 for (auto i : procsWithData) {
412 procsWithDataVector.push_back(i);
413 }
414 }
415
416 std::string filePrefix = name + "_D_";
417
418 NFilesIter nfi(nOutFiles, filePrefix, VisMF::GetGroupSets(), VisMF::GetSetBuf());
419
420 if (useSparseFPP) {
421 nfi.SetSparseFPP(procsWithDataVector);
422 } else {
423 nfi.SetDynamic();
424 }
425
426 const auto &fio = FAB::getFABio();
427
428 for ( ; nfi.ReadyToWrite(); ++nfi) {
429 for(MFIter mfi(fa); mfi.isValid(); ++mfi) {
430 const FAB &fab = fa[mfi];
431 {
432 std::stringstream hss;
433 fio.write_header(hss, fab, fab.nComp());
434 auto const tstr = hss.view();
435 nfi.Stream().write(tstr.data(), static_cast<std::streamsize>(tstr.size()));
436 nfi.Stream().flush();
437 }
438 auto const* fabdata = fab.dataPtr();
439#ifdef AMREX_USE_GPU
440 std::unique_ptr<FAB> hostfab;
441 if (fab.arena()->isManaged() || fab.arena()->isDevice()) {
442 hostfab = std::make_unique<FAB>(fab.box(), fab.nComp(), The_Pinned_Arena());
443 Gpu::dtoh_memcpy_async(hostfab->dataPtr(), fab.dataPtr(),
444 fab.size()*sizeof(typename FAB::value_type));
446 fabdata = hostfab->dataPtr();
447 }
448#endif
449 Long writeDataItems = fab.box().numPts() * fa.nComp();
450 Long writeDataSize = writeDataItems * data_bytes;
451 nfi.Stream().write((char *) fabdata, writeDataSize);
452 nfi.Stream().flush();
453 }
454 }
455
456 int coordinatorProc = ParallelDescriptor::IOProcessorNumber();
457 if (nfi.GetDynamic()) {
458 coordinatorProc = nfi.CoordinatorProc();
459 }
460
461 if (ParallelDescriptor::MyProc() == coordinatorProc) {
462 std::string header_file_name = name + "_H";
464 std::ofstream ofs;
465 ofs.rdbuf()->pubsetbuf(io_buffer.dataPtr(), io_buffer.size());
466 ofs.open(header_file_name.c_str(), std::ios::out | std::ios::trunc);
467 if (!ofs.good()) {
468 amrex::FileOpenFailed(header_file_name);
469 }
470
471 ofs << "amrex::FabArray<" << FAB::getClassName() << "> v1.0\n";
472 ofs << fa.nComp() << '\n';
473 ofs << fa.nGrowVect() << '\n';
474 fa.boxArray().writeOn(ofs);
475 ofs << '\n';
476
477 const DistributionMapping& dm = fa.DistributionMap();
478 int nfabs = fa.boxArray().size();
479 int nFiles = NFilesIter::ActualNFiles(nOutFiles);
480 int nprocs = ParallelDescriptor::NProcs();
481
482 Vector<Long> fabBytes(nfabs, 0);
483 std::map<int, Vector<int> > rankBoxOrder;
484 for (int i = 0; i < nfabs; ++i) {
485 std::stringstream hss;
486 FAB tmp(fa.fabbox(i), fa.nComp(), false);
487 fio.write_header(hss, tmp, tmp.nComp());
488 // Size includes header and data
489 fabBytes[i] = static_cast<std::streamoff>(hss.tellp()) + tmp.size() * data_bytes;
490 rankBoxOrder[dm[i]].push_back(i);
491 }
492
493 Vector<int> fileNumbers;
494 if (nfi.GetDynamic()) {
495 fileNumbers = nfi.FileNumbersWritten();
496 } else if (nfi.GetSparseFPP()) {
497 // if sparse, write to (file number = rank)
498 fileNumbers.resize(nprocs);
499 std::iota(fileNumbers.begin(), fileNumbers.end(), 0);
500 } else {
501 fileNumbers.resize(nprocs);
502 for (int i = 0; i < nprocs; ++i) {
503 fileNumbers[i] = NFilesIter::FileNumber(nFiles, i, VisMF::GetGroupSets());
504 }
505 }
506
507 Vector<VisMF::FabOnDisk> fod(nfabs);
508
509 const Vector< Vector<int> > &fileNumbersWriteOrder = nfi.FileNumbersWriteOrder();
510 for (auto const& rv : fileNumbersWriteOrder) {
511 Long currentOffset = 0;
512 for (auto rank : rv) {
513 auto rbo_it = rankBoxOrder.find(rank);
514 if (rbo_it != rankBoxOrder.end()) {
515 Vector<int> const& index = rbo_it->second;
516 int whichFileNumber = fileNumbers[rank];
517 std::string const& whichFileName =
518 VisMF::BaseName(NFilesIter::FileName(whichFileNumber, filePrefix));
519 for (int i : index) {
520 fod[i].m_name = whichFileName;
521 fod[i].m_head = currentOffset;
522 currentOffset += fabBytes[i];
523 }
524 }
525 }
526 }
527 ofs << fod;
528 }
529}
530
532namespace detail {
533template <typename FAB>
534void read_fab (FAB& fab, VisMF::FabOnDisk const& fod, std::string const& name)
535{
536 std::string fullname = VisMF::DirName(name);
537 fullname += fod.m_name;
539 std::ifstream ifs;
540 ifs.rdbuf()->pubsetbuf(io_buffer.dataPtr(), io_buffer.size());
541 ifs.open(fullname.c_str(), std::ios::in | std::ios::binary);
542 if (!ifs.good()) {
543 amrex::FileOpenFailed(fullname);
544 }
545 ifs.seekg(fod.m_head, std::ios::beg);
546 fab.readFrom(ifs);
547}
548}
550
568template <typename FAB>
569// This function does work for MultiFab, but let's disable it to avoid confusion.
570std::enable_if_t<std::is_same_v<FAB,IArrayBox>>
571Read (FabArray<FAB>& fa, const std::string& name)
572{
573 BL_PROFILE("Read(FabArray)");
574 AMREX_ASSERT(!name.empty() && name.back() != '/');
575
576 BoxArray ba;
577 int ncomp;
578 IntVect ngrow;
580 {
581 std::string header_file_name = name + "_H";
582 Vector<char> header_file_chars;
583 ParallelDescriptor::ReadAndBcastFile(header_file_name, header_file_chars);
584 std::string header_file_string(header_file_chars.data());
585 std::stringstream ifs(header_file_string, std::istringstream::in);
586
587 std::string type, version;
588 ifs >> type >> version;
589 AMREX_ASSERT(type == "amrex::FabArray<amrex::IArrayBox>" ||
590 type == "amrex::FabArray<amrex::FArrayBox>");
591 ifs >> ncomp;
592 ifs >> ngrow;
593 ba.readFrom(ifs);
594 ifs >> fod;
595 }
596
597 if (fa.empty()) {
598 fa.define(ba, DistributionMapping{ba}, ncomp, ngrow);
599 } else {
601 }
602
603#ifdef AMREX_USE_MPI
604 const int nopensperfile = VisMF::GetMFFileInStreams(); // # of concurrent readers per file
605 const int myproc = ParallelDescriptor::MyProc();
606 const int coordproc = ParallelDescriptor::IOProcessorNumber();
607
608 int nreqs = 0;
609 int allreadsindex = 0;
610 std::map<std::string, int> filenames; // <filename, allreadsindex>
611
612 const int nboxes = fa.size();
613 const auto& dm = fa.DistributionMap();
614 for (int i = 0; i < nboxes; ++i) {
615 if (myproc == dm[i]) {
616 ++nreqs;
617 }
618 if (myproc == coordproc) {
619 std::string const& fname = fod[i].m_name;
620 auto r =filenames.insert(std::make_pair(fname, allreadsindex));
621 if (r.second) {
622 ++allreadsindex;
623 }
624 }
625 }
626
627 const int readtag = ParallelDescriptor::SeqNum();
628 const int donetag = ParallelDescriptor::SeqNum();
629
630 if (myproc == coordproc) {
631 std::multiset<int> availablefiles; // [whichFile] supports multiple reads/file
632 Vector<std::map<int,std::map<Long,int> > > allreads; // [file]<proc,<seek,index>>
633
634 const auto nfiles = static_cast<int>(filenames.size());
635 for (int i = 0; i < nfiles; ++i) {
636 for (int j = 0; j < nopensperfile; ++j) {
637 availablefiles.insert(i);
638 }
639 }
640 allreads.resize(nfiles);
641 for (int i = 0; i < nboxes; ++i) {
642 const auto whichproc = dm[i];
643 const auto iseekpos = fod[i].m_head;
644 std::string const& fname = fod[i].m_name;
645 auto filenamesiter = filenames.find(fname);
646 if (filenamesiter != filenames.end()) {
647 const int fi = filenamesiter->second;
648 allreads[fi][whichproc].insert(std::make_pair(iseekpos,i));
649 } else {
650 amrex::Error("Error in amrex::Read: filename not found "+fname);
651 }
652 }
653
654 int totalioreqs = nboxes;
655 int reqspending = 0;
656 int iopfileindex;
657 std::deque<int> iopreads;
658 std::set<int> busyprocs;
659 while (totalioreqs > 0) {
660 auto afilesiter = availablefiles.begin();
661 while (afilesiter != availablefiles.end()) {
662 const int fi = *afilesiter;
663 if (allreads[fi].empty()) {
664 availablefiles.erase(fi);
665 afilesiter = availablefiles.begin();
666 continue;
667 }
668 auto whichread = allreads[fi].begin();
669 for ( ; whichread != allreads[fi].end(); ++whichread) {
670 const int tryproc = whichread->first;
671 if (!busyprocs.contains(tryproc)) { // not busy
672 busyprocs.insert(tryproc);
673 Vector<int> vreads;
674 vreads.reserve(whichread->second.size());
675 for (auto const& kv : whichread->second) {
676 vreads.push_back(kv.second);
677 }
678 if (tryproc == coordproc) {
679 iopfileindex = fi;
680 for (auto x : vreads) {
681 iopreads.push_back(x);
682 }
683 } else {
684 ParallelDescriptor::Send(vreads, tryproc, readtag);
685 ++reqspending;
686 }
687 availablefiles.erase(afilesiter);
688 afilesiter = availablefiles.begin();
689 break;
690 }
691 }
692 if (whichread == allreads[fi].end()) {
693 ++afilesiter;
694 } else {
695 allreads[fi].erase(whichread);
696 }
697 }
698
699 while (!iopreads.empty()) {
700 int i = iopreads.front();
701 detail::read_fab(fa[i], fod[i], name);
702 --totalioreqs;
703 iopreads.pop_front();
704 if (iopreads.empty()) {
705 availablefiles.insert(iopfileindex);
706 busyprocs.erase(coordproc);
707 }
708 int doneflag;
709 MPI_Status status;
710 ParallelDescriptor::IProbe(MPI_ANY_SOURCE, donetag, doneflag, status);
711 if (doneflag) {
712 break;
713 }
714 }
715
716 if (reqspending > 0) {
717 Vector<int> idone(2);
718 ParallelDescriptor::Message rmess = ParallelDescriptor::Recv(idone, MPI_ANY_SOURCE,
719 donetag);
720 const int i = idone[0];
721 const int procdone = rmess.pid();
722 totalioreqs -= idone[1];
723 --reqspending;
724 busyprocs.erase(procdone);
725 std::string const& fname = fod[i].m_name;
726 const int fi = filenames.find(fname)->second;
727 availablefiles.insert(fi);
728 }
729 }
730 } else {
731 Vector<int> recreads(nreqs, -1);
732 Vector<int> idone(2);
733 while (nreqs > 0) {
734 ParallelDescriptor::Message rmess = ParallelDescriptor::Recv(recreads, coordproc,
735 readtag);
736 const auto nrmess = static_cast<int>(rmess.count());
737 for (int ir = 0; ir < nrmess; ++ir) {
738 int i = recreads[ir];
739 detail::read_fab(fa[i], fod[i], name);
740 }
741 nreqs -= nrmess;
742 idone[0] = recreads[0];
743 idone[1] = nrmess;
744 ParallelDescriptor::Send(idone, coordproc, donetag);
745 }
746 }
747#else
748 for (MFIter mfi(fa); mfi.isValid(); ++mfi) {
749 detail::read_fab(fa[mfi], fod[mfi.index()], name);
750 }
751#endif
752}
753
754}
755
756#endif /*BL_VISMF_H*/
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_EXPORT
Definition AMReX_Extension.H:191
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1139
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:564
std::ostream & writeOn(std::ostream &) const
Output this BoxArray to a checkpoint file.
Definition AMReX_BoxArray.cpp:482
int readFrom(std::istream &is)
Initialize the BoxArray from the supplied istream. It is an error if the BoxArray has already been in...
Definition AMReX_BoxArray.cpp:468
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition AMReX_BoxArray.H:611
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
const Vector< int > & ProcessorMap() const noexcept
Returns a constant reference to the mapping of boxes in the underlying BoxArray to the CPU that holds...
Definition AMReX_DistributionMapping.cpp:49
A Fortran Array of REALs.
Definition AMReX_FArrayBox.H:231
IntVect nGrowVect() const noexcept
Definition AMReX_FabArrayBase.H:80
int size() const noexcept
Return the number of FABs in the FabArray.
Definition AMReX_FabArrayBase.H:110
const DistributionMapping & DistributionMap() const noexcept
Return constant reference to associated DistributionMapping.
Definition AMReX_FabArrayBase.H:131
bool empty() const noexcept
Definition AMReX_FabArrayBase.H:89
Box fabbox(int K) const noexcept
Return the Kth FABs Box in the FabArray. That is, the region the Kth fab is actually defined on.
Definition AMReX_FabArrayBase.cpp:217
int nComp() const noexcept
Return number of variables (aka components) associated with each point.
Definition AMReX_FabArrayBase.H:83
const BoxArray & boxArray() const noexcept
Return a constant reference to the BoxArray that defines the valid region associated with this FabArr...
Definition AMReX_FabArrayBase.H:95
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:351
void define(const BoxArray &bxs, const DistributionMapping &dm, int nvar, int ngrow, const MFInfo &info=MFInfo(), const FabFactory< FAB > &factory=DefaultFabFactory< FAB >())
Define this FabArray identically to that performed by the constructor having an analogous function si...
Definition AMReX_FabArray.H:2183
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:172
This class encapsulates writing to nfiles.
Definition AMReX_NFiles.H:27
Vector< int > FileNumbersWritten()
these are the file numbers to which each rank wrote [rank] a rank only writes to one file
Definition AMReX_NFiles.cpp:518
static int ActualNFiles(int nOutFiles)
this returns the actual number of files used the range [1, nProcs] is enforced
Definition AMReX_NFiles.H:130
bool GetDynamic() const
Definition AMReX_NFiles.H:54
bool GetSparseFPP() const
Definition AMReX_NFiles.H:65
std::fstream & Stream()
Definition AMReX_NFiles.H:98
void SetSparseFPP(const Vector< int > &ranksToWrite)
call this to use a file per process for a sparse set of writers. ranksToWrite.size() will set noutfil...
Definition AMReX_NFiles.cpp:115
bool ReadyToWrite(bool appendFirst=false)
if appendFirst is true, the first set for this iterator will open the files in append mode
Definition AMReX_NFiles.cpp:205
const Vector< Vector< int > > & FileNumbersWriteOrder() const
these are the order of ranks which wrote to each file [filenumber][ranks in order they wrote to filen...
Definition AMReX_NFiles.H:198
int FileNumber() const
Definition AMReX_NFiles.H:161
const std::string & FileName() const
Definition AMReX_NFiles.H:160
int CoordinatorProc() const
this is the processor coordinating dynamic set selection
Definition AMReX_NFiles.H:184
void SetDynamic(int deciderproc=-1)
call this to use dynamic set selection deciderproc defaults to nprocs - 1 if < 0
Definition AMReX_NFiles.cpp:70
Hold the description and status of communication data.
Definition AMReX_ParallelDescriptor.H:57
A Descriptor of the Real Type.
Definition AMReX_FabConv.H:105
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
T * dataPtr() noexcept
get access to the underlying data pointer
Definition AMReX_Vector.H:49
Long size() const noexcept
Definition AMReX_Vector.H:53
Definition AMReX_VisMFBuffer.H:13
static Long GetIOBufferSize()
Definition AMReX_VisMFBuffer.H:26
Vector< Setbuf_Char_Type > IO_Buffer
A simple character buffer for setbuf() usage.
Definition AMReX_VisMFBuffer.H:24
File I/O for FabArray<FArrayBox>. Wrapper class for reading/writing FabArray<FArrayBox> objects to di...
Definition AMReX_VisMF.H:36
IntVect nGrowVect() const
Definition AMReX_VisMF.cpp:496
static void SetGroupSets(bool groupsets)
Definition AMReX_VisMF.H:261
Real max(int fabIndex, int nComp) const
The max of the FAB (in valid region) at specified index and component.
Definition AMReX_VisMF.cpp:539
static void SetNOutFiles(int noutfiles, MPI_Comm comm=ParallelDescriptor::Communicator())
Definition AMReX_VisMF.cpp:157
static Long FileOffset(std::ostream &os)
The file offset of the passed ostream.
Definition AMReX_VisMF.cpp:580
static std::ifstream * OpenStream(const std::string &fileName)
Open the stream if it is not already open Close the stream if not persistent or forced Close all open...
Definition AMReX_VisMF.cpp:2257
static bool GetUseDynamicSetSelection()
Definition AMReX_VisMF.H:281
static bool GetSetBuf()
Definition AMReX_VisMF.H:263
const BoxArray & boxArray() const
The BoxArray of the on-disk FabArray<FArrayBox>.
Definition AMReX_VisMF.cpp:508
static void SetUseSingleWrite(bool usesinglewrite)
Definition AMReX_VisMF.H:270
static void SetHeaderVersion(VisMF::Header::Version version)
Definition AMReX_VisMF.H:257
static void Read(FabArray< FArrayBox > &mf, const std::string &name, const char *faHeader=nullptr, int coordinatorProc=ParallelDescriptor::IOProcessorNumber(), int allow_empty_mf=0)
Read a FabArray<FArrayBox> from disk written using VisMF::Write(). If the FabArray<FArrayBox> fafab h...
Definition AMReX_VisMF.cpp:1582
static void CloseAllStreams()
Definition AMReX_VisMF.cpp:2307
int size() const
Definition AMReX_VisMF.cpp:502
static int GetVerbose()
Definition AMReX_VisMF.H:253
static void SetNoFlushAfterWrite(bool nfaw)
Definition AMReX_VisMF.H:285
static VisMF::Header::Version GetHeaderVersion()
Definition AMReX_VisMF.H:256
static void SetCheckFilePositions(bool cfp)
Definition AMReX_VisMF.H:273
static void SetUseSingleRead(bool usesingleread)
Definition AMReX_VisMF.H:267
FArrayBox * readFAB(int idx, const std::string &mf_name)
Read the entire fab (all components).
Definition AMReX_VisMF.cpp:596
static bool GetCheckFilePositions()
Definition AMReX_VisMF.H:272
const FArrayBox & GetFab(int fabIndex, int compIndex) const
The FAB at the specified index and component. Reads it from disk if necessary. This reads only the sp...
Definition AMReX_VisMF.cpp:564
static void DeleteStream(const std::string &fileName)
Definition AMReX_VisMF.cpp:2296
static Long Write(const FabArray< FArrayBox > &mf, const std::string &name, VisMF::How how=NFiles, bool set_ghost=false)
Write a FabArray<FArrayBox> to disk in a "smart" way. Returns the total number of bytes written on th...
Definition AMReX_VisMF.cpp:980
static bool GetUseSynchronousReads()
Definition AMReX_VisMF.H:278
static int GetMFFileInStreams()
Definition AMReX_VisMF.H:250
static int GetNOutFiles()
Definition AMReX_VisMF.cpp:169
static bool GetNoFlushAfterWrite()
Definition AMReX_VisMF.H:284
static bool GetGroupSets()
Definition AMReX_VisMF.H:260
static std::string DirName(const std::string &filename)
Definition AMReX_VisMF.cpp:626
static Long WriteOnlyHeader(const FabArray< FArrayBox > &mf, const std::string &mf_name, VisMF::How how=NFiles)
Write only the header-file corresponding to FabArray<FArrayBox> to disk without the corresponding FAB...
Definition AMReX_VisMF.cpp:1182
VisMF(VisMF &&)=delete
int nGrow() const
The grow factor of the on-disk FabArray<FArrayBox>.
Definition AMReX_VisMF.cpp:490
static bool GetUsePersistentIFStreams()
Definition AMReX_VisMF.H:275
static void SetUseSynchronousReads(bool usepsr)
Definition AMReX_VisMF.H:279
static void SetVerbose(int v)
Definition AMReX_VisMF.H:254
VisMF & operator=(const VisMF &)=delete
static void CloseStream(const std::string &fileName, bool forceClose=false)
Definition AMReX_VisMF.cpp:2279
How
How we write out FabArray<FArrayBox>s. These are deprecated, we always use NFiles....
Definition AMReX_VisMF.H:43
@ NFiles
Definition AMReX_VisMF.H:43
@ OneFilePerCPU
Definition AMReX_VisMF.H:43
static void AsyncWrite(const FabArray< FArrayBox > &mf, const std::string &mf_name, bool valid_cells_only=false)
Definition AMReX_VisMF.cpp:2313
static void Finalize()
Definition AMReX_VisMF.cpp:150
void clear()
Delete()s all the FABs.
Definition AMReX_VisMF.cpp:2225
static void SetBarrierAfterLevel(bool bal)
Definition AMReX_VisMF.H:288
static bool Exist(const std::string &name)
Does FabArray exist?
Definition AMReX_VisMF.cpp:2087
static bool GetBarrierAfterLevel()
Definition AMReX_VisMF.H:287
static void SetMFFileInStreams(int nstreams, MPI_Comm comm=ParallelDescriptor::Communicator())
Definition AMReX_VisMF.cpp:163
static void ReadFAHeader(const std::string &fafabName, Vector< char > &header)
Read only the header of a FabArray, header will be resized here.
Definition AMReX_VisMF.cpp:2101
~VisMF()=default
static void Initialize()
Definition AMReX_VisMF.cpp:110
static bool GetUseSingleRead()
Definition AMReX_VisMF.H:266
static std::string BaseName(const std::string &filename)
Definition AMReX_VisMF.cpp:608
static void RemoveFiles(const std::string &name, bool verbose=false)
this will remove nfiles associated with name and the header
Definition AMReX_VisMF.cpp:1373
int nComp() const
The number of components in the on-disk FabArray<FArrayBox>.
Definition AMReX_VisMF.cpp:484
Real min(int fabIndex, int nComp) const
The min of the FAB (in valid region) at specified index and component.
Definition AMReX_VisMF.cpp:514
static bool Check(const std::string &name)
Check if the multifab is ok, false is returned if not ok.
Definition AMReX_VisMF.cpp:2112
static void SetUsePersistentIFStreams(bool usepifs)
Definition AMReX_VisMF.H:276
static void SetUseDynamicSetSelection(bool usedss)
Definition AMReX_VisMF.H:282
static bool GetUseSingleWrite()
Definition AMReX_VisMF.H:269
VisMF(const VisMF &)=delete
static void SetSetBuf(bool setbuf)
Definition AMReX_VisMF.H:264
static bool NoFabHeader(const VisMF::Header &hdr)
Definition AMReX_VisMF.cpp:2235
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
amrex_long Long
Definition AMReX_INT.H:30
Arena * The_Pinned_Arena()
Definition AMReX_Arena.cpp:860
int MyProc() noexcept
Definition AMReX_ParallelDescriptor.H:128
void ReadAndBcastFile(const std::string &filename, Vector< char > &charBuf, bool bExitOnError, const MPI_Comm &comm)
Definition AMReX_ParallelDescriptor.cpp:1495
int NProcs() noexcept
Definition AMReX_ParallelDescriptor.H:255
int IOProcessorNumber() noexcept
The MPI rank number of the I/O Processor (probably rank 0). This rank is usually used to write to std...
Definition AMReX_ParallelDescriptor.H:279
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
void dtoh_memcpy_async(void *p_h, const void *p_d, const std::size_t sz) noexcept
Definition AMReX_GpuDevice.H:435
MPI_Comm Communicator() noexcept
Definition AMReX_ParallelDescriptor.H:223
Message Send(const T *buf, size_t n, int dst_pid, int tag)
Definition AMReX_ParallelDescriptor.H:1193
int SeqNum() noexcept
Returns sequential message sequence numbers, usually used as tags for send/recv.
Definition AMReX_ParallelDescriptor.H:696
void IProbe(int, int, int &, MPI_Status &)
Definition AMReX_ParallelDescriptor.cpp:1222
Message Recv(T *, size_t n, int pid, int tag)
Definition AMReX_ParallelDescriptor.H:1235
int MPI_Comm
Definition AMReX_ccse-mpi.H:51
Definition AMReX_Amr.cpp:50
std::ostream & operator<<(std::ostream &os, AmrMesh const &amr_mesh)
Stream helper; forwards to the friend declared inside AmrMesh.
Definition AMReX_AmrMesh.cpp:1306
void FileOpenFailed(const std::string &file)
Output a message and abort when couldn't open the file.
Definition AMReX_Utility.cpp:137
std::enable_if_t< std::is_same_v< FAB, IArrayBox > > Write(const FabArray< FAB > &fa, const std::string &name)
Write iMultiFab/FabArray<IArrayBox>
Definition AMReX_VisMF.H:393
bool match(const BoxArray &x, const BoxArray &y)
Note that two BoxArrays that match are not necessarily equal.
Definition AMReX_BoxArray.cpp:1930
void Error(const std::string &msg)
Print out message to cerr and exit via amrex::Abort().
Definition AMReX.cpp:235
std::enable_if_t< std::is_same_v< FAB, IArrayBox > > Read(FabArray< FAB > &fa, const std::string &name)
Read iMultiFab/FabArray<IArrayBox>
Definition AMReX_VisMF.H:571
std::istream & operator>>(std::istream &is, BoxND< dim > &bx)
Read from istream.
Definition AMReX_Box.H:1825
__host__ __device__ Dim3 end(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2015
A structure containing info regarding an on-disk FAB.
Definition AMReX_VisMF.H:58
Long m_head
Offset to start of FAB in file.
Definition AMReX_VisMF.H:65
std::string m_name
The two data values in a FabOnDisk structure.
Definition AMReX_VisMF.H:64
FabOnDisk()=default
The default constructor – null out all fields.
An on-disk FabArray<FArrayBox> contains this info in a header file.
Definition AMReX_VisMF.H:69
Header & operator=(Header const &)=delete
Vector< Vector< Real > > m_min
The min()s of each component of FABs. [findex][comp].
Definition AMReX_VisMF.H:110
Header(Header const &)=delete
Vector< Real > m_famax
The max()s of each component of the FabArray. [comp].
Definition AMReX_VisMF.H:113
void CalculateMinMax(const FabArray< FArrayBox > &mf, int procToWrite=ParallelDescriptor::IOProcessorNumber(), MPI_Comm=ParallelDescriptor::Communicator())
Calculate the min and max arrays.
Definition AMReX_VisMF.cpp:746
Version
The versions of the FabArray<FArrayBox> Header code.
Definition AMReX_VisMF.H:71
@ NoFabHeader_v1
-— no fab headers, no fab mins or maxes
Definition AMReX_VisMF.H:76
@ NoFabHeaderMinMax_v1
Definition AMReX_VisMF.H:77
@ Undefined_v1
-— undefined
Definition AMReX_VisMF.H:72
@ Version_v1
Definition AMReX_VisMF.H:73
@ NoFabHeaderFAMinMax_v1
Definition AMReX_VisMF.H:79
int m_vers
The version of the Header.
Definition AMReX_VisMF.H:101
BoxArray m_ba
The BoxArray of the MF.
Definition AMReX_VisMF.H:105
Vector< Vector< Real > > m_max
The max()s of each component of FABs. [findex][comp].
Definition AMReX_VisMF.H:111
Header()
The default constructor.
How m_how
How the MF was written to disk.
Definition AMReX_VisMF.H:102
Header(Header &&rhs) noexcept=default
Vector< FabOnDisk > m_fod
FabOnDisk info for contained FABs.
Definition AMReX_VisMF.H:106
int m_ncomp
Number of components in MF.
Definition AMReX_VisMF.H:103
RealDescriptor m_writtenRD
Definition AMReX_VisMF.H:114
IntVect m_ngrow
The number of ghost cells in MF.
Definition AMReX_VisMF.H:104
Vector< Real > m_famin
The min()s of each component of the FabArray. [comp].
Definition AMReX_VisMF.H:112
This structure is used to store file ifstreams that remain open.
Definition AMReX_VisMF.H:131
~PersistentIFStream()
Definition AMReX_VisMF.cpp:2246
VisMF::IO_Buffer ioBuffer
Definition AMReX_VisMF.H:135
PersistentIFStream(PersistentIFStream &&)=delete
PersistentIFStream & operator=(PersistentIFStream const &)=delete
std::streampos currentPosition
Definition AMReX_VisMF.H:133
std::ifstream * pstr
Definition AMReX_VisMF.H:132
bool isOpen
Definition AMReX_VisMF.H:134
PersistentIFStream(PersistentIFStream const &)=delete
Definition AMReX_ccse-mpi.H:55