Block-Structured AMR Software Framework
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
Loading...
Searching...
No Matches
AMReX_FabArrayBase.H
Go to the documentation of this file.
1#ifndef BL_FABARRAYBASE_H_
2#define BL_FABARRAYBASE_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_BoxArray.H>
10#include <AMReX_Periodicity.H>
11#include <AMReX_Print.H>
12#include <AMReX_Arena.H>
13#include <AMReX_Gpu.H>
14
15#ifdef AMREX_USE_OMP
16#include <omp.h>
17#endif
18
19#include <ostream>
20#include <string>
21#include <utility>
22
23
24namespace amrex {
25
26class MFIter;
27class Geometry;
28class FArrayBox;
29template <typename FAB> class FabFactory;
30template <typename FAB> class FabArray;
31
32namespace EB2 { class IndexSpace; }
33
41{
42 friend class MFIter;
43
44 template <class FAB> friend void FillBoundary (Vector<FabArray<FAB>*> const& mf, const Periodicity& period);
45
46public:
47
48 FabArrayBase () = default;
49 ~FabArrayBase () = default;
50
52 const DistributionMapping& dm,
53 int nvar,
54 int ngrow);
55
57 const DistributionMapping& dm,
58 int nvar,
59 const IntVect& ngrow);
60
61 FabArrayBase (FabArrayBase&& rhs) noexcept = default;
62 FabArrayBase (const FabArrayBase& rhs) = default;
63 FabArrayBase& operator= (const FabArrayBase& rhs) = default;
65
66 void define (const BoxArray& bxs,
67 const DistributionMapping& dm,
68 int nvar,
69 int ngrow);
70
71 void define (const BoxArray& bxs,
72 const DistributionMapping& dm,
73 int nvar,
74 const IntVect& ngrow);
75
77 [[nodiscard]] int nGrow (int direction = 0) const noexcept { return n_grow[direction]; }
78
79 [[nodiscard]] IntVect nGrowVect () const noexcept { return n_grow; }
80
82 [[nodiscard]] int nComp () const noexcept { return n_comp; }
83
85 [[nodiscard]] IndexType ixType () const noexcept { return boxarray.ixType(); }
86
87 //Return whether this FabArray is empty
88 [[nodiscard]] bool empty () const noexcept { return boxarray.empty(); }
89
94 [[nodiscard]] const BoxArray& boxArray () const noexcept { return boxarray; }
95
100 [[nodiscard]] Box box (int K) const noexcept { return boxarray[K]; }
101
106 [[nodiscard]] Box fabbox (int K) const noexcept;
107
109 [[nodiscard]] int size () const noexcept { return static_cast<int>(boxarray.size()); }
110
112 [[nodiscard]] int local_size () const noexcept { return static_cast<int>(indexArray.size()); }
113
115 [[nodiscard]] const Vector<int> &IndexArray () const noexcept { return indexArray; }
116
118 [[nodiscard]] int localindex (int K) const noexcept {
119 auto low
120 = std::lower_bound(indexArray.begin(), indexArray.end(), K);
121 if (low != indexArray.end() && *low == K) {
122 return static_cast<int>(low - indexArray.begin());
123 }
124 else {
125 return -1;
126 }
127 }
128
130 [[nodiscard]] const DistributionMapping& DistributionMap () const noexcept { return distributionMap; }
131
135 [[nodiscard]] bool is_nodal () const noexcept;
139 [[nodiscard]] bool is_nodal (int dir) const noexcept;
143 [[nodiscard]] bool is_cell_centered () const noexcept;
144
145 void setMultiGhost(bool a_multi_ghost) {m_multi_ghost = a_multi_ghost;}
146
147 // These are provided for convenience to keep track of how many
148 // ghost cells are up to date. The number of filled ghost cells
149 // is updated by FillBoundary and ParallelCopy.
150 [[nodiscard]] IntVect nGrowFilled () const noexcept { return n_filled; }
151 void setNGrowFilled (IntVect const& ng) noexcept { n_filled = ng; }
152
154 [[nodiscard]] bool isFusingCandidate () const noexcept;
155
156 //
158 {
159 int size{0};
160 int maxsize{0};
161 Long maxuse{0};
162 Long nuse{0};
163 Long nbuild{0};
164 Long nerase{0};
165 Long bytes{0};
166 Long bytes_hwm{0};
167 std::string name;
168 explicit CacheStats (std::string name_)
169 : name(std::move(name_)) {;}
170 void recordBuild () noexcept {
171 ++size;
172 ++nbuild;
173 maxsize = std::max(maxsize, size);
174 }
175 void recordErase (Long n) noexcept {
176 // n: how many times the item to be deleted has been used.
177 --size;
178 ++nerase;
179 maxuse = std::max(maxuse, n);
180 }
181 void recordUse () noexcept { ++nuse; }
182 void print () const {
183 amrex::Print(Print::AllProcs) << "### " << name << " ###\n"
184 << " tot # of builds : " << nbuild << "\n"
185 << " tot # of erasures: " << nerase << "\n"
186 << " tot # of uses : " << nuse << "\n"
187 << " max cache size : " << maxsize << "\n"
188 << " max # of uses : " << maxuse << "\n";
189 }
190 };
191 //
194 {
199 CopyComTag () noexcept = default;
200 CopyComTag (const Box& db, const Box& sb, int didx, int sidx) noexcept
201 : dbox(db), sbox(sb), dstIndex(didx), srcIndex(sidx) {}
202 bool operator< (const CopyComTag& rhs) const noexcept {
203 return (srcIndex < rhs.srcIndex) || ((srcIndex == rhs.srcIndex) && (
204 (sbox.smallEnd() < rhs.sbox.smallEnd()
205 || ((sbox.smallEnd() == rhs.sbox.smallEnd()) && (
206 (dstIndex < rhs.dstIndex) || ((dstIndex == rhs.dstIndex) && (
207 (dbox.smallEnd() < rhs.dbox.smallEnd()))))))));
208 }
209 //
210 // Some typedefs & helper functions used throughout the code.
211 //
212 using CopyComTagsContainer = std::vector<CopyComTag>;
213
214 using MapOfCopyComTagContainers = std::map<int,CopyComTagsContainer>;
215 };
216 //
217 // Some useful typedefs.
218 //
221 //
223
229 struct BDKey {
230 BDKey () noexcept = default;
231 BDKey (const BoxArray::RefID& baid, const DistributionMapping::RefID& dmid) noexcept
232 : m_ba_id(baid), m_dm_id(dmid) {}
233 bool operator< (const BDKey& rhs) const noexcept {
234 return (m_ba_id < rhs.m_ba_id) ||
235 ((m_ba_id == rhs.m_ba_id) && (m_dm_id < rhs.m_dm_id));
236 }
237 bool operator== (const BDKey& rhs) const noexcept {
238 return m_ba_id == rhs.m_ba_id && m_dm_id == rhs.m_dm_id;
239 }
240 bool operator!= (const BDKey& rhs) const noexcept {
241 return m_ba_id != rhs.m_ba_id || m_dm_id != rhs.m_dm_id;
242 }
243 friend std::ostream& operator<< (std::ostream& os, const BDKey& id);
244 private:
247 };
248
249 [[nodiscard]] BDKey getBDKey () const noexcept {
251 }
252
253 void updateBDKey ();
254
255 //
267
268 //
271 {
272 int fromProc{0};
273 int toProc{0};
274 int fabIndex{0};
275 int fineIndex{0};
276 int srcComp{0};
277 int destComp{0};
278 int nComp{0};
279 int face{0};
281 int fillBoxId{0};
285 };
286
289
292
294 static void Initialize ();
295 static void Finalize ();
302
303 struct FPinfo
304 {
305 FPinfo (const FabArrayBase& srcfa,
306 const FabArrayBase& dstfa,
307 const Box& dstdomain,
308 const IntVect& dstng,
309 const BoxConverter& coarsener,
310 const Box& fdomain,
311 const Box& cdomain,
312 const EB2::IndexSpace* index_space);
313
314 [[nodiscard]] Long bytes () const;
315
319 std::unique_ptr<FabFactory<FArrayBox> > fact_crse_patch;
320 std::unique_ptr<FabFactory<FArrayBox> > fact_fine_patch;
321 //
326 std::unique_ptr<BoxConverter> m_coarsener;
327 //
328 Long m_nuse{0};
329 };
330
331 using FPinfoCache = std::multimap<BDKey,FabArrayBase::FPinfo*>;
332 using FPinfoCacheIter = FPinfoCache::iterator;
333
335
337
338 static const FPinfo& TheFPinfo (const FabArrayBase& srcfa,
339 const FabArrayBase& dstfa,
340 const IntVect& dstng,
341 const BoxConverter& coarsener,
342 const Geometry& fgeom,
343 const Geometry& cgeom,
344 const EB2::IndexSpace*);
345
346 void flushFPinfo (bool no_assertion=false) const;
347
348 //
350 struct CFinfo
351 {
352 CFinfo (const FabArrayBase& finefa,
353 const Geometry& finegm,
354 const IntVect& ng,
355 bool include_periodic,
356 bool include_physbndry);
357
358 [[nodiscard]] Long bytes () const;
359
360 static Box Domain (const Geometry& geom, const IntVect& ng,
361 bool include_periodic, bool include_physbndry);
362
366 //
372 //
373 Long m_nuse{0};
374 };
375
376 using CFinfoCache = std::multimap<BDKey,FabArrayBase::CFinfo*>;
377 using CFinfoCacheIter = CFinfoCache::iterator;
378
380
382
383 static const CFinfo& TheCFinfo (const FabArrayBase& finefa,
384 const Geometry& finegm,
385 const IntVect& ng,
386 bool include_periodic,
387 bool include_physbndry);
388
389 void flushCFinfo (bool no_assertion=false) const;
390
391 //
393 enum CpOp { COPY = 0, ADD = 1 };
394
395 const TileArray* getTileArray (const IntVect& tilesize) const;
396
397 // Memory Usage Tags
398 struct meminfo {
399 Long nbytes = 0L;
400 Long nbytes_hwm = 0L;
401 };
402 static std::map<std::string, meminfo> m_mem_usage;
403
404 static void updateMemUsage (std::string const& tag, Long nbytes, Arena const* ar);
405 static void printMemUsage ();
406 static Long queryMemUsage (const std::string& tag = std::string("All"));
407 static Long queryMemUsageHWM (const std::string& tag = std::string("All"));
408
409 static void pushRegionTag (const char* t);
410 static void pushRegionTag (std::string t);
411 static void popRegionTag ();
412
413 static AMREX_EXPORT std::vector<std::string> m_region_tag;
414 struct RegionTag {
415 RegionTag (const char* t) : tagged(true) { pushRegionTag(t); }
416 RegionTag (const std::string& t) : tagged(true) { pushRegionTag(t); }
417 RegionTag (RegionTag const&) = delete;
418 RegionTag (RegionTag && rhs) noexcept : tagged(rhs.tagged) { rhs.tagged = false; }
419 RegionTag& operator= (RegionTag const&) = delete;
421 ~RegionTag () { if (tagged) { popRegionTag(); } }
422 private:
423 bool tagged = false;
424 };
425
426//#ifndef AMREX_USE_GPU
427//protected:
428//#endif
429
430 void clear ();
431
437 const std::vector<bool>& OwnerShip () const noexcept { return ownership; }
438 bool isOwner (int li) const noexcept { return ownership[li]; }
439
440 //
441 // The data ...
442 //
446 std::vector<bool> ownership;
449 mutable BDKey m_bdkey;
450 IntVect n_filled; // Note that IntVect is zero by default.
451 bool m_multi_ghost = false;
452
453 //
454 // Tiling
455 //
456 // We use tile size as the key for the inner map.
457
458 using TAMap = std::map<std::pair<IntVect,IntVect>, TileArray>;
459 using TACache = std::map<BDKey, TAMap>;
460 //
463 //
464 void buildTileArray (const IntVect& tilesize, TileArray& ta) const;
465 //
467 bool no_assertion=false) const;
468 static void flushTileArrayCache ();
469
471 {
472 // The cache of local and send/recv per FillBoundary() or ParallelCopy().
473 bool m_threadsafe_loc = false;
474 bool m_threadsafe_rcv = false;
475 std::unique_ptr<CopyComTagsContainer> m_LocTags;
476 std::unique_ptr<MapOfCopyComTagContainers> m_SndTags;
477 std::unique_ptr<MapOfCopyComTagContainers> m_RcvTags;
478 };
479
480 void define_fb_metadata (CommMetaData& cmd, const IntVect& nghost, bool cross,
481 const Periodicity& period, bool multi_ghost) const;
482
483 //
485 struct FB
487 {
488 FB (const FabArrayBase& fa, const IntVect& nghost,
489 bool cross, const Periodicity& period,
490 bool enforce_periodicity_only, bool override_sync,
491 bool multi_ghost);
492
497 bool m_epo;
500 //
501 Long m_nuse{0};
502 bool m_multi_ghost = false;
503 //
504#if defined(__CUDACC__) && defined (AMREX_USE_CUDA)
505 CudaGraph<CopyMemory> m_localCopy;
506 CudaGraph<CopyMemory> m_copyToBuffer;
507 CudaGraph<CopyMemory> m_copyFromBuffer;
508#endif
509 //
510 [[nodiscard]] Long bytes () const;
511 private:
512 void define_fb (const FabArrayBase& fa);
513 void define_epo (const FabArrayBase& fa);
514 void define_os (const FabArrayBase& fa);
515 void tag_one_box (int krcv, BoxArray const& ba, DistributionMapping const& dm,
516 bool build_recv_tag);
517 };
518 //
519 using FBCache = std::multimap<BDKey,FabArrayBase::FB*>;
520 using FBCacheIter = FBCache::iterator;
521 //
524 //
525 const FB& getFB (const IntVect& nghost, const Periodicity& period,
526 bool cross=false, bool enforce_periodicity_only = false,
527 bool override_sync = false) const;
528 //
529 void flushFB (bool no_assertion=false) const;
530 static void flushFBCache ();
531
532 //
534 struct CPC
536 {
537 CPC (const FabArrayBase& dstfa, const IntVect& dstng,
538 const FabArrayBase& srcfa, const IntVect& srcng,
539 const Periodicity& period, bool to_ghost_cells_only = false);
540 CPC (const BoxArray& dstba, const DistributionMapping& dstdm,
541 const Vector<int>& dstidx, const IntVect& dstng,
542 const BoxArray& srcba, const DistributionMapping& srcdm,
543 const Vector<int>& srcidx, const IntVect& srcng,
544 const Periodicity& period, int myproc);
545 CPC (const BoxArray& ba, const IntVect& ng,
546 const DistributionMapping& dstdm, const DistributionMapping& srcdm);
547
548 [[nodiscard]] Long bytes () const;
549
555 bool m_tgco;
558 //
559 Long m_nuse{0};
560
561 private:
562 void define (const BoxArray& ba_dst, const DistributionMapping& dm_dst,
563 const Vector<int>& imap_dst,
564 const BoxArray& ba_src, const DistributionMapping& dm_src,
565 const Vector<int>& imap_src,
566 int MyProc = ParallelDescriptor::MyProc());
567 };
568
569 //
570 using CPCache = std::multimap<BDKey,FabArrayBase::CPC*>;
571 using CPCacheIter = CPCache::iterator;
572 //
575 //
576 const CPC& getCPC (const IntVect& dstng, const FabArrayBase& src, const IntVect& srcng,
577 const Periodicity& period, bool to_ghost_cells_only = false) const;
578 //
579 void flushCPC (bool no_assertion=false) const;
580 static void flushCPCache ();
581
582 //
584 struct RB90
586 {
587 RB90 (const FabArrayBase& fa, const IntVect& nghost, Box const& domain);
590 private:
591 void define (const FabArrayBase& fa);
592 };
593 //
594 using RB90Cache = std::multimap<BDKey,FabArrayBase::RB90*>;
595 using RB90CacheIter = RB90Cache::iterator;
596 //
598 //
599 const RB90& getRB90 (const IntVect& nghost, const Box& domain) const;
600 //
601 void flushRB90 (bool no_assertion=false) const;
602 static void flushRB90Cache ();
603
604 //
606 struct RB180
608 {
609 RB180 (const FabArrayBase& fa, const IntVect& nghost, Box const& domain);
612 private:
613 void define (const FabArrayBase& fa);
614 };
615 //
616 using RB180Cache = std::multimap<BDKey,FabArrayBase::RB180*>;
617 using RB180CacheIter = RB180Cache::iterator;
618 //
620 //
621 const RB180& getRB180 (const IntVect& nghost, const Box& domain) const;
622 //
623 void flushRB180 (bool no_assertion=false) const;
624 static void flushRB180Cache ();
625
626 //
628 struct PolarB
630 {
631 PolarB (const FabArrayBase& fa, const IntVect& nghost, Box const& domain);
634 private:
635 void define (const FabArrayBase& fa);
636 };
637 //
638 using PolarBCache = std::multimap<BDKey,FabArrayBase::PolarB*>;
639 using PolarBCacheIter = PolarBCache::iterator;
640 //
642 //
643 const PolarB& getPolarB (const IntVect& nghost, const Box& domain) const;
644 //
645 void flushPolarB (bool no_assertion=false) const;
646 static void flushPolarBCache ();
647
648#ifdef AMREX_USE_GPU
649 //
652 {
653 ParForInfo (const FabArrayBase& fa, const IntVect& nghost, int nthreads);
655
656 std::pair<int*,int*> const& getBlocks () const { return m_nblocks_x; }
657 BoxIndexer const* getBoxes () const { return m_boxes; }
658
659 ParForInfo () = delete;
660 ParForInfo (ParForInfo const&) = delete;
661 ParForInfo (ParForInfo &&) = delete;
662 void operator= (ParForInfo const&) = delete;
663 void operator= (ParForInfo &&) = delete;
664
668 std::pair<int*,int*> m_nblocks_x;
669 BoxIndexer* m_boxes = nullptr;
670 char* m_hp = nullptr;
671 char* m_dp = nullptr;
672 };
673
674 ParForInfo const& getParForInfo (const IntVect& nghost, int nthreads) const;
675
676 static std::multimap<BDKey,ParForInfo*> m_TheParForCache;
677
678 void flushParForInfo (bool no_assertion=false) const; // flushes its own cache
679 static void flushParForCache (); // flushes the entire cache
680
681#endif
682
683 //
685 static std::map<BDKey, int> m_BD_count;
686 //
688 void clearThisBD (bool no_assertion=false) const;
689 //
691 void addThisBD ();
692 //
694 {
699 Long num_build{0};
700
701 void recordBuild () noexcept {
703 ++num_build;
705 }
706 void recordDelete () noexcept {
708 }
709 void recordMaxNumBoxArrays (int n) noexcept {
711 }
712 void recordMaxNumBAUse (int n) noexcept {
713 max_num_ba_use = std::max(max_num_ba_use, n);
714 }
715 void print () const {
716 amrex::Print(Print::AllProcs) << "### FabArray ###\n"
717 << " tot # of builds : " << num_build << "\n"
718 << " max # of FabArrays : " << max_num_fabarrays << "\n"
719 << " max # of BoxArrays : " << max_num_boxarrays << "\n"
720 << " max # of BoxArray uses: " << max_num_ba_use << "\n";
721 }
722 };
724
726
727 [[nodiscard]] static bool getAllocSingleChunk () { return m_alloc_single_chunk; }
728};
729
730namespace detail {
732 : public Arena
733 {
734 public:
735 SingleChunkArena (Arena* a_arena, std::size_t a_size);
736 ~SingleChunkArena () override;
737
738 SingleChunkArena () = delete;
739 SingleChunkArena (const SingleChunkArena& rhs) = delete;
743
744 [[nodiscard]] void* alloc (std::size_t sz) override;
745 void free (void* pt) override;
746
747 // isDeviceAccessible and isHostAccessible can both be true.
748 [[nodiscard]] bool isDeviceAccessible () const override;
749 [[nodiscard]] bool isHostAccessible () const override;
750
751 [[nodiscard]] bool isManaged () const override;
752 [[nodiscard]] bool isDevice () const override;
753 [[nodiscard]] bool isPinned () const override;
754
755 [[nodiscard]] void* data () const noexcept { return (void*) m_root; }
756
757 private:
759 char* m_root = nullptr;
760 char* m_free = nullptr;
761 std::size_t m_size = 0;
762 };
763}
764
765[[nodiscard]] int nComp (FabArrayBase const& fa);
766[[nodiscard]] IntVect nGrowVect (FabArrayBase const& fa);
767[[nodiscard]] BoxArray const& boxArray (FabArrayBase const& fa);
768[[nodiscard]] DistributionMapping const& DistributionMap (FabArrayBase const& fa);
769
770#ifdef BL_USE_MPI
771bool CheckRcvStats (Vector<MPI_Status>& recv_stats, const Vector<std::size_t>& recv_size, int tag);
772#endif
773
774std::ostream& operator<< (std::ostream& os, const FabArrayBase::BDKey& id);
775
776}
777
778#endif
#define AMREX_EXPORT
Definition AMReX_Extension.H:191
A virtual base class for objects that manage their own dynamic memory allocation.
Definition AMReX_Arena.H:100
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:550
IndexType ixType() const noexcept
Return index type of this BoxArray.
Definition AMReX_BoxArray.H:837
RefID getRefID() const noexcept
Return a unique ID of the reference.
Definition AMReX_BoxArray.H:834
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition AMReX_BoxArray.H:597
bool empty() const noexcept
Return whether the BoxArray is empty.
Definition AMReX_BoxArray.H:603
Definition AMReX_Box.H:1193
AMREX_GPU_HOST_DEVICE const IntVectND< dim > & smallEnd() const &noexcept
Get the smallend of the BoxND.
Definition AMReX_Box.H:105
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:41
RefID getRefID() const noexcept
This gives a unique ID of the reference, which is different from dmID above.
Definition AMReX_DistributionMapping.H:379
Definition AMReX_EB2.H:26
Base class for FabArray.
Definition AMReX_FabArrayBase.H:41
const CPC & getCPC(const IntVect &dstng, const FabArrayBase &src, const IntVect &srcng, const Periodicity &period, bool to_ghost_cells_only=false) const
void flushTileArray(const IntVect &tilesize=IntVect::TheZeroVector(), bool no_assertion=false) const
IntVect nGrowVect() const noexcept
Definition AMReX_FabArrayBase.H:79
static AMREX_EXPORT IntVect mfiter_tile_size
Default tilesize in MFIter.
Definition AMReX_FabArrayBase.H:288
ParForInfo const & getParForInfo(const IntVect &nghost, int nthreads) const
static void Finalize()
void clearThisBD(bool no_assertion=false) const
clear BD count and caches associated with this BD, if no other is using this BD.
void flushPolarB(bool no_assertion=false) const
This flushes its own PolarB.
static AMREX_EXPORT bool m_alloc_single_chunk
Definition AMReX_FabArrayBase.H:725
std::multimap< BDKey, FabArrayBase::FPinfo * > FPinfoCache
Definition AMReX_FabArrayBase.H:331
Vector< int > indexArray
Definition AMReX_FabArrayBase.H:445
FabArrayBase()=default
void flushCPC(bool no_assertion=false) const
This flushes its own CPC.
const RB180 & getRB180(const IntVect &nghost, const Box &domain) const
static const FPinfo & TheFPinfo(const FabArrayBase &srcfa, const FabArrayBase &dstfa, const IntVect &dstng, const BoxConverter &coarsener, const Geometry &fgeom, const Geometry &cgeom, const EB2::IndexSpace *)
static std::map< std::string, meminfo > m_mem_usage
Definition AMReX_FabArrayBase.H:402
const Vector< int > & IndexArray() const noexcept
Return constant reference to indices in the FabArray that we have access.
Definition AMReX_FabArrayBase.H:115
void define_fb_metadata(CommMetaData &cmd, const IntVect &nghost, bool cross, const Periodicity &period, bool multi_ghost) const
FabArrayBase(const BoxArray &bxs, const DistributionMapping &dm, int nvar, int ngrow)
static bool getAllocSingleChunk()
Definition AMReX_FabArrayBase.H:727
BDKey m_bdkey
Definition AMReX_FabArrayBase.H:449
void addThisBD()
add the current BD into BD count database
static void printMemUsage()
void setNGrowFilled(IntVect const &ng) noexcept
Definition AMReX_FabArrayBase.H:151
RB90Cache::iterator RB90CacheIter
Definition AMReX_FabArrayBase.H:595
static CacheStats m_FPinfo_stats
Definition AMReX_FabArrayBase.H:336
void flushParForInfo(bool no_assertion=false) const
static std::multimap< BDKey, ParForInfo * > m_TheParForCache
Definition AMReX_FabArrayBase.H:676
void buildTileArray(const IntVect &tilesize, TileArray &ta) const
FBCache::iterator FBCacheIter
Definition AMReX_FabArrayBase.H:520
static void flushRB180Cache()
This flushes the entire cache.
FabArrayBase(FabArrayBase &&rhs) noexcept=default
int n_comp
Definition AMReX_FabArrayBase.H:448
static void pushRegionTag(const char *t)
~FabArrayBase()=default
static CacheStats m_FBC_stats
Definition AMReX_FabArrayBase.H:523
const PolarB & getPolarB(const IntVect &nghost, const Box &domain) const
bool isFusingCandidate() const noexcept
Is this a good candidate for kernel fusing?
static void flushCPCache()
This flusheds the entire cache.
bool is_cell_centered() const noexcept
This tests on whether the FabArray is cell-centered.
bool is_nodal() const noexcept
This tests on whether the FabArray is fully nodal.
bool isOwner(int li) const noexcept
Definition AMReX_FabArrayBase.H:438
static void flushTileArrayCache()
This flushes the entire cache.
IndexType ixType() const noexcept
Return index type.
Definition AMReX_FabArrayBase.H:85
static void pushRegionTag(std::string t)
static FBCache m_TheFBCache
Definition AMReX_FabArrayBase.H:522
static void updateMemUsage(std::string const &tag, Long nbytes, Arena const *ar)
static void flushPolarBCache()
This flushes the entire cache.
int size() const noexcept
Return the number of FABs in the FabArray.
Definition AMReX_FabArrayBase.H:109
void flushCFinfo(bool no_assertion=false) const
std::multimap< BDKey, FabArrayBase::RB180 * > RB180Cache
Definition AMReX_FabArrayBase.H:616
IntVect nGrowFilled() const noexcept
Definition AMReX_FabArrayBase.H:150
FabArrayBase & operator=(const FabArrayBase &rhs)=default
const TileArray * getTileArray(const IntVect &tilesize) const
void flushFPinfo(bool no_assertion=false) const
void define(const BoxArray &bxs, const DistributionMapping &dm, int nvar, int ngrow)
void setMultiGhost(bool a_multi_ghost)
Definition AMReX_FabArrayBase.H:145
const std::vector< bool > & OwnerShip() const noexcept
Return owenership of fabs. The concept of ownership only applies when UPC++ team is used....
Definition AMReX_FabArrayBase.H:437
static PolarBCache m_ThePolarBCache
Definition AMReX_FabArrayBase.H:641
static void popRegionTag()
std::map< std::pair< IntVect, IntVect >, TileArray > TAMap
Definition AMReX_FabArrayBase.H:458
int nGrow(int direction=0) const noexcept
Return the grow factor that defines the region of definition.
Definition AMReX_FabArrayBase.H:77
static CacheStats m_CFinfo_stats
Definition AMReX_FabArrayBase.H:381
static RB180Cache m_TheRB180Cache
Definition AMReX_FabArrayBase.H:619
CPCache::iterator CPCacheIter
Definition AMReX_FabArrayBase.H:571
CopyComTag::CopyComTagsContainer CopyComTagsContainer
Definition AMReX_FabArrayBase.H:219
static AMREX_EXPORT int MaxComp
The maximum number of components to copy() at a time.
Definition AMReX_FabArrayBase.H:291
std::multimap< BDKey, FabArrayBase::FB * > FBCache
Definition AMReX_FabArrayBase.H:519
static Long queryMemUsage(const std::string &tag=std::string("All"))
RB180Cache::iterator RB180CacheIter
Definition AMReX_FabArrayBase.H:617
BDKey getBDKey() const noexcept
Definition AMReX_FabArrayBase.H:249
static CPCache m_TheCPCache
Definition AMReX_FabArrayBase.H:573
CopyComTag::MapOfCopyComTagContainers MapOfCopyComTagContainers
Definition AMReX_FabArrayBase.H:220
static CacheStats m_TAC_stats
Definition AMReX_FabArrayBase.H:462
int localindex(int K) const noexcept
Return local index in the vector of FABs.
Definition AMReX_FabArrayBase.H:118
FabArrayBase(const FabArrayBase &rhs)=default
const DistributionMapping & DistributionMap() const noexcept
Return constant reference to associated DistributionMapping.
Definition AMReX_FabArrayBase.H:130
int local_size() const noexcept
Return the number of local FABs in the FabArray.
Definition AMReX_FabArrayBase.H:112
IntVect n_grow
Definition AMReX_FabArrayBase.H:447
void flushRB90(bool no_assertion=false) const
This flushes its own RB90.
static AMREX_EXPORT IntVect comm_tile_size
communication tile size
Definition AMReX_FabArrayBase.H:301
const RB90 & getRB90(const IntVect &nghost, const Box &domain) const
static const CFinfo & TheCFinfo(const FabArrayBase &finefa, const Geometry &finegm, const IntVect &ng, bool include_periodic, bool include_physbndry)
std::vector< bool > ownership
Definition AMReX_FabArrayBase.H:446
static CFinfoCache m_TheCrseFineCache
Definition AMReX_FabArrayBase.H:379
static Long queryMemUsageHWM(const std::string &tag=std::string("All"))
std::multimap< BDKey, FabArrayBase::CPC * > CPCache
Definition AMReX_FabArrayBase.H:570
FabArrayBase(const BoxArray &bxs, const DistributionMapping &dm, int nvar, const IntVect &ngrow)
std::multimap< BDKey, FabArrayBase::CFinfo * > CFinfoCache
Definition AMReX_FabArrayBase.H:376
bool empty() const noexcept
Definition AMReX_FabArrayBase.H:88
static void Initialize()
Initialize from ParmParse with "fabarray" prefix.
std::multimap< BDKey, FabArrayBase::RB90 * > RB90Cache
Definition AMReX_FabArrayBase.H:594
static void flushRB90Cache()
This flushes the entire cache.
static RB90Cache m_TheRB90Cache
Definition AMReX_FabArrayBase.H:597
void flushRB180(bool no_assertion=false) const
This flushes its own RB180.
static void flushFBCache()
This flushes the entire cache.
CpOp
parallel copy or add
Definition AMReX_FabArrayBase.H:393
@ ADD
Definition AMReX_FabArrayBase.H:393
@ COPY
Definition AMReX_FabArrayBase.H:393
bool m_multi_ghost
Definition AMReX_FabArrayBase.H:451
IntVect n_filled
Definition AMReX_FabArrayBase.H:450
Box box(int K) const noexcept
Return the Kth Box in the BoxArray. That is, the valid region of the Kth grid.
Definition AMReX_FabArrayBase.H:100
static FPinfoCache m_TheFillPatchCache
Definition AMReX_FabArrayBase.H:334
std::map< BDKey, TAMap > TACache
Definition AMReX_FabArrayBase.H:459
static CacheStats m_CPC_stats
Definition AMReX_FabArrayBase.H:574
void flushFB(bool no_assertion=false) const
This flushes its own FB.
static AMREX_EXPORT std::vector< std::string > m_region_tag
Definition AMReX_FabArrayBase.H:413
void define(const BoxArray &bxs, const DistributionMapping &dm, int nvar, const IntVect &ngrow)
static Long bytesOfMapOfCopyComTagContainers(const MapOfCopyComTagContainers &)
DistributionMapping distributionMap
Definition AMReX_FabArrayBase.H:444
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.
FPinfoCache::iterator FPinfoCacheIter
Definition AMReX_FabArrayBase.H:332
friend void FillBoundary(Vector< FabArray< FAB > * > const &mf, const Periodicity &period)
CFinfoCache::iterator CFinfoCacheIter
Definition AMReX_FabArrayBase.H:377
PolarBCache::iterator PolarBCacheIter
Definition AMReX_FabArrayBase.H:639
BoxArray boxarray
Definition AMReX_FabArrayBase.H:443
const FB & getFB(const IntVect &nghost, const Periodicity &period, bool cross=false, bool enforce_periodicity_only=false, bool override_sync=false) const
int nComp() const noexcept
Return number of variables (aka components) associated with each point.
Definition AMReX_FabArrayBase.H:82
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:94
std::multimap< BDKey, FabArrayBase::PolarB * > PolarBCache
Definition AMReX_FabArrayBase.H:638
static AMREX_EXPORT FabArrayStats m_FA_stats
Definition AMReX_FabArrayBase.H:723
static void flushParForCache()
static std::map< BDKey, int > m_BD_count
Keep track of how many FabArrays are built with the same BDKey.
Definition AMReX_FabArrayBase.H:685
static TACache m_TheTileArrayCache
Definition AMReX_FabArrayBase.H:461
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:344
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:73
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE constexpr IntVectND< dim > TheZeroVector() noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition AMReX_IntVect.H:670
Definition AMReX_MFIter.H:57
This provides length of period for periodic domains. 0 means it is not periodic in that direction....
Definition AMReX_Periodicity.H:17
This class provides the user with a few print options.
Definition AMReX_Print.H:35
static constexpr int AllProcs
Definition AMReX_Print.H:38
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:27
Long size() const noexcept
Definition AMReX_Vector.H:50
Definition AMReX_FabArrayBase.H:733
bool isPinned() const override
void * alloc(std::size_t sz) override
bool isHostAccessible() const override
std::size_t m_size
Definition AMReX_FabArrayBase.H:761
SingleChunkArena(Arena *a_arena, std::size_t a_size)
SingleChunkArena & operator=(const SingleChunkArena &rhs)=delete
bool isDeviceAccessible() const override
char * m_root
Definition AMReX_FabArrayBase.H:759
SingleChunkArena(const SingleChunkArena &rhs)=delete
void * data() const noexcept
Definition AMReX_FabArrayBase.H:755
void free(void *pt) override
A pure virtual function for deleting the arena pointed to by pt.
SingleChunkArena(SingleChunkArena &&rhs)=delete
bool isDevice() const override
char * m_free
Definition AMReX_FabArrayBase.H:760
DataAllocator m_dallocator
Definition AMReX_FabArrayBase.H:758
bool isManaged() const override
int MyProc() noexcept
return the rank number local to the current Parallel Context
Definition AMReX_ParallelDescriptor.H:125
Definition AMReX_Amr.cpp:49
int nComp(FabArrayBase const &fa)
DistributionMapping const & DistributionMap(FabArrayBase const &fa)
IntVect nGrowVect(FabArrayBase const &fa)
std::ostream & operator<<(std::ostream &os, AmrMesh const &amr_mesh)
Definition AMReX_AmrMesh.cpp:1236
BoxArray const & boxArray(FabArrayBase const &fa)
Definition AMReX_FabArrayCommI.H:896
Definition AMReX_BoxArray.H:276
Definition AMReX_BoxArray.H:822
Definition AMReX_Box.H:2027
Definition AMReX_DataAllocator.H:9
Definition AMReX_DistributionMapping.H:365
Definition AMReX_FabArrayBase.H:229
BoxArray::RefID m_ba_id
Definition AMReX_FabArrayBase.H:245
bool operator==(const BDKey &rhs) const noexcept
Definition AMReX_FabArrayBase.H:237
bool operator<(const BDKey &rhs) const noexcept
Definition AMReX_FabArrayBase.H:233
bool operator!=(const BDKey &rhs) const noexcept
Definition AMReX_FabArrayBase.H:240
BDKey() noexcept=default
DistributionMapping::RefID m_dm_id
Definition AMReX_FabArrayBase.H:246
friend std::ostream & operator<<(std::ostream &os, const BDKey &id)
coarse/fine boundary
Definition AMReX_FabArrayBase.H:351
Box m_fine_domain
Definition AMReX_FabArrayBase.H:368
IntVect m_ng
Definition AMReX_FabArrayBase.H:369
static Box Domain(const Geometry &geom, const IntVect &ng, bool include_periodic, bool include_physbndry)
CFinfo(const FabArrayBase &finefa, const Geometry &finegm, const IntVect &ng, bool include_periodic, bool include_physbndry)
bool m_include_periodic
Definition AMReX_FabArrayBase.H:370
bool m_include_physbndry
Definition AMReX_FabArrayBase.H:371
Long m_nuse
Definition AMReX_FabArrayBase.H:373
BoxArray ba_cfb
Definition AMReX_FabArrayBase.H:363
BDKey m_fine_bdk
Definition AMReX_FabArrayBase.H:367
Vector< int > fine_grid_idx
local array
Definition AMReX_FabArrayBase.H:365
DistributionMapping dm_cfb
Definition AMReX_FabArrayBase.H:364
parallel copy or add
Definition AMReX_FabArrayBase.H:536
BoxArray m_srcba
Definition AMReX_FabArrayBase.H:556
CPC(const BoxArray &dstba, const DistributionMapping &dstdm, const Vector< int > &dstidx, const IntVect &dstng, const BoxArray &srcba, const DistributionMapping &srcdm, const Vector< int > &srcidx, const IntVect &srcng, const Periodicity &period, int myproc)
IntVect m_srcng
Definition AMReX_FabArrayBase.H:552
bool m_tgco
Definition AMReX_FabArrayBase.H:555
BoxArray m_dstba
Definition AMReX_FabArrayBase.H:557
BDKey m_srcbdk
Definition AMReX_FabArrayBase.H:550
void define(const BoxArray &ba_dst, const DistributionMapping &dm_dst, const Vector< int > &imap_dst, const BoxArray &ba_src, const DistributionMapping &dm_src, const Vector< int > &imap_src, int MyProc=ParallelDescriptor::MyProc())
BDKey m_dstbdk
Definition AMReX_FabArrayBase.H:551
Long m_nuse
Definition AMReX_FabArrayBase.H:559
CPC(const BoxArray &ba, const IntVect &ng, const DistributionMapping &dstdm, const DistributionMapping &srcdm)
CPC(const FabArrayBase &dstfa, const IntVect &dstng, const FabArrayBase &srcfa, const IntVect &srcng, const Periodicity &period, bool to_ghost_cells_only=false)
Periodicity m_period
Definition AMReX_FabArrayBase.H:554
IntVect m_dstng
Definition AMReX_FabArrayBase.H:553
Definition AMReX_FabArrayBase.H:158
void recordUse() noexcept
Definition AMReX_FabArrayBase.H:181
CacheStats(std::string name_)
Definition AMReX_FabArrayBase.H:168
void recordErase(Long n) noexcept
Definition AMReX_FabArrayBase.H:175
std::string name
name of the cache
Definition AMReX_FabArrayBase.H:167
void recordBuild() noexcept
Definition AMReX_FabArrayBase.H:170
void print() const
Definition AMReX_FabArrayBase.H:182
Definition AMReX_FabArrayBase.H:471
bool m_threadsafe_loc
Definition AMReX_FabArrayBase.H:473
bool m_threadsafe_rcv
Definition AMReX_FabArrayBase.H:474
std::unique_ptr< MapOfCopyComTagContainers > m_RcvTags
Definition AMReX_FabArrayBase.H:477
std::unique_ptr< MapOfCopyComTagContainers > m_SndTags
Definition AMReX_FabArrayBase.H:476
std::unique_ptr< CopyComTagsContainer > m_LocTags
Definition AMReX_FabArrayBase.H:475
Used by a bunch of routines when communicating via MPI.
Definition AMReX_FabArrayBase.H:194
Box sbox
Definition AMReX_FabArrayBase.H:196
CopyComTag() noexcept=default
std::vector< CopyComTag > CopyComTagsContainer
Definition AMReX_FabArrayBase.H:212
bool operator<(const CopyComTag &rhs) const noexcept
Definition AMReX_FabArrayBase.H:202
int srcIndex
Definition AMReX_FabArrayBase.H:198
Box dbox
Definition AMReX_FabArrayBase.H:195
int dstIndex
Definition AMReX_FabArrayBase.H:197
std::map< int, CopyComTagsContainer > MapOfCopyComTagContainers
Definition AMReX_FabArrayBase.H:214
FillBoundary.
Definition AMReX_FabArrayBase.H:487
bool m_multi_ghost
Definition AMReX_FabArrayBase.H:502
Long m_nuse
Definition AMReX_FabArrayBase.H:501
void define_os(const FabArrayBase &fa)
IntVect m_ngrow
Definition AMReX_FabArrayBase.H:495
void define_epo(const FabArrayBase &fa)
bool m_cross
Definition AMReX_FabArrayBase.H:496
IndexType m_typ
Definition AMReX_FabArrayBase.H:493
void define_fb(const FabArrayBase &fa)
Periodicity m_period
Definition AMReX_FabArrayBase.H:499
void tag_one_box(int krcv, BoxArray const &ba, DistributionMapping const &dm, bool build_recv_tag)
bool m_epo
Definition AMReX_FabArrayBase.H:497
bool m_override_sync
Definition AMReX_FabArrayBase.H:498
IntVect m_crse_ratio
BoxArray in FabArrayBase may have crse_ratio.
Definition AMReX_FabArrayBase.H:494
FB(const FabArrayBase &fa, const IntVect &nghost, bool cross, const Periodicity &period, bool enforce_periodicity_only, bool override_sync, bool multi_ghost)
Definition AMReX_FabArrayBase.H:304
std::unique_ptr< FabFactory< FArrayBox > > fact_crse_patch
Definition AMReX_FabArrayBase.H:319
DistributionMapping dm_patch
Definition AMReX_FabArrayBase.H:318
BDKey m_srcbdk
Definition AMReX_FabArrayBase.H:322
BoxArray ba_crse_patch
Definition AMReX_FabArrayBase.H:316
std::unique_ptr< FabFactory< FArrayBox > > fact_fine_patch
Definition AMReX_FabArrayBase.H:320
std::unique_ptr< BoxConverter > m_coarsener
Definition AMReX_FabArrayBase.H:326
Box m_dstdomain
Definition AMReX_FabArrayBase.H:324
FPinfo(const FabArrayBase &srcfa, const FabArrayBase &dstfa, const Box &dstdomain, const IntVect &dstng, const BoxConverter &coarsener, const Box &fdomain, const Box &cdomain, const EB2::IndexSpace *index_space)
IntVect m_dstng
Definition AMReX_FabArrayBase.H:325
Long m_nuse
Definition AMReX_FabArrayBase.H:328
BDKey m_dstbdk
Definition AMReX_FabArrayBase.H:323
BoxArray ba_fine_patch
Definition AMReX_FabArrayBase.H:317
Definition AMReX_FabArrayBase.H:694
Long num_build
Definition AMReX_FabArrayBase.H:699
void recordDelete() noexcept
Definition AMReX_FabArrayBase.H:706
int max_num_boxarrays
Definition AMReX_FabArrayBase.H:697
void recordMaxNumBAUse(int n) noexcept
Definition AMReX_FabArrayBase.H:712
void recordMaxNumBoxArrays(int n) noexcept
Definition AMReX_FabArrayBase.H:709
void print() const
Definition AMReX_FabArrayBase.H:715
int max_num_fabarrays
Definition AMReX_FabArrayBase.H:696
int max_num_ba_use
Definition AMReX_FabArrayBase.H:698
int num_fabarrays
Definition AMReX_FabArrayBase.H:695
void recordBuild() noexcept
Definition AMReX_FabArrayBase.H:701
Used for collecting information used in communicating FABs.
Definition AMReX_FabArrayBase.H:271
int srcComp
Definition AMReX_FabArrayBase.H:276
int toProc
Definition AMReX_FabArrayBase.H:273
int fromProc
Definition AMReX_FabArrayBase.H:272
int fabArrayId
Definition AMReX_FabArrayBase.H:280
int procThatHasData
Definition AMReX_FabArrayBase.H:283
int procThatNeedsData
Definition AMReX_FabArrayBase.H:282
int destComp
Definition AMReX_FabArrayBase.H:277
int face
Definition AMReX_FabArrayBase.H:279
Box box
Definition AMReX_FabArrayBase.H:284
int fabIndex
Definition AMReX_FabArrayBase.H:274
int nComp
Definition AMReX_FabArrayBase.H:278
int fineIndex
Definition AMReX_FabArrayBase.H:275
int fillBoxId
Definition AMReX_FabArrayBase.H:281
For ParallelFor(FabArray)
Definition AMReX_FabArrayBase.H:652
ParForInfo(ParForInfo const &)=delete
BATransformer m_bat
Definition AMReX_FabArrayBase.H:665
ParForInfo(ParForInfo &&)=delete
BoxIndexer const * getBoxes() const
Definition AMReX_FabArrayBase.H:657
int m_nthreads
Definition AMReX_FabArrayBase.H:667
char * m_hp
Definition AMReX_FabArrayBase.H:670
void operator=(ParForInfo const &)=delete
std::pair< int *, int * > const & getBlocks() const
Definition AMReX_FabArrayBase.H:656
IntVect m_ng
Definition AMReX_FabArrayBase.H:666
std::pair< int *, int * > m_nblocks_x
Definition AMReX_FabArrayBase.H:668
char * m_dp
Definition AMReX_FabArrayBase.H:671
ParForInfo(const FabArrayBase &fa, const IntVect &nghost, int nthreads)
BoxIndexer * m_boxes
Definition AMReX_FabArrayBase.H:669
Fill polar boundary in spherical coordinates.
Definition AMReX_FabArrayBase.H:630
void define(const FabArrayBase &fa)
Box m_domain
Definition AMReX_FabArrayBase.H:633
PolarB(const FabArrayBase &fa, const IntVect &nghost, Box const &domain)
IntVect m_ngrow
Definition AMReX_FabArrayBase.H:632
Rotate Boundary by 180.
Definition AMReX_FabArrayBase.H:608
void define(const FabArrayBase &fa)
IntVect m_ngrow
Definition AMReX_FabArrayBase.H:610
Box m_domain
Definition AMReX_FabArrayBase.H:611
RB180(const FabArrayBase &fa, const IntVect &nghost, Box const &domain)
Rotate Boundary by 90.
Definition AMReX_FabArrayBase.H:586
RB90(const FabArrayBase &fa, const IntVect &nghost, Box const &domain)
void define(const FabArrayBase &fa)
Box m_domain
Definition AMReX_FabArrayBase.H:589
IntVect m_ngrow
Definition AMReX_FabArrayBase.H:588
Definition AMReX_FabArrayBase.H:414
~RegionTag()
Definition AMReX_FabArrayBase.H:421
bool tagged
Definition AMReX_FabArrayBase.H:423
RegionTag & operator=(RegionTag const &)=delete
RegionTag(RegionTag &&rhs) noexcept
Definition AMReX_FabArrayBase.H:418
RegionTag(const std::string &t)
Definition AMReX_FabArrayBase.H:416
RegionTag(RegionTag const &)=delete
RegionTag(const char *t)
Definition AMReX_FabArrayBase.H:415
Tiling.
Definition AMReX_FabArrayBase.H:258
Vector< int > indexMap
Definition AMReX_FabArrayBase.H:261
Vector< int > localIndexMap
Definition AMReX_FabArrayBase.H:262
Vector< int > localTileIndexMap
Definition AMReX_FabArrayBase.H:263
Vector< int > numLocalTiles
Definition AMReX_FabArrayBase.H:260
Long nuse
Definition AMReX_FabArrayBase.H:259
Vector< Box > tileArray
Definition AMReX_FabArrayBase.H:264
Definition AMReX_FabArrayBase.H:398
Long nbytes
Definition AMReX_FabArrayBase.H:399
Long nbytes_hwm
Definition AMReX_FabArrayBase.H:400