Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_BoxArray.H
Go to the documentation of this file.
1
2#ifndef BL_BOXARRAY_H
3#define BL_BOXARRAY_H
4#include <AMReX_Config.H>
5
6#include <AMReX_IndexType.H>
7#include <AMReX_BoxList.H>
8#include <AMReX_Array.H>
9#include <AMReX_Periodicity.H>
10#include <AMReX_Vector.H>
11
12#include <iosfwd>
13#include <atomic>
14#include <cstddef>
15#include <cstdint>
16#include <map>
17#include <memory>
18#include <unordered_map>
19
20namespace amrex
21{
22 class BoxArray;
23
25 [[nodiscard]] BoxArray boxComplement (const Box& b1in, const Box& b2);
26
28 [[nodiscard]] BoxArray complementIn (const Box& b, const BoxArray& ba);
29
31 [[nodiscard]] BoxArray intersect (const BoxArray& ba, const Box& b, int ng = 0);
32
33 [[nodiscard]] BoxArray intersect (const BoxArray& ba, const Box& b, const IntVect& ng);
34
36 [[nodiscard]] BoxArray intersect (const BoxArray& lhs, const BoxArray& rhs);
37
39 [[nodiscard]] BoxList intersect (const BoxArray& ba, const BoxList& bl);
40
41 [[nodiscard]] BoxArray convert (const BoxArray& ba, IndexType typ);
42 [[nodiscard]] BoxArray convert (const BoxArray& ba, const IntVect& typ);
43
44 [[nodiscard]] BoxArray coarsen (const BoxArray& ba, int ratio);
45 [[nodiscard]] BoxArray coarsen (const BoxArray& ba, const IntVect& ratio);
46
47 [[nodiscard]] BoxArray refine (const BoxArray& ba, int ratio);
48 [[nodiscard]] BoxArray refine (const BoxArray& ba, const IntVect& ratio);
49
51 [[nodiscard]] BoxList GetBndryCells (const BoxArray& ba, int ngrow);
52
54 void readBoxArray (BoxArray& ba, std::istream& s, bool b = false);
55
57 [[nodiscard]] bool match (const BoxArray& x, const BoxArray& y);
58
74 [[nodiscard]] BoxArray decompose (Box const& domain, int nboxes,
75 Array<bool,AMREX_SPACEDIM> const& decomp
76 = {AMREX_D_DECL(true,true,true)},
77 bool no_overlap = false);
78
80struct BARef
81{
82 BARef ();
83 explicit BARef (size_t size);
84 explicit BARef (const Box& b);
85 explicit BARef (const BoxList& bl);
86 explicit BARef (BoxList&& bl) noexcept;
87 explicit BARef (std::istream& is);
88 BARef (const BARef& rhs);
89 BARef (BARef&& rhs) = delete;
90 BARef& operator= (const BARef& rhs) = delete;
91 BARef& operator= (BARef&& rhs) = delete;
92
93 ~BARef ();
94
95 void define (const Box& bx);
96 void define (const BoxList& bl);
97 void define (BoxList&& bl) noexcept;
98 void define (std::istream& is, int& ndims);
100 void resize (Long n);
101#ifdef AMREX_MEM_PROFILING
102 void updateMemoryUsage_box (int s);
103 void updateMemoryUsage_hash (int s);
104#endif
105
106 [[nodiscard]] bool HasHashMap () const noexcept {
107 return has_hashmap.load(std::memory_order_acquire);
108 }
109
110 //
112 Vector<Box> m_abox;
113 //
115 mutable Box bbox;
116
117 mutable IntVect crsn;
118
119 using HashType = std::unordered_map< IntVect, std::vector<int>, IntVect::shift_hasher > ;
120 //using HashType = std::map< IntVect,std::vector<int> >;
121
122 mutable HashType hash;
123
124 mutable std::atomic<bool> has_hashmap = false;
125
126 static int numboxarrays;
127 static int numboxarrays_hwm;
128 static Long total_box_bytes;
129 static Long total_box_bytes_hwm;
130 static Long total_hash_bytes;
131 static Long total_hash_bytes_hwm;
132
133 static void Initialize ();
134 static void Finalize ();
135 static bool initialized;
136};
138
140struct BATnull
141{
142 [[nodiscard]] Box operator() (const Box& bx) const noexcept { return bx; }
143 [[nodiscard]] static constexpr Box coarsen (Box const& a_box) { return a_box; }
144 [[nodiscard]] static constexpr IntVect doiLo () { return IntVect::TheZeroVector(); }
145 [[nodiscard]] static constexpr IntVect doiHi () { return IntVect::TheZeroVector(); }
146 [[nodiscard]] static constexpr IndexType index_type () { return IndexType(); }
147 [[nodiscard]] static constexpr IntVect coarsen_ratio () { return IntVect::TheUnitVector(); }
148};
150
152struct BATindexType
153{
154 explicit BATindexType (IndexType a_typ) : m_typ(a_typ) {}
155 [[nodiscard]] Box operator() (const Box& bx) const noexcept { return amrex::convert(bx,m_typ); }
156 [[nodiscard]] static Box coarsen (Box const& a_box) noexcept { return a_box; }
157 [[nodiscard]] static constexpr IntVect doiLo () { return IntVect::TheZeroVector(); }
158 [[nodiscard]] IntVect doiHi () const noexcept { return m_typ.ixType(); }
159 [[nodiscard]] IndexType index_type () const noexcept { return m_typ; }
160 [[nodiscard]] static constexpr IntVect coarsen_ratio () { return IntVect::TheUnitVector(); }
161 IndexType m_typ;
162};
164
166struct BATcoarsenRatio
167{
168 explicit BATcoarsenRatio (IntVect const& a_crse_ratio) : m_crse_ratio(a_crse_ratio) {}
169 [[nodiscard]] Box operator() (const Box& bx) const noexcept { return amrex::coarsen(bx,m_crse_ratio); }
170 [[nodiscard]] Box coarsen (Box const& a_box) const noexcept { return amrex::coarsen(a_box,m_crse_ratio); }
171 [[nodiscard]] static constexpr IntVect doiLo () { return IntVect::TheZeroVector(); }
172 [[nodiscard]] static constexpr IntVect doiHi () { return IntVect::TheZeroVector(); }
173 [[nodiscard]] static constexpr IndexType index_type () { return IndexType(); }
174 [[nodiscard]] IntVect coarsen_ratio () const noexcept { return m_crse_ratio; }
175 IntVect m_crse_ratio;
176};
178
180struct BATindexType_coarsenRatio
181{
182 BATindexType_coarsenRatio (IndexType a_typ, IntVect const& a_crse_ratio)
183 : m_typ(a_typ), m_crse_ratio(a_crse_ratio) {}
184
185 [[nodiscard]] Box operator() (const Box& bx) const noexcept {
186 return amrex::convert(amrex::coarsen(bx,m_crse_ratio),m_typ);
187 }
188
189 [[nodiscard]] Box coarsen (Box const& a_box) const noexcept { return amrex::coarsen(a_box,m_crse_ratio); }
190
191 [[nodiscard]] static constexpr IntVect doiLo () { return IntVect::TheZeroVector(); }
192 [[nodiscard]] IntVect doiHi () const noexcept { return m_typ.ixType(); }
193
194 [[nodiscard]] IndexType index_type () const noexcept { return m_typ; }
195 [[nodiscard]] IntVect coarsen_ratio () const noexcept { return m_crse_ratio; }
196
197 IndexType m_typ;
198 IntVect m_crse_ratio;
199};
201
203struct BATbndryReg
204{
205 BATbndryReg (Orientation a_face, IndexType a_typ,
206 int a_in_rad, int a_out_rad, int a_extent_rad)
207 : m_face(a_face), m_typ(a_typ), m_crse_ratio(1)
208 {
209 m_loshft = IntVect(-a_extent_rad);
210 m_hishft = IntVect( a_extent_rad);
211 IntVect nodal = a_typ.ixType();
212 m_hishft += nodal;
213 const int d = a_face.coordDir();
214 if (nodal[d]) {
215 // InterpFaceRegister & SyncRegister in IAMR
216 if (m_face.isLow()) {
217 m_loshft[d] = 0;
218 m_hishft[d] = 0;
219 } else {
220 m_loshft[d] = 1;
221 m_hishft[d] = 1;
222 }
223 m_doilo = IntVect(0);
224 m_doihi = nodal;
225 } else {
226 // BndryRegister
227 if (m_face.isLow()) {
228 m_loshft[d] = nodal[d] - a_out_rad;
229 m_hishft[d] = nodal[d] + a_in_rad - 1;
230 } else {
231 m_loshft[d] = 1 - a_in_rad;
232 m_hishft[d] = a_out_rad;
233 }
234 m_doilo = IntVect(a_extent_rad);
235 m_doihi = IntVect(a_extent_rad);
236 m_doihi += nodal;
237 if (m_face.isLow()) { // domain of influence in index space
238 m_doilo[d] = a_out_rad;
239 m_doihi[d] = 0;
240 } else {
241 m_doilo[d] = 0;
242 m_doihi[d] = a_out_rad;
243 }
244 }
245 }
246
247 [[nodiscard]] Box operator() (const Box& a_bx) const noexcept {
248 IntVect lo = amrex::coarsen(a_bx.smallEnd(), m_crse_ratio);
249 IntVect hi = amrex::coarsen(a_bx.bigEnd(), m_crse_ratio);
250 const int d = m_face.coordDir();
251 if (m_face.isLow()) {
252 hi[d] = lo[d];
253 } else {
254 lo[d] = hi[d];
255 }
256 lo += m_loshft;
257 hi += m_hishft;
258 return Box(lo,hi,m_typ);
259 }
260
261 [[nodiscard]] Box coarsen (Box const& a_box) const noexcept { return amrex::coarsen(a_box,m_crse_ratio); }
262
263 [[nodiscard]] IntVect doiLo () const noexcept { return m_doilo; }
264 [[nodiscard]] IntVect doiHi () const noexcept { return m_doihi; }
265
266 [[nodiscard]] IndexType index_type () const noexcept { return m_typ; }
267 [[nodiscard]] IntVect coarsen_ratio () const noexcept { return m_crse_ratio; }
268
269 friend bool operator== (BATbndryReg const& a, BATbndryReg const& b) noexcept {
270 return a.m_face == b.m_face && a.m_typ == b.m_typ && a.m_crse_ratio == b.m_crse_ratio
271 && a.m_loshft == b.m_loshft && a.m_hishft == b.m_hishft
272 && a.m_doilo == b.m_doilo && a.m_doihi == b.m_doihi;
273 }
274
275 Orientation m_face;
276 IndexType m_typ;
277 IntVect m_crse_ratio;
278 IntVect m_loshft;
279 IntVect m_hishft;
280 IntVect m_doilo;
281 IntVect m_doihi;
282};
284
286struct BATransformer
287{
288 enum struct BATType : std::int8_t { null, indexType, coarsenRatio, indexType_coarsenRatio, bndryReg};
289
290 union BATOp {
291 BATOp () noexcept
292 : m_null() {}
293 BATOp (IndexType t) noexcept
294 : m_indexType(t) {}
295 BATOp (IntVect const& r) noexcept
296 : m_coarsenRatio(r) {}
297 BATOp (IndexType t, IntVect const& r) noexcept
298 : m_indexType_coarsenRatio(t,r) {}
299 BATOp (Orientation f, IndexType t, int in_rad, int out_rad, int extent_rad) noexcept
300 : m_bndryReg(f,t,in_rad,out_rad,extent_rad) {}
301 BATnull m_null;
302 BATindexType m_indexType;
303 BATcoarsenRatio m_coarsenRatio;
304 BATindexType_coarsenRatio m_indexType_coarsenRatio;
305 BATbndryReg m_bndryReg;
306 };
307
308 BATransformer () = default;
309
310 BATransformer (IndexType t)
311 : m_bat_type(t.cellCentered() ? BATType::null : BATType::indexType),
312 m_op (t.cellCentered() ? BATOp() : BATOp(t)) {}
313
314 BATransformer (Orientation f, IndexType t, int in_rad, int out_rad, int extent_rad)
315 : m_bat_type(BATType::bndryReg),
316 m_op(f,t,in_rad,out_rad,extent_rad) {}
317
318 [[nodiscard]] Box operator() (Box const& ab) const noexcept {
319 switch (m_bat_type)
320 {
321 case BATType::null:
322 return m_op.m_null(ab);
323 case BATType::indexType:
324 return m_op.m_indexType(ab);
325 case BATType::coarsenRatio:
326 return m_op.m_coarsenRatio(ab);
327 case BATType::indexType_coarsenRatio:
328 return m_op.m_indexType_coarsenRatio(ab);
329 default:
330 return m_op.m_bndryReg(ab);
331 }
332 }
333
334 [[nodiscard]] Box coarsen (Box const& a_box) const noexcept {
335 switch (m_bat_type)
336 {
337 case BATType::null:
338 return amrex::BATnull::coarsen(a_box);
339 case BATType::indexType:
340 return amrex::BATindexType::coarsen(a_box);
341 case BATType::coarsenRatio:
342 return m_op.m_coarsenRatio.coarsen(a_box);
343 case BATType::indexType_coarsenRatio:
344 return m_op.m_indexType_coarsenRatio.coarsen(a_box);
345 default:
346 return m_op.m_bndryReg.coarsen(a_box);
347 }
348 }
349
350 [[nodiscard]] IntVect doiLo () const noexcept {
351 switch (m_bat_type)
352 {
353 case BATType::null:
354 return amrex::BATnull::doiLo();
355 case BATType::indexType:
356 return amrex::BATindexType::doiLo();
357 case BATType::coarsenRatio:
358 return amrex::BATcoarsenRatio::doiLo();
359 case BATType::indexType_coarsenRatio:
360 return amrex::BATindexType_coarsenRatio::doiLo();
361 default:
362 return m_op.m_bndryReg.doiLo();
363 }
364 }
365
366 [[nodiscard]] IntVect doiHi () const noexcept {
367 switch (m_bat_type)
368 {
369 case BATType::null:
370 return amrex::BATnull::doiHi();
371 case BATType::indexType:
372 return m_op.m_indexType.doiHi();
373 case BATType::coarsenRatio:
374 return amrex::BATcoarsenRatio::doiHi();
375 case BATType::indexType_coarsenRatio:
376 return m_op.m_indexType_coarsenRatio.doiHi();
377 default:
378 return m_op.m_bndryReg.doiHi();
379 }
380 }
381
382 [[nodiscard]] IndexType index_type () const noexcept {
383 switch (m_bat_type)
384 {
385 case BATType::null:
386 return amrex::BATnull::index_type();
387 case BATType::indexType:
388 return m_op.m_indexType.index_type();
389 case BATType::coarsenRatio:
390 return amrex::BATcoarsenRatio::index_type();
391 case BATType::indexType_coarsenRatio:
392 return m_op.m_indexType_coarsenRatio.index_type();
393 default:
394 return m_op.m_bndryReg.index_type();
395 }
396 }
397
398 [[nodiscard]] IntVect coarsen_ratio () const noexcept {
399 switch (m_bat_type)
400 {
401 case BATType::null:
402 return amrex::BATnull::coarsen_ratio();
403 case BATType::indexType:
404 return amrex::BATindexType::coarsen_ratio();
405 case BATType::coarsenRatio:
406 return m_op.m_coarsenRatio.coarsen_ratio();
407 case BATType::indexType_coarsenRatio:
408 return m_op.m_indexType_coarsenRatio.coarsen_ratio();
409 default:
410 return m_op.m_bndryReg.coarsen_ratio();
411 }
412 }
413
414 [[nodiscard]] bool is_null () const noexcept {
415 return m_bat_type == BATType::null;
416 }
417
418 [[nodiscard]] bool is_simple () const noexcept {
419 return m_bat_type != BATType::bndryReg;
420 }
421
422 void set_coarsen_ratio (IntVect const& a_ratio) noexcept {
423 switch (m_bat_type)
424 {
425 case BATType::null:
426 {
427 if (a_ratio == IntVect::TheUnitVector()) {
428 return;
429 } else {
430 m_bat_type = BATType::coarsenRatio;
431 m_op.m_coarsenRatio.m_crse_ratio = a_ratio;
432 return;
433 }
434 }
435 case BATType::indexType:
436 {
437 if (a_ratio == IntVect::TheUnitVector()) {
438 return;
439 } else {
440 m_bat_type = BATType::indexType_coarsenRatio;
441 auto t = m_op.m_indexType.m_typ;
442 m_op.m_indexType_coarsenRatio.m_typ = t;
443 m_op.m_indexType_coarsenRatio.m_crse_ratio = a_ratio;
444 return;
445 }
446 }
447 case BATType::coarsenRatio:
448 {
449 if (a_ratio == IntVect::TheUnitVector()) {
450 m_bat_type = BATType::null;
451 return;
452 } else {
453 m_op.m_coarsenRatio.m_crse_ratio = a_ratio;
454 return;
455 }
456 }
457 case BATType::indexType_coarsenRatio:
458 {
459 if (a_ratio == IntVect::TheUnitVector()) {
460 m_bat_type = BATType::indexType;
461 auto t = m_op.m_indexType_coarsenRatio.m_typ;
462 m_op.m_indexType.m_typ = t;
463 return;
464 } else {
465 m_op.m_indexType_coarsenRatio.m_crse_ratio = a_ratio;
466 return;
467 }
468 }
469 default:
470 {
471 m_op.m_bndryReg.m_crse_ratio = a_ratio;
472 return;
473 }
474 }
475 }
476
477 void set_index_type (IndexType typ) noexcept {
478 switch (m_bat_type)
479 {
480 case BATType::null:
481 {
482 if (typ.cellCentered()) {
483 return;
484 } else {
485 m_bat_type = BATType::indexType;
486 m_op.m_indexType.m_typ = typ;
487 return;
488 }
489 }
490 case BATType::indexType:
491 {
492 if (typ.cellCentered()) {
493 m_bat_type = BATType::null;
494 return;
495 } else {
496 m_op.m_indexType.m_typ = typ;
497 return;
498 }
499 }
500 case BATType::coarsenRatio:
501 {
502 if (typ.cellCentered()) {
503 return;
504 } else {
505 m_bat_type = BATType::indexType_coarsenRatio;
506 auto r = m_op.m_coarsenRatio.m_crse_ratio;
507 m_op.m_indexType_coarsenRatio.m_typ = typ;
508 m_op.m_indexType_coarsenRatio.m_crse_ratio = r;
509 return;
510 }
511 }
512 case BATType::indexType_coarsenRatio:
513 {
514 if (typ.cellCentered()) {
515 m_bat_type = BATType::coarsenRatio;
516 auto r = m_op.m_indexType_coarsenRatio.m_crse_ratio;
517 m_op.m_coarsenRatio.m_crse_ratio = r;
518 return;
519 } else {
520 m_op.m_indexType_coarsenRatio.m_typ = typ;
521 return;
522 }
523 }
524 default:
525 {
526 m_op.m_bndryReg.m_typ = typ;
527 return;
528 }
529 }
530 }
531
532 friend bool operator== (BATransformer const& a, BATransformer const& b) noexcept {
533 if (a.m_bat_type != BATType::bndryReg && b.m_bat_type != BATType::bndryReg) {
534 return a.index_type() == b.index_type()
535 && a.coarsen_ratio() == b.coarsen_ratio();
536 } else if (a.m_bat_type == BATType::bndryReg && b.m_bat_type == BATType::bndryReg) {
537 return a.m_op.m_bndryReg == b.m_op.m_bndryReg;
538 } else {
539 return false;
540 }
541 }
542
543 BATType m_bat_type{BATType::null};
544 BATOp m_op;
545};
547
548// for backward compatibility
549using BndryBATransformer = BATransformer;
550
551class MFIter;
552class AmrMesh;
553class FabArrayBase;
554
564{
565public:
566
567 BoxArray () noexcept;
568 BoxArray (const BoxArray& rhs) = default;
569 BoxArray (BoxArray&& rhs) noexcept = default;
570 BoxArray& operator= (BoxArray const& rhs) = default;
571 BoxArray& operator= (BoxArray&& rhs) noexcept = default;
572 ~BoxArray() noexcept = default;
573
575 explicit BoxArray (const Box& bx);
576
578 explicit BoxArray (size_t n);
579
581 BoxArray (const Box* bxvec,
582 int nbox);
583
585 explicit BoxArray (const BoxList& bl);
586 explicit BoxArray (BoxList&& bl) noexcept;
587
588 BoxArray (const BoxArray& rhs, const BATransformer& trans);
589
590 BoxArray (BoxList&& bl, IntVect const& max_grid_size);
591
596 void define (const Box& bx);
601 void define (const BoxList& bl);
602 void define (BoxList&& bl) noexcept;
603
605 void clear ();
606
608 void resize (Long len);
609
611 [[nodiscard]] Long size () const noexcept { return m_ref->m_abox.size(); }
612
614 [[nodiscard]] Long capacity () const noexcept { return static_cast<Long>(m_ref->m_abox.capacity()); }
615
617 [[nodiscard]] bool empty () const noexcept { return m_ref->m_abox.empty(); }
618
620 [[nodiscard]] Long numPts() const noexcept;
621
623 [[nodiscard]] double d_numPts () const noexcept;
630 int readFrom (std::istream& is);
631
633 std::ostream& writeOn (std::ostream&) const;
634
636 [[nodiscard]] bool operator== (const BoxArray& rhs) const noexcept;
637
639 [[nodiscard]] bool operator!= (const BoxArray& rhs) const noexcept;
640
641 [[nodiscard]] bool operator== (const Vector<Box>& bv) const noexcept;
642 [[nodiscard]] bool operator!= (const Vector<Box>& bv) const noexcept;
643
645 [[nodiscard]] bool CellEqual (const BoxArray& rhs) const noexcept;
646
648 BoxArray& maxSize (int block_size);
649
650 BoxArray& maxSize (const IntVect& block_size);
651
655 BoxArray& minmaxSize (const IntVect& min_size, const IntVect& max_size);
656
658 BoxArray& refine (int refinement_ratio);
659
661 BoxArray& refine (const IntVect& iv);
662
664 BoxArray& coarsen (int refinement_ratio);
665
667 [[nodiscard]] bool coarsenable (int refinement_ratio, int min_width=1) const;
668 [[nodiscard]] bool coarsenable (const IntVect& refinement_ratio, int min_width=1) const;
669 [[nodiscard]] bool coarsenable (const IntVect& refinement_ratio, const IntVect& min_width) const;
670
672 BoxArray& coarsen (const IntVect& iv);
673
675 BoxArray& growcoarsen (int n, const IntVect& iv);
676 BoxArray& growcoarsen (IntVect const& ngrow, const IntVect& iv);
677
679 BoxArray& grow (int n);
680
682 BoxArray& grow (const IntVect& iv);
687 BoxArray& grow (int idir, int n_cell);
692 BoxArray& growLo (int idir, int n_cell);
697 BoxArray& growHi (int idir, int n_cell);
707 BoxArray& surroundingNodes (int dir);
708
711
713 BoxArray& enclosedCells (int dir);
714
717
718 BoxArray& convert (const IntVect& iv);
719
721 BoxArray& convert (Box (*fp)(const Box&));
722
724 BoxArray& shift (int dir, int nzones);
725
727 BoxArray& shift (const IntVect &iv);
728
730 void set (int i, const Box& ibox);
731
733 [[nodiscard]] Box operator[] (int index) const noexcept {
734 return m_bat(m_ref->m_abox[index]);
735 }
736
738 [[nodiscard]] Box operator[] (const MFIter& mfi) const noexcept;
739
741 [[nodiscard]] Box get (int index) const noexcept { return operator[](index); }
742
744 [[nodiscard]] Box getCellCenteredBox (int index) const noexcept {
745 return m_bat.coarsen(m_ref->m_abox[index]);
746 }
747
752 [[nodiscard]] bool ok () const;
753
755 [[nodiscard]] bool isDisjoint () const;
756
758 [[nodiscard]] BoxList boxList () const;
759
761 [[nodiscard]] bool contains (const IntVect& v) const;
762
767 [[nodiscard]]
768 bool contains (const Box& b, bool assume_disjoint_ba = false,
769 const IntVect& ng = IntVect(0)) const;
770
774 [[nodiscard]]
775 bool contains (const BoxArray& ba, bool assume_disjoint_ba = false,
776 const IntVect& ng = IntVect(0)) const;
777
785 [[nodiscard]]
786 bool contains (const BoxArray& ba, Periodicity const& period) const;
787
789 [[nodiscard]] Box minimalBox () const;
790 [[nodiscard]] Box minimalBox (Long& npts_avg_box) const;
791
796 [[nodiscard]]
797 bool intersects (const Box& b, int ng = 0) const;
798
799 [[nodiscard]]
800 bool intersects (const Box& b, const IntVect& ng) const;
801
803 [[nodiscard]]
804 std::vector< std::pair<int,Box> > intersections (const Box& bx) const;
805
807 [[nodiscard]]
808 std::vector< std::pair<int,Box> > intersections (const Box& bx, bool first_only, int ng) const;
809
810 [[nodiscard]]
811 std::vector< std::pair<int,Box> > intersections (const Box& bx, bool first_only, const IntVect& ng) const;
812
814 void intersections (const Box& bx, std::vector< std::pair<int,Box> >& isects) const;
815
817 void intersections (const Box& bx, std::vector< std::pair<int,Box> >& isects,
818 bool first_only, int ng) const;
819
820 void intersections (const Box& bx, std::vector< std::pair<int,Box> >& isects,
821 bool first_only, const IntVect& ng) const;
822
824 [[nodiscard]] BoxList complementIn (const Box& b) const;
825 void complementIn (BoxList& bl, const Box& b) const;
826
828 [[nodiscard]] BoxList complementIn (const Box& b, const Periodicity& period) const;
829
831 void clear_hash_bin () const;
832
834 void removeOverlap (bool simplify=true);
835
837 [[nodiscard]] static bool SameRefs (const BoxArray& lhs, const BoxArray& rhs) { return lhs.m_ref == rhs.m_ref; }
838
839 struct RefID {
840 RefID () noexcept = default;
841 explicit RefID (BARef* data_) noexcept : data(data_) {}
842 bool operator< (const RefID& rhs) const noexcept { return std::less<>()(data,rhs.data); }
843 bool operator== (const RefID& rhs) const noexcept { return data == rhs.data; }
844 bool operator!= (const RefID& rhs) const noexcept { return data != rhs.data; }
845 friend std::ostream& operator<< (std::ostream& os, const RefID& id);
846 private:
847 BARef* data{nullptr};
848 };
849
851 [[nodiscard]] RefID getRefID () const noexcept { return RefID { m_ref.get() }; }
852
854 [[nodiscard]] IndexType ixType () const noexcept { return m_bat.index_type(); }
855
857 [[nodiscard]] IntVect crseRatio () const noexcept { return m_bat.coarsen_ratio(); }
858
859 static void Initialize ();
860 static void Finalize ();
861 static bool initialized;
862
864 void uniqify ();
865
866 [[nodiscard]] BoxList const& simplified_list () const; // For regular AMR grids only!!!
867 [[nodiscard]] BoxArray simplified () const; // For regular AMR grids only!!!
868
869 [[nodiscard]] BATransformer const& transformer () const;
870
871 [[nodiscard]] std::weak_ptr<BARef> getWeakRef () const;
872 [[nodiscard]] std::shared_ptr<BARef> const& getSharedRef () const;
873 std::shared_ptr<BARef>& getSharedRef ();
874
875 friend class AmrMesh;
876 friend class FabArrayBase;
877
878private:
880 void type_update ();
881
882 [[nodiscard]] BARef::HashType& getHashMap () const;
883
884 [[nodiscard]] IntVect getDoiLo () const noexcept;
885 [[nodiscard]] IntVect getDoiHi () const noexcept;
886
887 BATransformer m_bat;
889 std::shared_ptr<BARef> m_ref;
890 mutable std::shared_ptr<BoxList> m_simplified_list;
891};
892
894std::ostream& operator<< (std::ostream& os, const BoxArray& ba);
895
896std::ostream& operator<< (std::ostream& os, const BoxArray::RefID& id);
897
898}
899
900#endif /*BL_BOXARRAY_H*/
int idir
Definition AMReX_HypreMLABecLap.cpp:1143
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
Definition AMReX_AmrMesh.H:69
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:564
IndexType ixType() const noexcept
Return index type of this BoxArray.
Definition AMReX_BoxArray.H:854
Box operator[](int index) const noexcept
Return element index of this BoxArray.
Definition AMReX_BoxArray.H:733
Long capacity() const noexcept
Return the number of boxes that can be held in the current allocated storage.
Definition AMReX_BoxArray.H:614
RefID getRefID() const noexcept
Return a unique ID of the reference.
Definition AMReX_BoxArray.H:851
bool intersects(const Box &b, int ng=0) const
True if the Box intersects with this BoxArray(+ghostcells). The Box must have the same IndexType as t...
Definition AMReX_BoxArray.cpp:1169
BoxArray & grow(int n)
Grow each Box in the BoxArray by the specified amount.
Definition AMReX_BoxArray.cpp:706
static void Finalize()
Definition AMReX_BoxArray.cpp:282
BoxArray() noexcept
Definition AMReX_BoxArray.cpp:287
std::vector< std::pair< int, Box > > intersections(const Box &bx) const
Return intersections of Box and BoxArray.
Definition AMReX_BoxArray.cpp:1186
std::ostream & writeOn(std::ostream &) const
Output this BoxArray to a checkpoint file.
Definition AMReX_BoxArray.cpp:482
bool contains(const IntVect &v) const
True if the IntVect is within any of the Boxes in this BoxArray.
Definition AMReX_BoxArray.cpp:977
static void Initialize()
Definition AMReX_BoxArray.cpp:271
BoxList const & simplified_list() const
Definition AMReX_BoxArray.cpp:1636
BoxList boxList() const
Create a BoxList from this BoxArray.
Definition AMReX_BoxArray.cpp:949
BoxArray simplified() const
Definition AMReX_BoxArray.cpp:1647
IntVect crseRatio() const noexcept
Return crse ratio of this BoxArray.
Definition AMReX_BoxArray.H:857
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
void define(const Box &bx)
Initialize the BoxArray from a single box. It is an error if the BoxArray has already been initialize...
Definition AMReX_BoxArray.cpp:352
Long numPts() const noexcept
Returns the total number of cells contained in all boxes in the BoxArray.
Definition AMReX_BoxArray.cpp:394
void clear_hash_bin() const
Clear out the internal hash table used by intersections.
Definition AMReX_BoxArray.cpp:1440
Box minimalBox() const
Return smallest Box that contains all Boxes in this BoxArray.
Definition AMReX_BoxArray.cpp:1067
BoxArray & minmaxSize(const IntVect &min_size, const IntVect &max_size)
Definition AMReX_BoxArray.cpp:577
BoxArray & refine(int refinement_ratio)
Refine each Box in the BoxArray to the specified ratio.
Definition AMReX_BoxArray.cpp:593
std::weak_ptr< BARef > getWeakRef() const
Definition AMReX_BoxArray.cpp:1659
bool coarsenable(int refinement_ratio, int min_width=1) const
Coarsen each Box in the BoxArray to the specified ratio.
Definition AMReX_BoxArray.cpp:615
BoxArray & growLo(int idir, int n_cell)
Grow each Box in the BoxArray on the low end by n_cell cells in the idir direction.
Definition AMReX_BoxArray.cpp:752
void resize(Long len)
Resize the BoxArray. See Vector<T>::resize() for the gory details.
Definition AMReX_BoxArray.cpp:387
BoxArray & enclosedCells()
Apply Box::enclosedCells() to each Box in the BoxArray.
Definition AMReX_BoxArray.cpp:799
std::shared_ptr< BARef > const & getSharedRef() const
Definition AMReX_BoxArray.cpp:1665
Box getCellCenteredBox(int index) const noexcept
Return cell-centered box at element index of this BoxArray.
Definition AMReX_BoxArray.H:744
double d_numPts() const noexcept
Returns the total number of cells (in double type) contained in all boxes in the BoxArray.
Definition AMReX_BoxArray.cpp:431
BATransformer const & transformer() const
Definition AMReX_BoxArray.cpp:1653
static bool initialized
Definition AMReX_BoxArray.H:861
static bool SameRefs(const BoxArray &lhs, const BoxArray &rhs)
whether two BoxArrays share the same data
Definition AMReX_BoxArray.H:837
BoxArray & growHi(int idir, int n_cell)
Grow each Box in the BoxArray on the high end by n_cell cells in the idir direction.
Definition AMReX_BoxArray.cpp:768
BoxList complementIn(const Box &b) const
Return box - boxarray.
Definition AMReX_BoxArray.cpp:1314
void removeOverlap(bool simplify=true)
Change the BoxArray to one with no overlap and then simplify it (see the simplify function in BoxList...
Definition AMReX_BoxArray.cpp:1456
BoxArray & coarsen(int refinement_ratio)
Coarsen each Box in the BoxArray to the specified ratio.
Definition AMReX_BoxArray.cpp:672
BoxArray & shift(int dir, int nzones)
Apply Box::shift(int,int) to each Box in the BoxArray.
Definition AMReX_BoxArray.cpp:847
BoxArray & surroundingNodes()
Apply surroundingNodes(Box) to each Box in BoxArray. See the documentation of Box for details.
Definition AMReX_BoxArray.cpp:784
BoxArray(BoxArray &&rhs) noexcept=default
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition AMReX_BoxArray.H:611
Box get(int index) const noexcept
Return element index of this BoxArray.
Definition AMReX_BoxArray.H:741
~BoxArray() noexcept=default
bool CellEqual(const BoxArray &rhs) const noexcept
Are the BoxArrays equal after conversion to cell-centered.
Definition AMReX_BoxArray.cpp:546
void set(int i, const Box &ibox)
Set element i in this BoxArray to Box ibox.
Definition AMReX_BoxArray.cpp:878
BoxArray & growcoarsen(int n, const IntVect &iv)
Grow and then coarsen each Box in the BoxArray.
Definition AMReX_BoxArray.cpp:685
BoxArray & convert(IndexType typ)
Apply Box::convert(IndexType) to each Box in the BoxArray.
Definition AMReX_BoxArray.cpp:813
bool ok() const
Return true if Box is valid and they all have the same IndexType. Is true by default if the BoxArray ...
Definition AMReX_BoxArray.cpp:894
BoxArray & operator=(BoxArray const &rhs)=default
void uniqify()
Make ourselves unique.
Definition AMReX_BoxArray.cpp:1613
BoxArray(const BoxArray &rhs)=default
bool empty() const noexcept
Return whether the BoxArray is empty.
Definition AMReX_BoxArray.H:617
bool isDisjoint() const
Return true if set of intersecting Boxes in BoxArray is null.
Definition AMReX_BoxArray.cpp:920
BoxArray & maxSize(int block_size)
Forces each Box in BoxArray to have sides <= block_size.
Definition AMReX_BoxArray.cpp:553
void clear()
Remove all Boxes from the BoxArray.
Definition AMReX_BoxArray.cpp:379
A class for managing a List of Boxes that share a common IndexType. This class implements operations ...
Definition AMReX_BoxList.H:52
Base class for FabArray.
Definition AMReX_FabArrayBase.H:42
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
This provides length of period for periodic domains. 0 means it is not periodic in that direction....
Definition AMReX_Periodicity.H:17
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
amrex_long Long
Definition AMReX_INT.H:30
__host__ __device__ BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) coarsening ratio.
Definition AMReX_Box.H:1409
__host__ __device__ BoxND< dim > refine(const BoxND< dim > &b, int ref_ratio) noexcept
Definition AMReX_Box.H:1459
Definition AMReX_Amr.cpp:50
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1558
BoxArray intersect(const BoxArray &ba, const Box &b, int ng)
Make a BoxArray from the intersection of Box b and BoxArray(+ghostcells).
Definition AMReX_BoxArray.cpp:1717
AMReX * Initialize(MPI_Comm mpi_comm, std::ostream &a_osout=std::cout, std::ostream &a_oserr=std::cerr, ErrorHandler a_errhandler=nullptr, int a_device_id=-1)
Definition AMReX.cpp:343
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
BoxList GetBndryCells(const BoxArray &ba, int ngrow)
Find the ghost cells of a given BoxArray.
Definition AMReX_BoxArray.cpp:1842
IndexTypeND< 3 > IndexType
IndexType is an alias for amrex::IndexTypeND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:36
void Finalize(AMReX *pamrex)
Definition AMReX.cpp:824
BoxArray decompose(Box const &domain, int nboxes, Array< bool, 3 > const &decomp, bool no_overlap)
Decompose domain box into BoxArray.
Definition AMReX_BoxArray.cpp:1943
bool match(const BoxArray &x, const BoxArray &y)
Note that two BoxArrays that match are not necessarily equal.
Definition AMReX_BoxArray.cpp:1930
BoxArray complementIn(const Box &b, const BoxArray &ba)
Make a BoxArray from the complement of BoxArray ba in Box b.
Definition AMReX_BoxArray.cpp:1710
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
void readBoxArray(BoxArray &ba, std::istream &is, bool bReadSpecial)
Read a BoxArray from a stream. If b is true, read in a special way.
Definition AMReX_BoxArray.cpp:1900
BATransformer BndryBATransformer
Definition AMReX_BoxArray.H:549
BoxArray boxComplement(const Box &b1in, const Box &b2)
Make a BoxArray from the the complement of b2 in b1in.
Definition AMReX_BoxArray.cpp:1703
Definition AMReX_BoxArray.H:839
friend std::ostream & operator<<(std::ostream &os, const RefID &id)
Definition AMReX_BoxArray.cpp:2126
RefID() noexcept=default
bool operator<(const RefID &rhs) const noexcept
Definition AMReX_BoxArray.H:842
bool operator!=(const RefID &rhs) const noexcept
Definition AMReX_BoxArray.H:844
bool operator==(const RefID &rhs) const noexcept
Definition AMReX_BoxArray.H:843