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 <cstddef>
14#include <map>
15#include <memory>
16#include <unordered_map>
17
18namespace amrex
19{
20 class BoxArray;
21
23 [[nodiscard]] BoxArray boxComplement (const Box& b1in, const Box& b2);
24
26 [[nodiscard]] BoxArray complementIn (const Box& b, const BoxArray& ba);
27
29 [[nodiscard]] BoxArray intersect (const BoxArray& ba, const Box& b, int ng = 0);
30
31 [[nodiscard]] BoxArray intersect (const BoxArray& ba, const Box& b, const IntVect& ng);
32
34 [[nodiscard]] BoxArray intersect (const BoxArray& lhs, const BoxArray& rhs);
35
37 [[nodiscard]] BoxList intersect (const BoxArray& ba, const BoxList& bl);
38
39 [[nodiscard]] BoxArray convert (const BoxArray& ba, IndexType typ);
40 [[nodiscard]] BoxArray convert (const BoxArray& ba, const IntVect& typ);
41
42 [[nodiscard]] BoxArray coarsen (const BoxArray& ba, int ratio);
43 [[nodiscard]] BoxArray coarsen (const BoxArray& ba, const IntVect& ratio);
44
45 [[nodiscard]] BoxArray refine (const BoxArray& ba, int ratio);
46 [[nodiscard]] BoxArray refine (const BoxArray& ba, const IntVect& ratio);
47
49 [[nodiscard]] BoxList GetBndryCells (const BoxArray& ba, int ngrow);
50
52 void readBoxArray (BoxArray& ba, std::istream& s, bool b = false);
53
55 [[nodiscard]] bool match (const BoxArray& x, const BoxArray& y);
56
72 [[nodiscard]] BoxArray decompose (Box const& domain, int nboxes,
73 Array<bool,AMREX_SPACEDIM> const& decomp
74 = {AMREX_D_DECL(true,true,true)},
75 bool no_overlap = false);
76
78struct BARef
79{
80 BARef ();
81 explicit BARef (size_t size);
82 explicit BARef (const Box& b);
83 explicit BARef (const BoxList& bl);
84 explicit BARef (BoxList&& bl) noexcept;
85 explicit BARef (std::istream& is);
86 BARef (const BARef& rhs);
87 BARef (BARef&& rhs) = delete;
88 BARef& operator= (const BARef& rhs) = delete;
89 BARef& operator= (BARef&& rhs) = delete;
90
91 ~BARef ();
92
93 void define (const Box& bx);
94 void define (const BoxList& bl);
95 void define (BoxList&& bl) noexcept;
96 void define (std::istream& is, int& ndims);
98 void resize (Long n);
99#ifdef AMREX_MEM_PROFILING
100 void updateMemoryUsage_box (int s);
101 void updateMemoryUsage_hash (int s);
102#endif
103
104 [[nodiscard]] inline bool HasHashMap () const {
105 bool r;
106#ifdef AMREX_USE_OMP
107#pragma omp atomic read
108#endif
109 r = has_hashmap;
110 return r;
111 }
112
113 //
115 Vector<Box> m_abox;
116 //
118 mutable Box bbox;
119
120 mutable IntVect crsn;
121
122 using HashType = std::unordered_map< IntVect, std::vector<int>, IntVect::shift_hasher > ;
123 //using HashType = std::map< IntVect,std::vector<int> >;
124
125 mutable HashType hash;
126
127 mutable bool has_hashmap = false;
128
129 static int numboxarrays;
130 static int numboxarrays_hwm;
131 static Long total_box_bytes;
132 static Long total_box_bytes_hwm;
133 static Long total_hash_bytes;
134 static Long total_hash_bytes_hwm;
135
136 static void Initialize ();
137 static void Finalize ();
138 static bool initialized;
139};
141
143struct BATnull
144{
145 [[nodiscard]] Box operator() (const Box& bx) const noexcept { return bx; }
146 [[nodiscard]] static constexpr Box coarsen (Box const& a_box) { return a_box; }
147 [[nodiscard]] static constexpr IntVect doiLo () { return IntVect::TheZeroVector(); }
148 [[nodiscard]] static constexpr IntVect doiHi () { return IntVect::TheZeroVector(); }
149 [[nodiscard]] static constexpr IndexType index_type () { return IndexType(); }
150 [[nodiscard]] static constexpr IntVect coarsen_ratio () { return IntVect::TheUnitVector(); }
151};
153
155struct BATindexType
156{
157 explicit BATindexType (IndexType a_typ) : m_typ(a_typ) {}
158 [[nodiscard]] Box operator() (const Box& bx) const noexcept { return amrex::convert(bx,m_typ); }
159 [[nodiscard]] static Box coarsen (Box const& a_box) noexcept { return a_box; }
160 [[nodiscard]] static constexpr IntVect doiLo () { return IntVect::TheZeroVector(); }
161 [[nodiscard]] IntVect doiHi () const noexcept { return m_typ.ixType(); }
162 [[nodiscard]] IndexType index_type () const noexcept { return m_typ; }
163 [[nodiscard]] static constexpr IntVect coarsen_ratio () { return IntVect::TheUnitVector(); }
164 IndexType m_typ;
165};
167
169struct BATcoarsenRatio
170{
171 explicit BATcoarsenRatio (IntVect const& a_crse_ratio) : m_crse_ratio(a_crse_ratio) {}
172 [[nodiscard]] Box operator() (const Box& bx) const noexcept { return amrex::coarsen(bx,m_crse_ratio); }
173 [[nodiscard]] Box coarsen (Box const& a_box) const noexcept { return amrex::coarsen(a_box,m_crse_ratio); }
174 [[nodiscard]] static constexpr IntVect doiLo () { return IntVect::TheZeroVector(); }
175 [[nodiscard]] static constexpr IntVect doiHi () { return IntVect::TheZeroVector(); }
176 [[nodiscard]] static constexpr IndexType index_type () { return IndexType(); }
177 [[nodiscard]] IntVect coarsen_ratio () const noexcept { return m_crse_ratio; }
178 IntVect m_crse_ratio;
179};
181
183struct BATindexType_coarsenRatio
184{
185 BATindexType_coarsenRatio (IndexType a_typ, IntVect const& a_crse_ratio)
186 : m_typ(a_typ), m_crse_ratio(a_crse_ratio) {}
187
188 [[nodiscard]] Box operator() (const Box& bx) const noexcept {
189 return amrex::convert(amrex::coarsen(bx,m_crse_ratio),m_typ);
190 }
191
192 [[nodiscard]] Box coarsen (Box const& a_box) const noexcept { return amrex::coarsen(a_box,m_crse_ratio); }
193
194 [[nodiscard]] static constexpr IntVect doiLo () { return IntVect::TheZeroVector(); }
195 [[nodiscard]] IntVect doiHi () const noexcept { return m_typ.ixType(); }
196
197 [[nodiscard]] IndexType index_type () const noexcept { return m_typ; }
198 [[nodiscard]] IntVect coarsen_ratio () const noexcept { return m_crse_ratio; }
199
200 IndexType m_typ;
201 IntVect m_crse_ratio;
202};
204
206struct BATbndryReg
207{
208 BATbndryReg (Orientation a_face, IndexType a_typ,
209 int a_in_rad, int a_out_rad, int a_extent_rad)
210 : m_face(a_face), m_typ(a_typ), m_crse_ratio(1)
211 {
212 m_loshft = IntVect(-a_extent_rad);
213 m_hishft = IntVect( a_extent_rad);
214 IntVect nodal = a_typ.ixType();
215 m_hishft += nodal;
216 const int d = a_face.coordDir();
217 if (nodal[d]) {
218 // InterpFaceRegister & SyncRegister in IAMR
219 if (m_face.isLow()) {
220 m_loshft[d] = 0;
221 m_hishft[d] = 0;
222 } else {
223 m_loshft[d] = 1;
224 m_hishft[d] = 1;
225 }
226 m_doilo = IntVect(0);
227 m_doihi = nodal;
228 } else {
229 // BndryRegister
230 if (m_face.isLow()) {
231 m_loshft[d] = nodal[d] - a_out_rad;
232 m_hishft[d] = nodal[d] + a_in_rad - 1;
233 } else {
234 m_loshft[d] = 1 - a_in_rad;
235 m_hishft[d] = a_out_rad;
236 }
237 m_doilo = IntVect(a_extent_rad);
238 m_doihi = IntVect(a_extent_rad);
239 m_doihi += nodal;
240 if (m_face.isLow()) { // domain of influence in index space
241 m_doilo[d] = a_out_rad;
242 m_doihi[d] = 0;
243 } else {
244 m_doilo[d] = 0;
245 m_doihi[d] = a_out_rad;
246 }
247 }
248 }
249
250 [[nodiscard]] Box operator() (const Box& a_bx) const noexcept {
251 IntVect lo = amrex::coarsen(a_bx.smallEnd(), m_crse_ratio);
252 IntVect hi = amrex::coarsen(a_bx.bigEnd(), m_crse_ratio);
253 const int d = m_face.coordDir();
254 if (m_face.isLow()) {
255 hi[d] = lo[d];
256 } else {
257 lo[d] = hi[d];
258 }
259 lo += m_loshft;
260 hi += m_hishft;
261 return Box(lo,hi,m_typ);
262 }
263
264 [[nodiscard]] Box coarsen (Box const& a_box) const noexcept { return amrex::coarsen(a_box,m_crse_ratio); }
265
266 [[nodiscard]] IntVect doiLo () const noexcept { return m_doilo; }
267 [[nodiscard]] IntVect doiHi () const noexcept { return m_doihi; }
268
269 [[nodiscard]] IndexType index_type () const noexcept { return m_typ; }
270 [[nodiscard]] IntVect coarsen_ratio () const noexcept { return m_crse_ratio; }
271
272 friend bool operator== (BATbndryReg const& a, BATbndryReg const& b) noexcept {
273 return a.m_face == b.m_face && a.m_typ == b.m_typ && a.m_crse_ratio == b.m_crse_ratio
274 && a.m_loshft == b.m_loshft && a.m_hishft == b.m_hishft
275 && a.m_doilo == b.m_doilo && a.m_doihi == b.m_doihi;
276 }
277
278 Orientation m_face;
279 IndexType m_typ;
280 IntVect m_crse_ratio;
281 IntVect m_loshft;
282 IntVect m_hishft;
283 IntVect m_doilo;
284 IntVect m_doihi;
285};
287
289struct BATransformer
290{
291 enum struct BATType { null, indexType, coarsenRatio, indexType_coarsenRatio, bndryReg};
292
293 union BATOp {
294 BATOp () noexcept
295 : m_null() {}
296 BATOp (IndexType t) noexcept
297 : m_indexType(t) {}
298 BATOp (IntVect const& r) noexcept
299 : m_coarsenRatio(r) {}
300 BATOp (IndexType t, IntVect const& r) noexcept
301 : m_indexType_coarsenRatio(t,r) {}
302 BATOp (Orientation f, IndexType t, int in_rad, int out_rad, int extent_rad) noexcept
303 : m_bndryReg(f,t,in_rad,out_rad,extent_rad) {}
304 BATnull m_null;
305 BATindexType m_indexType;
306 BATcoarsenRatio m_coarsenRatio;
307 BATindexType_coarsenRatio m_indexType_coarsenRatio;
308 BATbndryReg m_bndryReg;
309 };
310
311 BATransformer () = default;
312
313 BATransformer (IndexType t)
314 : m_bat_type(t.cellCentered() ? BATType::null : BATType::indexType),
315 m_op (t.cellCentered() ? BATOp() : BATOp(t)) {}
316
317 BATransformer (Orientation f, IndexType t, int in_rad, int out_rad, int extent_rad)
318 : m_bat_type(BATType::bndryReg),
319 m_op(f,t,in_rad,out_rad,extent_rad) {}
320
321 [[nodiscard]] Box operator() (Box const& ab) const noexcept {
322 switch (m_bat_type)
323 {
324 case BATType::null:
325 return m_op.m_null(ab);
326 case BATType::indexType:
327 return m_op.m_indexType(ab);
328 case BATType::coarsenRatio:
329 return m_op.m_coarsenRatio(ab);
330 case BATType::indexType_coarsenRatio:
331 return m_op.m_indexType_coarsenRatio(ab);
332 default:
333 return m_op.m_bndryReg(ab);
334 }
335 }
336
337 [[nodiscard]] Box coarsen (Box const& a_box) const noexcept {
338 switch (m_bat_type)
339 {
340 case BATType::null:
341 return amrex::BATnull::coarsen(a_box);
342 case BATType::indexType:
343 return amrex::BATindexType::coarsen(a_box);
344 case BATType::coarsenRatio:
345 return m_op.m_coarsenRatio.coarsen(a_box);
346 case BATType::indexType_coarsenRatio:
347 return m_op.m_indexType_coarsenRatio.coarsen(a_box);
348 default:
349 return m_op.m_bndryReg.coarsen(a_box);
350 }
351 }
352
353 [[nodiscard]] IntVect doiLo () const noexcept {
354 switch (m_bat_type)
355 {
356 case BATType::null:
357 return amrex::BATnull::doiLo();
358 case BATType::indexType:
359 return amrex::BATindexType::doiLo();
360 case BATType::coarsenRatio:
361 return amrex::BATcoarsenRatio::doiLo();
362 case BATType::indexType_coarsenRatio:
363 return amrex::BATindexType_coarsenRatio::doiLo();
364 default:
365 return m_op.m_bndryReg.doiLo();
366 }
367 }
368
369 [[nodiscard]] IntVect doiHi () const noexcept {
370 switch (m_bat_type)
371 {
372 case BATType::null:
373 return amrex::BATnull::doiHi();
374 case BATType::indexType:
375 return m_op.m_indexType.doiHi();
376 case BATType::coarsenRatio:
377 return amrex::BATcoarsenRatio::doiHi();
378 case BATType::indexType_coarsenRatio:
379 return m_op.m_indexType_coarsenRatio.doiHi();
380 default:
381 return m_op.m_bndryReg.doiHi();
382 }
383 }
384
385 [[nodiscard]] IndexType index_type () const noexcept {
386 switch (m_bat_type)
387 {
388 case BATType::null:
389 return amrex::BATnull::index_type();
390 case BATType::indexType:
391 return m_op.m_indexType.index_type();
392 case BATType::coarsenRatio:
393 return amrex::BATcoarsenRatio::index_type();
394 case BATType::indexType_coarsenRatio:
395 return m_op.m_indexType_coarsenRatio.index_type();
396 default:
397 return m_op.m_bndryReg.index_type();
398 }
399 }
400
401 [[nodiscard]] IntVect coarsen_ratio () const noexcept {
402 switch (m_bat_type)
403 {
404 case BATType::null:
405 return amrex::BATnull::coarsen_ratio();
406 case BATType::indexType:
407 return amrex::BATindexType::coarsen_ratio();
408 case BATType::coarsenRatio:
409 return m_op.m_coarsenRatio.coarsen_ratio();
410 case BATType::indexType_coarsenRatio:
411 return m_op.m_indexType_coarsenRatio.coarsen_ratio();
412 default:
413 return m_op.m_bndryReg.coarsen_ratio();
414 }
415 }
416
417 [[nodiscard]] bool is_null () const noexcept {
418 return m_bat_type == BATType::null;
419 }
420
421 [[nodiscard]] bool is_simple () const noexcept {
422 return m_bat_type != BATType::bndryReg;
423 }
424
425 void set_coarsen_ratio (IntVect const& a_ratio) noexcept {
426 switch (m_bat_type)
427 {
428 case BATType::null:
429 {
430 if (a_ratio == IntVect::TheUnitVector()) {
431 return;
432 } else {
433 m_bat_type = BATType::coarsenRatio;
434 m_op.m_coarsenRatio.m_crse_ratio = a_ratio;
435 return;
436 }
437 }
438 case BATType::indexType:
439 {
440 if (a_ratio == IntVect::TheUnitVector()) {
441 return;
442 } else {
443 m_bat_type = BATType::indexType_coarsenRatio;
444 auto t = m_op.m_indexType.m_typ;
445 m_op.m_indexType_coarsenRatio.m_typ = t;
446 m_op.m_indexType_coarsenRatio.m_crse_ratio = a_ratio;
447 return;
448 }
449 }
450 case BATType::coarsenRatio:
451 {
452 if (a_ratio == IntVect::TheUnitVector()) {
453 m_bat_type = BATType::null;
454 return;
455 } else {
456 m_op.m_coarsenRatio.m_crse_ratio = a_ratio;
457 return;
458 }
459 }
460 case BATType::indexType_coarsenRatio:
461 {
462 if (a_ratio == IntVect::TheUnitVector()) {
463 m_bat_type = BATType::indexType;
464 auto t = m_op.m_indexType_coarsenRatio.m_typ;
465 m_op.m_indexType.m_typ = t;
466 return;
467 } else {
468 m_op.m_indexType_coarsenRatio.m_crse_ratio = a_ratio;
469 return;
470 }
471 }
472 default:
473 {
474 m_op.m_bndryReg.m_crse_ratio = a_ratio;
475 return;
476 }
477 }
478 }
479
480 void set_index_type (IndexType typ) noexcept {
481 switch (m_bat_type)
482 {
483 case BATType::null:
484 {
485 if (typ.cellCentered()) {
486 return;
487 } else {
488 m_bat_type = BATType::indexType;
489 m_op.m_indexType.m_typ = typ;
490 return;
491 }
492 }
493 case BATType::indexType:
494 {
495 if (typ.cellCentered()) {
496 m_bat_type = BATType::null;
497 return;
498 } else {
499 m_op.m_indexType.m_typ = typ;
500 return;
501 }
502 }
503 case BATType::coarsenRatio:
504 {
505 if (typ.cellCentered()) {
506 return;
507 } else {
508 m_bat_type = BATType::indexType_coarsenRatio;
509 auto r = m_op.m_coarsenRatio.m_crse_ratio;
510 m_op.m_indexType_coarsenRatio.m_typ = typ;
511 m_op.m_indexType_coarsenRatio.m_crse_ratio = r;
512 return;
513 }
514 }
515 case BATType::indexType_coarsenRatio:
516 {
517 if (typ.cellCentered()) {
518 m_bat_type = BATType::coarsenRatio;
519 auto r = m_op.m_indexType_coarsenRatio.m_crse_ratio;
520 m_op.m_coarsenRatio.m_crse_ratio = r;
521 return;
522 } else {
523 m_op.m_indexType_coarsenRatio.m_typ = typ;
524 return;
525 }
526 }
527 default:
528 {
529 m_op.m_bndryReg.m_typ = typ;
530 return;
531 }
532 }
533 }
534
535 friend bool operator== (BATransformer const& a, BATransformer const& b) noexcept {
536 if (a.m_bat_type != BATType::bndryReg && b.m_bat_type != BATType::bndryReg) {
537 return a.index_type() == b.index_type()
538 && a.coarsen_ratio() == b.coarsen_ratio();
539 } else if (a.m_bat_type == BATType::bndryReg && b.m_bat_type == BATType::bndryReg) {
540 return a.m_op.m_bndryReg == b.m_op.m_bndryReg;
541 } else {
542 return false;
543 }
544 }
545
546 BATType m_bat_type{BATType::null};
547 BATOp m_op;
548};
550
551// for backward compatibility
552using BndryBATransformer = BATransformer;
553
554class MFIter;
555class AmrMesh;
556class FabArrayBase;
557
567{
568public:
569
570 BoxArray () noexcept;
571 BoxArray (const BoxArray& rhs) = default;
572 BoxArray (BoxArray&& rhs) noexcept = default;
573 BoxArray& operator= (BoxArray const& rhs) = default;
574 BoxArray& operator= (BoxArray&& rhs) noexcept = default;
575 ~BoxArray() noexcept = default;
576
578 explicit BoxArray (const Box& bx);
579
581 explicit BoxArray (size_t n);
582
584 BoxArray (const Box* bxvec,
585 int nbox);
586
588 explicit BoxArray (const BoxList& bl);
589 explicit BoxArray (BoxList&& bl) noexcept;
590
591 BoxArray (const BoxArray& rhs, const BATransformer& trans);
592
593 BoxArray (BoxList&& bl, IntVect const& max_grid_size);
594
599 void define (const Box& bx);
604 void define (const BoxList& bl);
605 void define (BoxList&& bl) noexcept;
606
608 void clear ();
609
611 void resize (Long len);
612
614 [[nodiscard]] Long size () const noexcept { return m_ref->m_abox.size(); }
615
617 [[nodiscard]] Long capacity () const noexcept { return static_cast<Long>(m_ref->m_abox.capacity()); }
618
620 [[nodiscard]] bool empty () const noexcept { return m_ref->m_abox.empty(); }
621
623 [[nodiscard]] Long numPts() const noexcept;
624
626 [[nodiscard]] double d_numPts () const noexcept;
633 int readFrom (std::istream& is);
634
636 std::ostream& writeOn (std::ostream&) const;
637
639 [[nodiscard]] bool operator== (const BoxArray& rhs) const noexcept;
640
642 [[nodiscard]] bool operator!= (const BoxArray& rhs) const noexcept;
643
644 [[nodiscard]] bool operator== (const Vector<Box>& bv) const noexcept;
645 [[nodiscard]] bool operator!= (const Vector<Box>& bv) const noexcept;
646
648 [[nodiscard]] bool CellEqual (const BoxArray& rhs) const noexcept;
649
651 BoxArray& maxSize (int block_size);
652
653 BoxArray& maxSize (const IntVect& block_size);
654
658 BoxArray& minmaxSize (const IntVect& min_size, const IntVect& max_size);
659
661 BoxArray& refine (int refinement_ratio);
662
664 BoxArray& refine (const IntVect& iv);
665
667 BoxArray& coarsen (int refinement_ratio);
668
670 [[nodiscard]] bool coarsenable (int refinement_ratio, int min_width=1) const;
671 [[nodiscard]] bool coarsenable (const IntVect& refinement_ratio, int min_width=1) const;
672 [[nodiscard]] bool coarsenable (const IntVect& refinement_ratio, const IntVect& min_width) const;
673
675 BoxArray& coarsen (const IntVect& iv);
676
678 BoxArray& growcoarsen (int n, const IntVect& iv);
679 BoxArray& growcoarsen (IntVect const& ngrow, const IntVect& iv);
680
682 BoxArray& grow (int n);
683
685 BoxArray& grow (const IntVect& iv);
690 BoxArray& grow (int idir, int n_cell);
695 BoxArray& growLo (int idir, int n_cell);
700 BoxArray& growHi (int idir, int n_cell);
710 BoxArray& surroundingNodes (int dir);
711
714
716 BoxArray& enclosedCells (int dir);
717
720
721 BoxArray& convert (const IntVect& iv);
722
724 BoxArray& convert (Box (*fp)(const Box&));
725
727 BoxArray& shift (int dir, int nzones);
728
730 BoxArray& shift (const IntVect &iv);
731
733 void set (int i, const Box& ibox);
734
736 [[nodiscard]] Box operator[] (int index) const noexcept {
737 return m_bat(m_ref->m_abox[index]);
738 }
739
741 [[nodiscard]] Box operator[] (const MFIter& mfi) const noexcept;
742
744 [[nodiscard]] Box get (int index) const noexcept { return operator[](index); }
745
747 [[nodiscard]] Box getCellCenteredBox (int index) const noexcept {
748 return m_bat.coarsen(m_ref->m_abox[index]);
749 }
750
755 [[nodiscard]] bool ok () const;
756
758 [[nodiscard]] bool isDisjoint () const;
759
761 [[nodiscard]] BoxList boxList () const;
762
764 [[nodiscard]] bool contains (const IntVect& v) const;
765
770 [[nodiscard]]
771 bool contains (const Box& b, bool assume_disjoint_ba = false,
772 const IntVect& ng = IntVect(0)) const;
773
777 [[nodiscard]]
778 bool contains (const BoxArray& ba, bool assume_disjoint_ba = false,
779 const IntVect& ng = IntVect(0)) const;
780
788 [[nodiscard]]
789 bool contains (const BoxArray& ba, Periodicity const& period) const;
790
792 [[nodiscard]] Box minimalBox () const;
793 [[nodiscard]] Box minimalBox (Long& npts_avg_box) const;
794
799 [[nodiscard]]
800 bool intersects (const Box& b, int ng = 0) const;
801
802 [[nodiscard]]
803 bool intersects (const Box& b, const IntVect& ng) const;
804
806 [[nodiscard]]
807 std::vector< std::pair<int,Box> > intersections (const Box& bx) const;
808
810 [[nodiscard]]
811 std::vector< std::pair<int,Box> > intersections (const Box& bx, bool first_only, int ng) const;
812
813 [[nodiscard]]
814 std::vector< std::pair<int,Box> > intersections (const Box& bx, bool first_only, const IntVect& ng) const;
815
817 void intersections (const Box& bx, std::vector< std::pair<int,Box> >& isects) const;
818
820 void intersections (const Box& bx, std::vector< std::pair<int,Box> >& isects,
821 bool first_only, int ng) const;
822
823 void intersections (const Box& bx, std::vector< std::pair<int,Box> >& isects,
824 bool first_only, const IntVect& ng) const;
825
827 [[nodiscard]] BoxList complementIn (const Box& b) const;
828 void complementIn (BoxList& bl, const Box& b) const;
829
831 [[nodiscard]] BoxList complementIn (const Box& b, const Periodicity& period) const;
832
834 void clear_hash_bin () const;
835
837 void removeOverlap (bool simplify=true);
838
840 [[nodiscard]] static bool SameRefs (const BoxArray& lhs, const BoxArray& rhs) { return lhs.m_ref == rhs.m_ref; }
841
842 struct RefID {
843 RefID () noexcept = default;
844 explicit RefID (BARef* data_) noexcept : data(data_) {}
845 bool operator< (const RefID& rhs) const noexcept { return std::less<>()(data,rhs.data); }
846 bool operator== (const RefID& rhs) const noexcept { return data == rhs.data; }
847 bool operator!= (const RefID& rhs) const noexcept { return data != rhs.data; }
848 friend std::ostream& operator<< (std::ostream& os, const RefID& id);
849 private:
850 BARef* data{nullptr};
851 };
852
854 [[nodiscard]] RefID getRefID () const noexcept { return RefID { m_ref.get() }; }
855
857 [[nodiscard]] IndexType ixType () const noexcept { return m_bat.index_type(); }
858
860 [[nodiscard]] IntVect crseRatio () const noexcept { return m_bat.coarsen_ratio(); }
861
862 static void Initialize ();
863 static void Finalize ();
864 static bool initialized;
865
867 void uniqify ();
868
869 [[nodiscard]] BoxList const& simplified_list () const; // For regular AMR grids only!!!
870 [[nodiscard]] BoxArray simplified () const; // For regular AMR grids only!!!
871
872 [[nodiscard]] BATransformer const& transformer () const;
873
874 [[nodiscard]] std::weak_ptr<BARef> getWeakRef () const;
875 [[nodiscard]] std::shared_ptr<BARef> const& getSharedRef () const;
876 std::shared_ptr<BARef>& getSharedRef ();
877
878 friend class AmrMesh;
879 friend class FabArrayBase;
880
881private:
883 void type_update ();
884
885 [[nodiscard]] BARef::HashType& getHashMap () const;
886
887 [[nodiscard]] IntVect getDoiLo () const noexcept;
888 [[nodiscard]] IntVect getDoiHi () const noexcept;
889
890 BATransformer m_bat;
892 std::shared_ptr<BARef> m_ref;
893 mutable std::shared_ptr<BoxList> m_simplified_list;
894};
895
897std::ostream& operator<< (std::ostream& os, const BoxArray& ba);
898
899std::ostream& operator<< (std::ostream& os, const BoxArray::RefID& id);
900
901}
902
903#endif /*BL_BOXARRAY_H*/
int idir
Definition AMReX_HypreMLABecLap.cpp:1093
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
Definition AMReX_AmrMesh.H:63
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:567
IndexType ixType() const noexcept
Return index type of this BoxArray.
Definition AMReX_BoxArray.H:857
Box operator[](int index) const noexcept
Return element index of this BoxArray.
Definition AMReX_BoxArray.H:736
Long capacity() const noexcept
Return the number of boxes that can be held in the current allocated storage.
Definition AMReX_BoxArray.H:617
RefID getRefID() const noexcept
Return a unique ID of the reference.
Definition AMReX_BoxArray.H:854
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
BATransformer m_bat
Definition AMReX_BoxArray.H:890
BoxList const & simplified_list() const
Definition AMReX_BoxArray.cpp:1633
BoxList boxList() const
Create a BoxList from this BoxArray.
Definition AMReX_BoxArray.cpp:949
BoxArray simplified() const
Definition AMReX_BoxArray.cpp:1644
IntVect crseRatio() const noexcept
Return crse ratio of this BoxArray.
Definition AMReX_BoxArray.H:860
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
IntVect getDoiLo() const noexcept
Definition AMReX_BoxArray.cpp:1542
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
std::shared_ptr< BoxList > m_simplified_list
Definition AMReX_BoxArray.H:893
BARef::HashType & getHashMap() const
IntVect getDoiHi() const noexcept
Definition AMReX_BoxArray.cpp:1548
void clear_hash_bin() const
Clear out the internal hash table used by intersections.
Definition AMReX_BoxArray.cpp:1433
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:1656
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:1662
Box getCellCenteredBox(int index) const noexcept
Return cell-centered box at element index of this BoxArray.
Definition AMReX_BoxArray.H:747
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:1650
static bool initialized
Definition AMReX_BoxArray.H:864
static bool SameRefs(const BoxArray &lhs, const BoxArray &rhs)
whether two BoxArrays share the same data
Definition AMReX_BoxArray.H:840
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 type_update()
Update BoxArray index type according the box type, and then convert boxes to cell-centered.
Definition AMReX_BoxArray.cpp:1528
std::shared_ptr< BARef > m_ref
The data – a reference-counted pointer to a Ref.
Definition AMReX_BoxArray.H:892
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:1449
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:614
Box get(int index) const noexcept
Return element index of this BoxArray.
Definition AMReX_BoxArray.H:744
~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:1610
BoxArray(const BoxArray &rhs)=default
bool empty() const noexcept
Return whether the BoxArray is empty.
Definition AMReX_BoxArray.H:620
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
__host__ __device__ BoxND & coarsen(int ref_ratio) noexcept
Coarsen BoxND by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo/rati...
Definition AMReX_Box.H:722
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:63
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:49
__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:1714
bool initialized
Definition AMReX_DistributionMapping.cpp:32
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:332
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:27
BoxList GetBndryCells(const BoxArray &ba, int ngrow)
Find the ghost cells of a given BoxArray.
Definition AMReX_BoxArray.cpp:1839
IndexTypeND< 3 > IndexType
IndexType is an alias for amrex::IndexTypeND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
void Finalize(AMReX *pamrex)
Definition AMReX.cpp:792
BoxArray decompose(Box const &domain, int nboxes, Array< bool, 3 > const &decomp, bool no_overlap)
Decompose domain box into BoxArray.
Definition AMReX_BoxArray.cpp:1940
bool match(const BoxArray &x, const BoxArray &y)
Note that two BoxArrays that match are not necessarily equal.
Definition AMReX_BoxArray.cpp:1927
BoxArray complementIn(const Box &b, const BoxArray &ba)
Make a BoxArray from the complement of BoxArray ba in Box b.
Definition AMReX_BoxArray.cpp:1707
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
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:1897
BATransformer BndryBATransformer
Definition AMReX_BoxArray.H:552
BoxArray boxComplement(const Box &b1in, const Box &b2)
Make a BoxArray from the the complement of b2 in b1in.
Definition AMReX_BoxArray.cpp:1700
Definition AMReX_BoxArray.H:842
friend std::ostream & operator<<(std::ostream &os, const RefID &id)
Definition AMReX_BoxArray.cpp:2123
RefID() noexcept=default
bool operator<(const RefID &rhs) const noexcept
Definition AMReX_BoxArray.H:845
bool operator!=(const RefID &rhs) const noexcept
Definition AMReX_BoxArray.H:847
bool operator==(const RefID &rhs) const noexcept
Definition AMReX_BoxArray.H:846
BARef * data
Definition AMReX_BoxArray.H:850