Line data Source code
1 : #ifndef BL_FABARRAYBASE_H_
2 : #define BL_FABARRAYBASE_H_
3 : #include <AMReX_Config.H>
4 :
5 : #include <AMReX_BoxArray.H>
6 : #include <AMReX_DataAllocator.H>
7 : #include <AMReX_DistributionMapping.H>
8 : #include <AMReX_ParallelDescriptor.H>
9 : #include <AMReX_ParallelReduce.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 :
24 : namespace amrex {
25 :
26 : class MFIter;
27 : class Geometry;
28 : class FArrayBox;
29 : template <typename FAB> class FabFactory;
30 : template <typename FAB> class FabArray;
31 :
32 : namespace EB2 { class IndexSpace; }
33 :
34 : /**
35 : * \brief Base class for FabArray.
36 : *
37 : * Not to be confused with FArrayBox or `FAB` shorthands.
38 : * Can be read as FArrayBox-like Array Base.
39 : */
40 11120 : class FabArrayBase
41 : {
42 : friend class MFIter;
43 :
44 : template <class FAB> friend void FillBoundary (Vector<FabArray<FAB>*> const& mf, const Periodicity& period);
45 :
46 : public:
47 :
48 403505 : FabArrayBase () = default;
49 453172 : ~FabArrayBase () = default;
50 :
51 : FabArrayBase (const BoxArray& bxs,
52 : const DistributionMapping& dm,
53 : int nvar,
54 : int ngrow);
55 :
56 : FabArrayBase (const BoxArray& bxs,
57 : const DistributionMapping& dm,
58 : int nvar,
59 : const IntVect& ngrow);
60 :
61 0 : FabArrayBase (FabArrayBase&& rhs) noexcept = default;
62 : FabArrayBase (const FabArrayBase& rhs) = default;
63 : FabArrayBase& operator= (const FabArrayBase& rhs) = default;
64 : FabArrayBase& operator= (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 :
76 : //! Return the grow factor that defines the region of definition.
77 5795380 : [[nodiscard]] int nGrow (int direction = 0) const noexcept { return n_grow[direction]; }
78 :
79 796 : [[nodiscard]] IntVect nGrowVect () const noexcept { return n_grow; }
80 :
81 : //! Return number of variables (aka components) associated with each point.
82 326736 : [[nodiscard]] int nComp () const noexcept { return n_comp; }
83 :
84 : //! Return index type.
85 330 : [[nodiscard]] IndexType ixType () const noexcept { return boxarray.ixType(); }
86 :
87 : //Return whether this FabArray is empty
88 12952 : [[nodiscard]] bool empty () const noexcept { return boxarray.empty(); }
89 :
90 : /**
91 : * \brief Return a constant reference to the BoxArray that defines the
92 : * valid region associated with this FabArray.
93 : */
94 97992 : [[nodiscard]] const BoxArray& boxArray () const noexcept { return boxarray; }
95 :
96 : /**
97 : * \brief Return the Kth Box in the BoxArray.
98 : * That is, the valid region of the Kth grid.
99 : */
100 246523 : [[nodiscard]] Box box (int K) const noexcept { return boxarray[K]; }
101 :
102 : /**
103 : * \brief Return the Kth FABs Box in the FabArray.
104 : * That is, the region the Kth fab is actually defined on.
105 : */
106 : [[nodiscard]] Box fabbox (int K) const noexcept;
107 :
108 : //! Return the number of FABs in the FabArray.
109 : [[nodiscard]] int size () const noexcept { return static_cast<int>(boxarray.size()); }
110 :
111 : //! Return the number of local FABs in the FabArray.
112 365725 : [[nodiscard]] int local_size () const noexcept { return static_cast<int>(indexArray.size()); }
113 :
114 : //! Return constant reference to indices in the FabArray that we have access.
115 : [[nodiscard]] const Vector<int> &IndexArray () const noexcept { return indexArray; }
116 :
117 : //! Return local index in the vector of FABs.
118 27148 : [[nodiscard]] int localindex (int K) const noexcept {
119 : auto low
120 27148 : = std::lower_bound(indexArray.begin(), indexArray.end(), K);
121 27148 : if (low != indexArray.end() && *low == K) {
122 27148 : return static_cast<int>(low - indexArray.begin());
123 : }
124 : else {
125 0 : return -1;
126 : }
127 : }
128 :
129 : //! Return constant reference to associated DistributionMapping.
130 91018 : [[nodiscard]] const DistributionMapping& DistributionMap () const noexcept { return distributionMap; }
131 :
132 : /**
133 : * \brief This tests on whether the FabArray is fully nodal.
134 : */
135 : [[nodiscard]] bool is_nodal () const noexcept;
136 : /**
137 : * \brief This tests on whether the FabArray is nodal in direction dir.
138 : */
139 : [[nodiscard]] bool is_nodal (int dir) const noexcept;
140 : /**
141 : * \brief This tests on whether the FabArray is cell-centered.
142 : */
143 : [[nodiscard]] bool is_cell_centered () const noexcept;
144 :
145 233097 : 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 :
153 : //! Is this a good candidate for kernel fusing?
154 : [[nodiscard]] bool isFusingCandidate () const noexcept;
155 :
156 : //
157 : struct CacheStats
158 : {
159 : int size{0}; //!< current size: nbuild - nerase
160 : int maxsize{0}; //!< highest water mark of size
161 : Long maxuse{0}; //!< max # of uses of a cached item
162 : Long nuse{0}; //!< # of uses of the whole cache
163 : Long nbuild{0}; //!< # of build operations
164 : Long nerase{0}; //!< # of erase operations
165 : Long bytes{0};
166 : Long bytes_hwm{0};
167 : std::string name; //!< name of the cache
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 : //
192 : //! Used by a bunch of routines when communicating via MPI.
193 : struct CopyComTag
194 : {
195 : Box dbox;
196 : Box sbox;
197 : int dstIndex;
198 : int srcIndex;
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 : //
219 : using CopyComTagsContainer = CopyComTag::CopyComTagsContainer;
220 : using MapOfCopyComTagContainers = CopyComTag::MapOfCopyComTagContainers;
221 : //
222 : static Long bytesOfMapOfCopyComTagContainers (const MapOfCopyComTagContainers&);
223 :
224 : /**
225 : * Key for unique combination of BoxArray and DistributionMapping
226 : * Note both BoxArray and DistributionMapping are reference counted.
227 : * Objects with the same references have the same key.
228 : */
229 : struct BDKey {
230 403505 : BDKey () noexcept = default;
231 0 : BDKey (const BoxArray::RefID& baid, const DistributionMapping::RefID& dmid) noexcept
232 0 : : 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 0 : bool operator!= (const BDKey& rhs) const noexcept {
241 0 : 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:
245 : BoxArray::RefID m_ba_id;
246 : DistributionMapping::RefID m_dm_id;
247 : };
248 :
249 0 : [[nodiscard]] BDKey getBDKey () const noexcept {
250 0 : return {boxarray.getRefID(), distributionMap.getRefID()};
251 : }
252 :
253 : void updateBDKey ();
254 :
255 : //
256 : //! Tiling
257 : struct TileArray
258 : {
259 : Long nuse{-1};
260 : Vector<int> numLocalTiles;
261 : Vector<int> indexMap;
262 : Vector<int> localIndexMap;
263 : Vector<int> localTileIndexMap;
264 : Vector<Box> tileArray;
265 : [[nodiscard]] Long bytes () const;
266 : };
267 :
268 : //
269 : //! Used for collecting information used in communicating FABs.
270 : struct FabComTag
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};
280 : int fabArrayId{0};
281 : int fillBoxId{0};
282 : int procThatNeedsData{0};
283 : int procThatHasData{0};
284 : Box box;
285 : };
286 :
287 : //! Default tilesize in MFIter
288 : static AMREX_EXPORT IntVect mfiter_tile_size;
289 :
290 : //! The maximum number of components to copy() at a time.
291 : static AMREX_EXPORT int MaxComp;
292 :
293 : //! Initialize from ParmParse with "fabarray" prefix.
294 : static void Initialize ();
295 : static void Finalize ();
296 : /**
297 : * To maximize thread efficiency we now can decompose things like
298 : * intersections among boxes into smaller tiles. This sets
299 : * their maximum size.
300 : */
301 : static AMREX_EXPORT IntVect comm_tile_size; //!< communication tile size
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 :
316 : BoxArray ba_crse_patch;
317 : BoxArray ba_fine_patch;
318 : DistributionMapping dm_patch;
319 : std::unique_ptr<FabFactory<FArrayBox> > fact_crse_patch;
320 : std::unique_ptr<FabFactory<FArrayBox> > fact_fine_patch;
321 : //
322 : BDKey m_srcbdk;
323 : BDKey m_dstbdk;
324 : Box m_dstdomain;
325 : IntVect m_dstng;
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 :
334 : static FPinfoCache m_TheFillPatchCache;
335 :
336 : static CacheStats m_FPinfo_stats;
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 : //
349 : //! coarse/fine boundary
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 :
363 : BoxArray ba_cfb;
364 : DistributionMapping dm_cfb;
365 : Vector<int> fine_grid_idx; //!< local array
366 : //
367 : BDKey m_fine_bdk;
368 : Box m_fine_domain;
369 : IntVect m_ng;
370 : bool m_include_periodic;
371 : bool m_include_physbndry;
372 : //
373 : Long m_nuse{0};
374 : };
375 :
376 : using CFinfoCache = std::multimap<BDKey,FabArrayBase::CFinfo*>;
377 : using CFinfoCacheIter = CFinfoCache::iterator;
378 :
379 : static CFinfoCache m_TheCrseFineCache;
380 :
381 : static CacheStats m_CFinfo_stats;
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 : //
392 : //! parallel copy or add
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;
420 : RegionTag& operator= (RegionTag &&) = 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 :
432 : /**
433 : * \brief Return owenership of fabs. The concept of ownership only applies when UPC++
434 : * team is used. In that case, each fab is shared by team workers, with one
435 : * taking the ownership.
436 : */
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 : //
443 : mutable BoxArray boxarray;
444 : DistributionMapping distributionMap;
445 : Vector<int> indexArray;
446 : std::vector<bool> ownership;
447 : IntVect n_grow;
448 : int n_comp;
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 : //
461 : static TACache m_TheTileArrayCache;
462 : static CacheStats m_TAC_stats;
463 : //
464 : void buildTileArray (const IntVect& tilesize, TileArray& ta) const;
465 : //
466 : void flushTileArray (const IntVect& tilesize = IntVect::TheZeroVector(),
467 : bool no_assertion=false) const;
468 : static void flushTileArrayCache (); //!< This flushes the entire cache.
469 :
470 : struct CommMetaData
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 : //
484 : //! FillBoundary
485 : struct FB
486 : : CommMetaData
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 :
493 : IndexType m_typ;
494 : IntVect m_crse_ratio; //!< BoxArray in FabArrayBase may have crse_ratio.
495 : IntVect m_ngrow;
496 : bool m_cross;
497 : bool m_epo;
498 : bool m_override_sync;
499 : Periodicity m_period;
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 : //
522 : static FBCache m_TheFBCache;
523 : static CacheStats m_FBC_stats;
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; //!< This flushes its own FB.
530 : static void flushFBCache (); //!< This flushes the entire cache.
531 :
532 : //
533 : //! parallel copy or add
534 : struct CPC
535 : : CommMetaData
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 :
550 : BDKey m_srcbdk;
551 : BDKey m_dstbdk;
552 : IntVect m_srcng;
553 : IntVect m_dstng;
554 : Periodicity m_period;
555 : bool m_tgco;
556 : BoxArray m_srcba;
557 : BoxArray m_dstba;
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 : //
573 : static CPCache m_TheCPCache;
574 : static CacheStats m_CPC_stats;
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; //!< This flushes its own CPC.
580 : static void flushCPCache (); //!< This flusheds the entire cache.
581 :
582 : //
583 : //! Rotate Boundary by 90
584 : struct RB90
585 : : CommMetaData
586 : {
587 : RB90 (const FabArrayBase& fa, const IntVect& nghost, Box const& domain);
588 : IntVect m_ngrow;
589 : Box m_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 : //
597 : static RB90Cache m_TheRB90Cache;
598 : //
599 : const RB90& getRB90 (const IntVect& nghost, const Box& domain) const;
600 : //
601 : void flushRB90 (bool no_assertion=false) const; //!< This flushes its own RB90.
602 : static void flushRB90Cache (); //!< This flushes the entire cache.
603 :
604 : //
605 : //! Rotate Boundary by 180
606 : struct RB180
607 : : CommMetaData
608 : {
609 : RB180 (const FabArrayBase& fa, const IntVect& nghost, Box const& domain);
610 : IntVect m_ngrow;
611 : Box m_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 : //
619 : static RB180Cache m_TheRB180Cache;
620 : //
621 : const RB180& getRB180 (const IntVect& nghost, const Box& domain) const;
622 : //
623 : void flushRB180 (bool no_assertion=false) const; //!< This flushes its own RB180.
624 : static void flushRB180Cache (); //!< This flushes the entire cache.
625 :
626 : //
627 : //! Fill polar boundary in spherical coordinates.
628 : struct PolarB
629 : : CommMetaData
630 : {
631 : PolarB (const FabArrayBase& fa, const IntVect& nghost, Box const& domain);
632 : IntVect m_ngrow;
633 : Box m_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 : //
641 : static PolarBCache m_ThePolarBCache;
642 : //
643 : const PolarB& getPolarB (const IntVect& nghost, const Box& domain) const;
644 : //
645 : void flushPolarB (bool no_assertion=false) const; //!< This flushes its own PolarB.
646 : static void flushPolarBCache (); //!< This flushes the entire cache.
647 :
648 : #ifdef AMREX_USE_GPU
649 : //
650 : //! For ParallelFor(FabArray)
651 : struct ParForInfo
652 : {
653 : ParForInfo (const FabArrayBase& fa, const IntVect& nghost, int nthreads);
654 : ~ParForInfo ();
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 :
665 : BATransformer m_bat;
666 : IntVect m_ng;
667 : int m_nthreads;
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 : //
684 : //! Keep track of how many FabArrays are built with the same BDKey.
685 : static std::map<BDKey, int> m_BD_count;
686 : //
687 : //! clear BD count and caches associated with this BD, if no other is using this BD.
688 : void clearThisBD (bool no_assertion=false) const;
689 : //
690 : //! add the current BD into BD count database
691 : void addThisBD ();
692 : //
693 : struct FabArrayStats
694 : {
695 : int num_fabarrays{0};
696 : int max_num_fabarrays{0};
697 : int max_num_boxarrays{0};
698 : int max_num_ba_use{1};
699 : Long num_build{0};
700 :
701 403505 : void recordBuild () noexcept {
702 403505 : ++num_fabarrays;
703 403505 : ++num_build;
704 403505 : max_num_fabarrays = std::max(max_num_fabarrays, num_fabarrays);
705 403505 : }
706 452872 : void recordDelete () noexcept {
707 452872 : --num_fabarrays;
708 452872 : }
709 : void recordMaxNumBoxArrays (int n) noexcept {
710 : max_num_boxarrays = std::max(max_num_boxarrays, n);
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 : };
723 : static AMREX_EXPORT FabArrayStats m_FA_stats;
724 :
725 : static AMREX_EXPORT bool m_alloc_single_chunk;
726 :
727 63839 : [[nodiscard]] static bool getAllocSingleChunk () { return m_alloc_single_chunk; }
728 : };
729 :
730 : namespace detail {
731 : class SingleChunkArena final
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;
740 : SingleChunkArena (SingleChunkArena&& rhs) = delete;
741 : SingleChunkArena& operator= (const SingleChunkArena& rhs) = delete;
742 : SingleChunkArena& operator= (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:
758 : DataAllocator m_dallocator;
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
771 : bool CheckRcvStats (Vector<MPI_Status>& recv_stats, const Vector<std::size_t>& recv_size, int tag);
772 : #endif
773 :
774 : std::ostream& operator<< (std::ostream& os, const FabArrayBase::BDKey& id);
775 :
776 : }
777 :
778 : #endif
|