Block-Structured AMR Software Framework
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 <unordered_map>
16 
17 namespace amrex
18 {
19  class BoxArray;
20 
22  [[nodiscard]] BoxArray boxComplement (const Box& b1in, const Box& b2);
23 
25  [[nodiscard]] BoxArray complementIn (const Box& b, const BoxArray& ba);
26 
28  [[nodiscard]] BoxArray intersect (const BoxArray& ba, const Box& b, int ng = 0);
29 
30  [[nodiscard]] BoxArray intersect (const BoxArray& ba, const Box& b, const IntVect& ng);
31 
33  [[nodiscard]] BoxArray intersect (const BoxArray& lhs, const BoxArray& rhs);
34 
36  [[nodiscard]] BoxList intersect (const BoxArray& ba, const BoxList& bl);
37 
38  [[nodiscard]] BoxArray convert (const BoxArray& ba, IndexType typ);
39  [[nodiscard]] BoxArray convert (const BoxArray& ba, const IntVect& typ);
40 
41  [[nodiscard]] BoxArray coarsen (const BoxArray& ba, int ratio);
42  [[nodiscard]] BoxArray coarsen (const BoxArray& ba, const IntVect& ratio);
43 
44  [[nodiscard]] BoxArray refine (const BoxArray& ba, int ratio);
45  [[nodiscard]] BoxArray refine (const BoxArray& ba, const IntVect& ratio);
46 
48  [[nodiscard]] BoxList GetBndryCells (const BoxArray& ba, int ngrow);
49 
51  void readBoxArray (BoxArray& ba, std::istream& s, bool b = false);
52 
54  [[nodiscard]] bool match (const BoxArray& x, const BoxArray& y);
55 
56 struct BARef
57 {
58  BARef ();
59  explicit BARef (size_t size);
60  explicit BARef (const Box& b);
61  explicit BARef (const BoxList& bl);
62  explicit BARef (BoxList&& bl) noexcept;
63  explicit BARef (std::istream& is);
64  BARef (const BARef& rhs);
65  BARef (BARef&& rhs) = delete;
66  BARef& operator= (const BARef& rhs) = delete;
67  BARef& operator= (BARef&& rhs) = delete;
68 
69  ~BARef ();
70 
71  void define (const Box& bx);
72  void define (const BoxList& bl);
73  void define (BoxList&& bl) noexcept;
74  void define (std::istream& is, int& ndims);
76  void resize (Long n);
77 #ifdef AMREX_MEM_PROFILING
78  void updateMemoryUsage_box (int s);
79  void updateMemoryUsage_hash (int s);
80 #endif
81 
82  [[nodiscard]] inline bool HasHashMap () const {
83  bool r;
84 #ifdef AMREX_USE_OMP
85 #pragma omp atomic read
86 #endif
87  r = has_hashmap;
88  return r;
89  }
90 
91  //
94  //
96  mutable Box bbox;
97 
98  mutable IntVect crsn;
99 
100  using HashType = std::unordered_map< IntVect, std::vector<int>, IntVect::shift_hasher > ;
101  //using HashType = std::map< IntVect,std::vector<int> >;
102 
103  mutable HashType hash;
104 
105  mutable bool has_hashmap = false;
106 
107  static int numboxarrays;
108  static int numboxarrays_hwm;
109  static Long total_box_bytes;
110  static Long total_box_bytes_hwm;
111  static Long total_hash_bytes;
112  static Long total_hash_bytes_hwm;
113 
114  static void Initialize ();
115  static void Finalize ();
116  static bool initialized;
117 };
118 
119 struct BATnull
120 {
121  [[nodiscard]] Box operator() (const Box& bx) const noexcept { return bx; }
122  [[nodiscard]] static constexpr Box coarsen (Box const& a_box) { return a_box; }
123  [[nodiscard]] static constexpr IntVect doiLo () { return IntVect::TheZeroVector(); }
124  [[nodiscard]] static constexpr IntVect doiHi () { return IntVect::TheZeroVector(); }
125  [[nodiscard]] static constexpr IndexType index_type () { return IndexType(); }
126  [[nodiscard]] static constexpr IntVect coarsen_ratio () { return IntVect::TheUnitVector(); }
127 };
128 
130 {
131  explicit BATindexType (IndexType a_typ) : m_typ(a_typ) {}
132  [[nodiscard]] Box operator() (const Box& bx) const noexcept { return amrex::convert(bx,m_typ); }
133  [[nodiscard]] static Box coarsen (Box const& a_box) noexcept { return a_box; }
134  [[nodiscard]] static constexpr IntVect doiLo () { return IntVect::TheZeroVector(); }
135  [[nodiscard]] IntVect doiHi () const noexcept { return m_typ.ixType(); }
136  [[nodiscard]] IndexType index_type () const noexcept { return m_typ; }
137  [[nodiscard]] static constexpr IntVect coarsen_ratio () { return IntVect::TheUnitVector(); }
139 };
140 
142 {
143  explicit BATcoarsenRatio (IntVect const& a_crse_ratio) : m_crse_ratio(a_crse_ratio) {}
144  [[nodiscard]] Box operator() (const Box& bx) const noexcept { return amrex::coarsen(bx,m_crse_ratio); }
145  [[nodiscard]] Box coarsen (Box const& a_box) const noexcept { return amrex::coarsen(a_box,m_crse_ratio); }
146  [[nodiscard]] static constexpr IntVect doiLo () { return IntVect::TheZeroVector(); }
147  [[nodiscard]] static constexpr IntVect doiHi () { return IntVect::TheZeroVector(); }
148  [[nodiscard]] static constexpr IndexType index_type () { return IndexType(); }
149  [[nodiscard]] IntVect coarsen_ratio () const noexcept { return m_crse_ratio; }
151 };
152 
154 {
155  BATindexType_coarsenRatio (IndexType a_typ, IntVect const& a_crse_ratio)
156  : m_typ(a_typ), m_crse_ratio(a_crse_ratio) {}
157 
158  [[nodiscard]] Box operator() (const Box& bx) const noexcept {
160  }
161 
162  [[nodiscard]] Box coarsen (Box const& a_box) const noexcept { return amrex::coarsen(a_box,m_crse_ratio); }
163 
164  [[nodiscard]] static constexpr IntVect doiLo () { return IntVect::TheZeroVector(); }
165  [[nodiscard]] IntVect doiHi () const noexcept { return m_typ.ixType(); }
166 
167  [[nodiscard]] IndexType index_type () const noexcept { return m_typ; }
168  [[nodiscard]] IntVect coarsen_ratio () const noexcept { return m_crse_ratio; }
169 
172 };
173 
175 {
177  int a_in_rad, int a_out_rad, int a_extent_rad)
178  : m_face(a_face), m_typ(a_typ), m_crse_ratio(1)
179  {
180  m_loshft = IntVect(-a_extent_rad);
181  m_hishft = IntVect( a_extent_rad);
182  IntVect nodal = a_typ.ixType();
183  m_hishft += nodal;
184  const int d = a_face.coordDir();
185  if (nodal[d]) {
186  // InterpFaceRegister & SyncRegister in IAMR
187  if (m_face.isLow()) {
188  m_loshft[d] = 0;
189  m_hishft[d] = 0;
190  } else {
191  m_loshft[d] = 1;
192  m_hishft[d] = 1;
193  }
194  m_doilo = IntVect(0);
195  m_doihi = nodal;
196  } else {
197  // BndryRegister
198  if (m_face.isLow()) {
199  m_loshft[d] = nodal[d] - a_out_rad;
200  m_hishft[d] = nodal[d] + a_in_rad - 1;
201  } else {
202  m_loshft[d] = 1 - a_in_rad;
203  m_hishft[d] = a_out_rad;
204  }
205  m_doilo = IntVect(a_extent_rad);
206  m_doihi = IntVect(a_extent_rad);
207  m_doihi += nodal;
208  if (m_face.isLow()) { // domain of influence in index space
209  m_doilo[d] = a_out_rad;
210  m_doihi[d] = 0;
211  } else {
212  m_doilo[d] = 0;
213  m_doihi[d] = a_out_rad;
214  }
215  }
216  }
217 
218  [[nodiscard]] Box operator() (const Box& a_bx) const noexcept {
219  IntVect lo = amrex::coarsen(a_bx.smallEnd(), m_crse_ratio);
220  IntVect hi = amrex::coarsen(a_bx.bigEnd(), m_crse_ratio);
221  const int d = m_face.coordDir();
222  if (m_face.isLow()) {
223  hi[d] = lo[d];
224  } else {
225  lo[d] = hi[d];
226  }
227  lo += m_loshft;
228  hi += m_hishft;
229  return Box(lo,hi,m_typ);
230  }
231 
232  [[nodiscard]] Box coarsen (Box const& a_box) const noexcept { return amrex::coarsen(a_box,m_crse_ratio); }
233 
234  [[nodiscard]] IntVect doiLo () const noexcept { return m_doilo; }
235  [[nodiscard]] IntVect doiHi () const noexcept { return m_doihi; }
236 
237  [[nodiscard]] IndexType index_type () const noexcept { return m_typ; }
238  [[nodiscard]] IntVect coarsen_ratio () const noexcept { return m_crse_ratio; }
239 
240  friend bool operator== (BATbndryReg const& a, BATbndryReg const& b) noexcept {
241  return a.m_face == b.m_face && a.m_typ == b.m_typ && a.m_crse_ratio == b.m_crse_ratio
242  && a.m_loshft == b.m_loshft && a.m_hishft == b.m_hishft
243  && a.m_doilo == b.m_doilo && a.m_doihi == b.m_doihi;
244  }
245 
253 };
254 
256 {
258 
259  union BATOp {
260  BATOp () noexcept
261  : m_null() {}
262  BATOp (IndexType t) noexcept
263  : m_indexType(t) {}
264  BATOp (IntVect const& r) noexcept
265  : m_coarsenRatio(r) {}
266  BATOp (IndexType t, IntVect const& r) noexcept
267  : m_indexType_coarsenRatio(t,r) {}
268  BATOp (Orientation f, IndexType t, int in_rad, int out_rad, int extent_rad) noexcept
269  : m_bndryReg(f,t,in_rad,out_rad,extent_rad) {}
275  };
276 
277  BATransformer () = default;
278 
280  : m_bat_type(t.cellCentered() ? BATType::null : BATType::indexType),
281  m_op (t.cellCentered() ? BATOp() : BATOp(t)) {}
282 
283  BATransformer (Orientation f, IndexType t, int in_rad, int out_rad, int extent_rad)
284  : m_bat_type(BATType::bndryReg),
285  m_op(f,t,in_rad,out_rad,extent_rad) {}
286 
287  [[nodiscard]] Box operator() (Box const& ab) const noexcept {
288  switch (m_bat_type)
289  {
290  case BATType::null:
291  return m_op.m_null(ab);
292  case BATType::indexType:
293  return m_op.m_indexType(ab);
295  return m_op.m_coarsenRatio(ab);
297  return m_op.m_indexType_coarsenRatio(ab);
298  default:
299  return m_op.m_bndryReg(ab);
300  }
301  }
302 
303  [[nodiscard]] Box coarsen (Box const& a_box) const noexcept {
304  switch (m_bat_type)
305  {
306  case BATType::null:
307  return amrex::BATnull::coarsen(a_box);
308  case BATType::indexType:
309  return amrex::BATindexType::coarsen(a_box);
311  return m_op.m_coarsenRatio.coarsen(a_box);
313  return m_op.m_indexType_coarsenRatio.coarsen(a_box);
314  default:
315  return m_op.m_bndryReg.coarsen(a_box);
316  }
317  }
318 
319  [[nodiscard]] IntVect doiLo () const noexcept {
320  switch (m_bat_type)
321  {
322  case BATType::null:
323  return amrex::BATnull::doiLo();
324  case BATType::indexType:
330  default:
331  return m_op.m_bndryReg.doiLo();
332  }
333  }
334 
335  [[nodiscard]] IntVect doiHi () const noexcept {
336  switch (m_bat_type)
337  {
338  case BATType::null:
339  return amrex::BATnull::doiHi();
340  case BATType::indexType:
341  return m_op.m_indexType.doiHi();
346  default:
347  return m_op.m_bndryReg.doiHi();
348  }
349  }
350 
351  [[nodiscard]] IndexType index_type () const noexcept {
352  switch (m_bat_type)
353  {
354  case BATType::null:
356  case BATType::indexType:
357  return m_op.m_indexType.index_type();
362  default:
363  return m_op.m_bndryReg.index_type();
364  }
365  }
366 
367  [[nodiscard]] IntVect coarsen_ratio () const noexcept {
368  switch (m_bat_type)
369  {
370  case BATType::null:
372  case BATType::indexType:
378  default:
379  return m_op.m_bndryReg.coarsen_ratio();
380  }
381  }
382 
383  [[nodiscard]] bool is_null () const noexcept {
384  return m_bat_type == BATType::null;
385  }
386 
387  [[nodiscard]] bool is_simple () const noexcept {
388  return m_bat_type != BATType::bndryReg;
389  }
390 
391  void set_coarsen_ratio (IntVect const& a_ratio) noexcept {
392  switch (m_bat_type)
393  {
394  case BATType::null:
395  {
396  if (a_ratio == IntVect::TheUnitVector()) {
397  return;
398  } else {
400  m_op.m_coarsenRatio.m_crse_ratio = a_ratio;
401  return;
402  }
403  }
404  case BATType::indexType:
405  {
406  if (a_ratio == IntVect::TheUnitVector()) {
407  return;
408  } else {
410  auto t = m_op.m_indexType.m_typ;
413  return;
414  }
415  }
417  {
418  if (a_ratio == IntVect::TheUnitVector()) {
420  return;
421  } else {
422  m_op.m_coarsenRatio.m_crse_ratio = a_ratio;
423  return;
424  }
425  }
427  {
428  if (a_ratio == IntVect::TheUnitVector()) {
431  m_op.m_indexType.m_typ = t;
432  return;
433  } else {
435  return;
436  }
437  }
438  default:
439  {
440  m_op.m_bndryReg.m_crse_ratio = a_ratio;
441  return;
442  }
443  }
444  }
445 
446  void set_index_type (IndexType typ) noexcept {
447  switch (m_bat_type)
448  {
449  case BATType::null:
450  {
451  if (typ.cellCentered()) {
452  return;
453  } else {
455  m_op.m_indexType.m_typ = typ;
456  return;
457  }
458  }
459  case BATType::indexType:
460  {
461  if (typ.cellCentered()) {
463  return;
464  } else {
465  m_op.m_indexType.m_typ = typ;
466  return;
467  }
468  }
470  {
471  if (typ.cellCentered()) {
472  return;
473  } else {
478  return;
479  }
480  }
482  {
483  if (typ.cellCentered()) {
487  return;
488  } else {
490  return;
491  }
492  }
493  default:
494  {
495  m_op.m_bndryReg.m_typ = typ;
496  return;
497  }
498  }
499  }
500 
501  friend bool operator== (BATransformer const& a, BATransformer const& b) noexcept {
502  if (a.m_bat_type != BATType::bndryReg && b.m_bat_type != BATType::bndryReg) {
503  return a.index_type() == b.index_type()
504  && a.coarsen_ratio() == b.coarsen_ratio();
505  } else if (a.m_bat_type == BATType::bndryReg && b.m_bat_type == BATType::bndryReg) {
506  return a.m_op.m_bndryReg == b.m_op.m_bndryReg;
507  } else {
508  return false;
509  }
510  }
511 
514 };
515 
516 // for backward compatibility
518 
519 class MFIter;
520 class AmrMesh;
521 class FabArrayBase;
522 
529 class BoxArray
530 {
531 public:
532 
533  BoxArray () noexcept;
534  BoxArray (const BoxArray& rhs) = default;
535  BoxArray (BoxArray&& rhs) noexcept = default;
536  BoxArray& operator= (BoxArray const& rhs) = default;
537  BoxArray& operator= (BoxArray&& rhs) noexcept = default;
538  ~BoxArray() noexcept = default;
539 
541  explicit BoxArray (const Box& bx);
542 
544  explicit BoxArray (size_t n);
545 
547  BoxArray (const Box* bxvec,
548  int nbox);
549 
551  explicit BoxArray (const BoxList& bl);
552  explicit BoxArray (BoxList&& bl) noexcept;
553 
554  BoxArray (const BoxArray& rhs, const BATransformer& trans);
555 
556  BoxArray (BoxList&& bl, IntVect const& max_grid_size);
557 
562  void define (const Box& bx);
567  void define (const BoxList& bl);
568  void define (BoxList&& bl) noexcept;
569 
571  void clear ();
572 
574  void resize (Long len);
575 
577  [[nodiscard]] Long size () const noexcept { return m_ref->m_abox.size(); }
578 
580  [[nodiscard]] Long capacity () const noexcept { return static_cast<Long>(m_ref->m_abox.capacity()); }
581 
583  [[nodiscard]] bool empty () const noexcept { return m_ref->m_abox.empty(); }
584 
586  [[nodiscard]] Long numPts() const noexcept;
587 
589  [[nodiscard]] double d_numPts () const noexcept;
596  int readFrom (std::istream& is);
597 
599  std::ostream& writeOn (std::ostream&) const;
600 
602  [[nodiscard]] bool operator== (const BoxArray& rhs) const noexcept;
603 
605  [[nodiscard]] bool operator!= (const BoxArray& rhs) const noexcept;
606 
607  [[nodiscard]] bool operator== (const Vector<Box>& bv) const noexcept;
608  [[nodiscard]] bool operator!= (const Vector<Box>& bv) const noexcept;
609 
611  [[nodiscard]] bool CellEqual (const BoxArray& rhs) const noexcept;
612 
614  BoxArray& maxSize (int block_size);
615 
616  BoxArray& maxSize (const IntVect& block_size);
617 
621  BoxArray& minmaxSize (const IntVect& min_size, const IntVect& max_size);
622 
624  BoxArray& refine (int refinement_ratio);
625 
627  BoxArray& refine (const IntVect& iv);
628 
630  BoxArray& coarsen (int refinement_ratio);
631 
633  [[nodiscard]] bool coarsenable (int refinement_ratio, int min_width=1) const;
634  [[nodiscard]] bool coarsenable (const IntVect& refinement_ratio, int min_width=1) const;
635  [[nodiscard]] bool coarsenable (const IntVect& refinement_ratio, const IntVect& min_width) const;
636 
638  BoxArray& coarsen (const IntVect& iv);
639 
641  BoxArray& growcoarsen (int n, const IntVect& iv);
642  BoxArray& growcoarsen (IntVect const& ngrow, const IntVect& iv);
643 
645  BoxArray& grow (int n);
646 
648  BoxArray& grow (const IntVect& iv);
653  BoxArray& grow (int idir, int n_cell);
658  BoxArray& growLo (int idir, int n_cell);
663  BoxArray& growHi (int idir, int n_cell);
673  BoxArray& surroundingNodes (int dir);
674 
677 
679  BoxArray& enclosedCells (int dir);
680 
682  BoxArray& convert (IndexType typ);
683 
684  BoxArray& convert (const IntVect& iv);
685 
687  BoxArray& convert (Box (*fp)(const Box&));
688 
690  BoxArray& shift (int dir, int nzones);
691 
693  BoxArray& shift (const IntVect &iv);
694 
696  void set (int i, const Box& ibox);
697 
699  [[nodiscard]] Box operator[] (int index) const noexcept {
700  return m_bat(m_ref->m_abox[index]);
701  }
702 
704  [[nodiscard]] Box operator[] (const MFIter& mfi) const noexcept;
705 
707  [[nodiscard]] Box get (int index) const noexcept { return operator[](index); }
708 
710  [[nodiscard]] Box getCellCenteredBox (int index) const noexcept {
711  return m_bat.coarsen(m_ref->m_abox[index]);
712  }
713 
718  [[nodiscard]] bool ok () const;
719 
721  [[nodiscard]] bool isDisjoint () const;
722 
724  [[nodiscard]] BoxList boxList () const;
725 
727  [[nodiscard]] bool contains (const IntVect& v) const;
728 
733  [[nodiscard]]
734  bool contains (const Box& b, bool assume_disjoint_ba = false,
735  const IntVect& ng = IntVect(0)) const;
736 
740  [[nodiscard]]
741  bool contains (const BoxArray& ba, bool assume_disjoint_ba = false,
742  const IntVect& ng = IntVect(0)) const;
743 
751  [[nodiscard]]
752  bool contains (const BoxArray& ba, Periodicity const& period) const;
753 
755  [[nodiscard]] Box minimalBox () const;
756  [[nodiscard]] Box minimalBox (Long& npts_avg_box) const;
757 
762  [[nodiscard]]
763  bool intersects (const Box& b, int ng = 0) const;
764 
765  [[nodiscard]]
766  bool intersects (const Box& b, const IntVect& ng) const;
767 
769  [[nodiscard]]
770  std::vector< std::pair<int,Box> > intersections (const Box& bx) const;
771 
773  [[nodiscard]]
774  std::vector< std::pair<int,Box> > intersections (const Box& bx, bool first_only, int ng) const;
775 
776  [[nodiscard]]
777  std::vector< std::pair<int,Box> > intersections (const Box& bx, bool first_only, const IntVect& ng) const;
778 
780  void intersections (const Box& bx, std::vector< std::pair<int,Box> >& isects) const;
781 
783  void intersections (const Box& bx, std::vector< std::pair<int,Box> >& isects,
784  bool first_only, int ng) const;
785 
786  void intersections (const Box& bx, std::vector< std::pair<int,Box> >& isects,
787  bool first_only, const IntVect& ng) const;
788 
790  [[nodiscard]] BoxList complementIn (const Box& b) const;
791  void complementIn (BoxList& bl, const Box& b) const;
792 
794  void clear_hash_bin () const;
795 
797  void removeOverlap (bool simplify=true);
798 
800  [[nodiscard]] static bool SameRefs (const BoxArray& lhs, const BoxArray& rhs) { return lhs.m_ref == rhs.m_ref; }
801 
802  struct RefID {
803  RefID () noexcept = default;
804  explicit RefID (BARef* data_) noexcept : data(data_) {}
805  bool operator< (const RefID& rhs) const noexcept { return std::less<>()(data,rhs.data); }
806  bool operator== (const RefID& rhs) const noexcept { return data == rhs.data; }
807  bool operator!= (const RefID& rhs) const noexcept { return data != rhs.data; }
808  friend std::ostream& operator<< (std::ostream& os, const RefID& id);
809  private:
810  BARef* data{nullptr};
811  };
812 
814  [[nodiscard]] RefID getRefID () const noexcept { return RefID { m_ref.get() }; }
815 
817  [[nodiscard]] IndexType ixType () const noexcept { return m_bat.index_type(); }
818 
820  [[nodiscard]] IntVect crseRatio () const noexcept { return m_bat.coarsen_ratio(); }
821 
822  static void Initialize ();
823  static void Finalize ();
824  static bool initialized;
825 
827  void uniqify ();
828 
829  [[nodiscard]] BoxList const& simplified_list () const; // For regular AMR grids only!!!
830  [[nodiscard]] BoxArray simplified () const; // For regular AMR grids only!!!
831 
832  [[nodiscard]] BATransformer const& transformer () const;
833 
834  friend class AmrMesh;
835  friend class FabArrayBase;
836 
837 private:
839  void type_update ();
840 
841  [[nodiscard]] BARef::HashType& getHashMap () const;
842 
843  [[nodiscard]] IntVect getDoiLo () const noexcept;
844  [[nodiscard]] IntVect getDoiHi () const noexcept;
845 
848  std::shared_ptr<BARef> m_ref;
849  mutable std::shared_ptr<BoxList> m_simplified_list;
850 };
851 
853 std::ostream& operator<< (std::ostream& os, const BoxArray& ba);
854 
855 std::ostream& operator<< (std::ostream& os, const BoxArray::RefID& id);
856 
857 }
858 
859 #endif /*BL_BOXARRAY_H*/
int idir
Definition: AMReX_HypreMLABecLap.cpp:1093
Definition: AMReX_AmrMesh.H:62
A collection of Boxes stored in an Array.
Definition: AMReX_BoxArray.H:530
IndexType ixType() const noexcept
Return index type of this BoxArray.
Definition: AMReX_BoxArray.H:817
Box operator[](int index) const noexcept
Return element index of this BoxArray.
Definition: AMReX_BoxArray.H:699
Long capacity() const noexcept
Return the number of boxes that can be held in the current allocated storage.
Definition: AMReX_BoxArray.H:580
RefID getRefID() const noexcept
Return a unique ID of the reference.
Definition: AMReX_BoxArray.H:814
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:1159
BoxArray & grow(int n)
Grow each Box in the BoxArray by the specified amount.
Definition: AMReX_BoxArray.cpp:696
static void Finalize()
Definition: AMReX_BoxArray.cpp:272
BoxArray() noexcept
Definition: AMReX_BoxArray.cpp:277
std::vector< std::pair< int, Box > > intersections(const Box &bx) const
Return intersections of Box and BoxArray.
Definition: AMReX_BoxArray.cpp:1176
std::ostream & writeOn(std::ostream &) const
Output this BoxArray to a checkpoint file.
Definition: AMReX_BoxArray.cpp:472
bool contains(const IntVect &v) const
True if the IntVect is within any of the Boxes in this BoxArray.
Definition: AMReX_BoxArray.cpp:967
static void Initialize()
Definition: AMReX_BoxArray.cpp:261
BATransformer m_bat
Definition: AMReX_BoxArray.H:846
BoxList const & simplified_list() const
Definition: AMReX_BoxArray.cpp:1601
BoxList boxList() const
Create a BoxList from this BoxArray.
Definition: AMReX_BoxArray.cpp:939
BoxArray simplified() const
Definition: AMReX_BoxArray.cpp:1612
IntVect crseRatio() const noexcept
Return crse ratio of this BoxArray.
Definition: AMReX_BoxArray.H:820
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:458
IntVect getDoiLo() const noexcept
Definition: AMReX_BoxArray.cpp:1512
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:342
Long numPts() const noexcept
Returns the total number of cells contained in all boxes in the BoxArray.
Definition: AMReX_BoxArray.cpp:384
std::shared_ptr< BoxList > m_simplified_list
Definition: AMReX_BoxArray.H:849
BARef::HashType & getHashMap() const
Definition: AMReX_BoxArray.cpp:1524
IntVect getDoiHi() const noexcept
Definition: AMReX_BoxArray.cpp:1518
void clear_hash_bin() const
Clear out the internal hash table used by intersections.
Definition: AMReX_BoxArray.cpp:1403
Box minimalBox() const
Return smallest Box that contains all Boxes in this BoxArray.
Definition: AMReX_BoxArray.cpp:1057
BoxArray & minmaxSize(const IntVect &min_size, const IntVect &max_size)
Definition: AMReX_BoxArray.cpp:567
BoxArray & refine(int refinement_ratio)
Refine each Box in the BoxArray to the specified ratio.
Definition: AMReX_BoxArray.cpp:583
BoxArray & operator=(BoxArray const &rhs)=default
bool coarsenable(int refinement_ratio, int min_width=1) const
Coarsen each Box in the BoxArray to the specified ratio.
Definition: AMReX_BoxArray.cpp:605
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:742
void resize(Long len)
Resize the BoxArray. See Vector<T>::resize() for the gory details.
Definition: AMReX_BoxArray.cpp:377
BoxArray & enclosedCells()
Apply Box::enclosedCells() to each Box in the BoxArray.
Definition: AMReX_BoxArray.cpp:789
Box getCellCenteredBox(int index) const noexcept
Return cell-centered box at element index of this BoxArray.
Definition: AMReX_BoxArray.H:710
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:421
BATransformer const & transformer() const
Definition: AMReX_BoxArray.cpp:1618
static bool initialized
Definition: AMReX_BoxArray.H:824
static bool SameRefs(const BoxArray &lhs, const BoxArray &rhs)
whether two BoxArrays share the same data
Definition: AMReX_BoxArray.H:800
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:758
BoxList complementIn(const Box &b) const
Return box - boxarray.
Definition: AMReX_BoxArray.cpp:1304
void type_update()
Update BoxArray index type according the box type, and then convert boxes to cell-centered.
Definition: AMReX_BoxArray.cpp:1498
std::shared_ptr< BARef > m_ref
The data – a reference-counted pointer to a Ref.
Definition: AMReX_BoxArray.H:848
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:1419
BoxArray & coarsen(int refinement_ratio)
Coarsen each Box in the BoxArray to the specified ratio.
Definition: AMReX_BoxArray.cpp:662
BoxArray & shift(int dir, int nzones)
Apply Box::shift(int,int) to each Box in the BoxArray.
Definition: AMReX_BoxArray.cpp:837
BoxArray & surroundingNodes()
Apply surroundingNodes(Box) to each Box in BoxArray. See the documentation of Box for details.
Definition: AMReX_BoxArray.cpp:774
BoxArray(BoxArray &&rhs) noexcept=default
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition: AMReX_BoxArray.H:577
Box get(int index) const noexcept
Return element index of this BoxArray.
Definition: AMReX_BoxArray.H:707
~BoxArray() noexcept=default
bool CellEqual(const BoxArray &rhs) const noexcept
Are the BoxArrays equal after conversion to cell-centered.
Definition: AMReX_BoxArray.cpp:536
void set(int i, const Box &ibox)
Set element i in this BoxArray to Box ibox.
Definition: AMReX_BoxArray.cpp:868
BoxArray & growcoarsen(int n, const IntVect &iv)
Grow and then coarsen each Box in the BoxArray.
Definition: AMReX_BoxArray.cpp:675
BoxArray & convert(IndexType typ)
Apply Box::convert(IndexType) to each Box in the BoxArray.
Definition: AMReX_BoxArray.cpp:803
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:884
void uniqify()
Make ourselves unique.
Definition: AMReX_BoxArray.cpp:1578
BoxArray(const BoxArray &rhs)=default
bool empty() const noexcept
Return whether the BoxArray is empty.
Definition: AMReX_BoxArray.H:583
bool isDisjoint() const
Return true if set of intersecting Boxes in BoxArray is null.
Definition: AMReX_BoxArray.cpp:910
BoxArray & maxSize(int block_size)
Forces each Box in BoxArray to have sides <= block_size.
Definition: AMReX_BoxArray.cpp:543
void clear()
Remove all Boxes from the BoxArray.
Definition: AMReX_BoxArray.cpp:369
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:41
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE CellIndex ixType(int dir) const noexcept
Returns the CellIndex in direction dir.
Definition: AMReX_IndexType.H:116
AMREX_GPU_HOST_DEVICE static constexpr AMREX_FORCE_INLINE IntVectND< dim > TheUnitVector() noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition: AMReX_IntVect.H:682
AMREX_GPU_HOST_DEVICE static constexpr AMREX_FORCE_INLINE IntVectND< dim > TheZeroVector() noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition: AMReX_IntVect.H:672
Definition: AMReX_MFIter.H:57
Encapsulation of the Orientation of the Faces of a Box.
Definition: AMReX_Orientation.H:29
AMREX_GPU_HOST_DEVICE int coordDir() const noexcept
Returns the coordinate direction.
Definition: AMReX_Orientation.H:83
AMREX_GPU_HOST_DEVICE bool isLow() const noexcept
Returns true if Orientation is low.
Definition: AMReX_Orientation.H:89
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:27
AMREX_EXPORT int max_grid_size
Definition: AMReX_EB2.cpp:23
AMREX_GPU_HOST_DEVICE Long size(T const &b) noexcept
integer version
Definition: AMReX_GpuRange.H:26
static int f(amrex::Real t, N_Vector y_data, N_Vector y_rhs, void *user_data)
Definition: AMReX_SundialsIntegrator.H:44
Definition: AMReX_Amr.cpp:49
BoxND< AMREX_SPACEDIM > Box
Definition: AMReX_BaseFwd.H:27
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:1664
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Returns a BoxND with different type.
Definition: AMReX_Box.H:1435
BoxList GetBndryCells(const BoxArray &ba, int ngrow)
Find the ghost cells of a given BoxArray.
Definition: AMReX_BoxArray.cpp:1789
IntVectND< AMREX_SPACEDIM > IntVect
Definition: AMReX_BaseFwd.H:30
bool match(const BoxArray &x, const BoxArray &y)
Note that two BoxArrays that match are not necessarily equal.
Definition: AMReX_BoxArray.cpp:1877
BoxArray complementIn(const Box &b, const BoxArray &ba)
Make a BoxArray from the complement of BoxArray ba in Box b.
Definition: AMReX_BoxArray.cpp:1657
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > refine(const BoxND< dim > &b, int ref_ratio) noexcept
Definition: AMReX_Box.H:1342
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo/rati...
Definition: AMReX_Box.H:1304
IndexTypeND< AMREX_SPACEDIM > IndexType
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:1847
BoxArray boxComplement(const Box &b1in, const Box &b2)
Make a BoxArray from the the complement of b2 in b1in.
Definition: AMReX_BoxArray.cpp:1650
Definition: AMReX_BoxArray.H:57
static Long total_box_bytes_hwm
Definition: AMReX_BoxArray.H:110
~BARef()
Definition: AMReX_BoxArray.cpp:85
BARef & operator=(const BARef &rhs)=delete
static int numboxarrays
Definition: AMReX_BoxArray.H:107
static Long total_hash_bytes
Definition: AMReX_BoxArray.H:111
BARef(BARef &&rhs)=delete
IntVect crsn
Definition: AMReX_BoxArray.H:98
std::unordered_map< IntVect, std::vector< int >, IntVect::shift_hasher > HashType
Definition: AMReX_BoxArray.H:100
static Long total_hash_bytes_hwm
Definition: AMReX_BoxArray.H:112
void resize(Long n)
Definition: AMReX_BoxArray.cpp:180
static bool initialized
Definition: AMReX_BoxArray.H:116
HashType hash
Definition: AMReX_BoxArray.H:103
BARef()
Definition: AMReX_BoxArray.cpp:35
Box bbox
Box hash stuff.
Definition: AMReX_BoxArray.H:96
static int numboxarrays_hwm
Definition: AMReX_BoxArray.H:108
bool has_hashmap
Definition: AMReX_BoxArray.H:105
Vector< Box > m_abox
The data.
Definition: AMReX_BoxArray.H:93
static Long total_box_bytes
Definition: AMReX_BoxArray.H:109
static void Initialize()
Definition: AMReX_BoxArray.cpp:231
static void Finalize()
Definition: AMReX_BoxArray.cpp:255
void define(const Box &bx)
Definition: AMReX_BoxArray.cpp:143
bool HasHashMap() const
Definition: AMReX_BoxArray.H:82
Definition: AMReX_BoxArray.H:175
IndexType index_type() const noexcept
Definition: AMReX_BoxArray.H:237
Box operator()(const Box &a_bx) const noexcept
Definition: AMReX_BoxArray.H:218
BATbndryReg(Orientation a_face, IndexType a_typ, int a_in_rad, int a_out_rad, int a_extent_rad)
Definition: AMReX_BoxArray.H:176
friend bool operator==(BATbndryReg const &a, BATbndryReg const &b) noexcept
Definition: AMReX_BoxArray.H:240
IntVect doiHi() const noexcept
Definition: AMReX_BoxArray.H:235
IntVect coarsen_ratio() const noexcept
Definition: AMReX_BoxArray.H:238
IntVect doiLo() const noexcept
Definition: AMReX_BoxArray.H:234
Orientation m_face
Definition: AMReX_BoxArray.H:246
IntVect m_loshft
Definition: AMReX_BoxArray.H:249
IndexType m_typ
Definition: AMReX_BoxArray.H:247
Box coarsen(Box const &a_box) const noexcept
Definition: AMReX_BoxArray.H:232
IntVect m_doihi
Definition: AMReX_BoxArray.H:252
IntVect m_hishft
Definition: AMReX_BoxArray.H:250
IntVect m_doilo
Definition: AMReX_BoxArray.H:251
IntVect m_crse_ratio
Definition: AMReX_BoxArray.H:248
Definition: AMReX_BoxArray.H:142
BATcoarsenRatio(IntVect const &a_crse_ratio)
Definition: AMReX_BoxArray.H:143
static constexpr IntVect doiHi()
Definition: AMReX_BoxArray.H:147
Box coarsen(Box const &a_box) const noexcept
Definition: AMReX_BoxArray.H:145
static constexpr IndexType index_type()
Definition: AMReX_BoxArray.H:148
static constexpr IntVect doiLo()
Definition: AMReX_BoxArray.H:146
IntVect m_crse_ratio
Definition: AMReX_BoxArray.H:150
Box operator()(const Box &bx) const noexcept
Definition: AMReX_BoxArray.H:144
IntVect coarsen_ratio() const noexcept
Definition: AMReX_BoxArray.H:149
Definition: AMReX_BoxArray.H:154
BATindexType_coarsenRatio(IndexType a_typ, IntVect const &a_crse_ratio)
Definition: AMReX_BoxArray.H:155
static constexpr IntVect doiLo()
Definition: AMReX_BoxArray.H:164
IntVect doiHi() const noexcept
Definition: AMReX_BoxArray.H:165
Box operator()(const Box &bx) const noexcept
Definition: AMReX_BoxArray.H:158
IntVect m_crse_ratio
Definition: AMReX_BoxArray.H:171
IntVect coarsen_ratio() const noexcept
Definition: AMReX_BoxArray.H:168
IndexType index_type() const noexcept
Definition: AMReX_BoxArray.H:167
Box coarsen(Box const &a_box) const noexcept
Definition: AMReX_BoxArray.H:162
IndexType m_typ
Definition: AMReX_BoxArray.H:170
Definition: AMReX_BoxArray.H:130
IntVect doiHi() const noexcept
Definition: AMReX_BoxArray.H:135
Box operator()(const Box &bx) const noexcept
Definition: AMReX_BoxArray.H:132
IndexType index_type() const noexcept
Definition: AMReX_BoxArray.H:136
static Box coarsen(Box const &a_box) noexcept
Definition: AMReX_BoxArray.H:133
BATindexType(IndexType a_typ)
Definition: AMReX_BoxArray.H:131
IndexType m_typ
Definition: AMReX_BoxArray.H:138
static constexpr IntVect coarsen_ratio()
Definition: AMReX_BoxArray.H:137
static constexpr IntVect doiLo()
Definition: AMReX_BoxArray.H:134
Definition: AMReX_BoxArray.H:120
static constexpr Box coarsen(Box const &a_box)
Definition: AMReX_BoxArray.H:122
static constexpr IndexType index_type()
Definition: AMReX_BoxArray.H:125
static constexpr IntVect doiHi()
Definition: AMReX_BoxArray.H:124
Box operator()(const Box &bx) const noexcept
Definition: AMReX_BoxArray.H:121
static constexpr IntVect coarsen_ratio()
Definition: AMReX_BoxArray.H:126
static constexpr IntVect doiLo()
Definition: AMReX_BoxArray.H:123
Definition: AMReX_BoxArray.H:256
BATransformer(IndexType t)
Definition: AMReX_BoxArray.H:279
void set_index_type(IndexType typ) noexcept
Definition: AMReX_BoxArray.H:446
BATransformer(Orientation f, IndexType t, int in_rad, int out_rad, int extent_rad)
Definition: AMReX_BoxArray.H:283
void set_coarsen_ratio(IntVect const &a_ratio) noexcept
Definition: AMReX_BoxArray.H:391
IntVect doiLo() const noexcept
Definition: AMReX_BoxArray.H:319
IntVect coarsen_ratio() const noexcept
Definition: AMReX_BoxArray.H:367
bool is_null() const noexcept
Definition: AMReX_BoxArray.H:383
Box operator()(Box const &ab) const noexcept
Definition: AMReX_BoxArray.H:287
BATType
Definition: AMReX_BoxArray.H:257
Box coarsen(Box const &a_box) const noexcept
Definition: AMReX_BoxArray.H:303
friend bool operator==(BATransformer const &a, BATransformer const &b) noexcept
Definition: AMReX_BoxArray.H:501
bool is_simple() const noexcept
Definition: AMReX_BoxArray.H:387
IndexType index_type() const noexcept
Definition: AMReX_BoxArray.H:351
IntVect doiHi() const noexcept
Definition: AMReX_BoxArray.H:335
BATType m_bat_type
Definition: AMReX_BoxArray.H:512
BATOp m_op
Definition: AMReX_BoxArray.H:513
Definition: AMReX_BoxArray.H:802
friend std::ostream & operator<<(std::ostream &os, const RefID &id)
Definition: AMReX_BoxArray.cpp:1891
RefID() noexcept=default
bool operator<(const RefID &rhs) const noexcept
Definition: AMReX_BoxArray.H:805
bool operator!=(const RefID &rhs) const noexcept
Definition: AMReX_BoxArray.H:807
bool operator==(const RefID &rhs) const noexcept
Definition: AMReX_BoxArray.H:806
BARef * data
Definition: AMReX_BoxArray.H:810
Definition: AMReX_BoxArray.H:259
BATindexType m_indexType
Definition: AMReX_BoxArray.H:271
BATOp(IndexType t) noexcept
Definition: AMReX_BoxArray.H:262
BATOp() noexcept
Definition: AMReX_BoxArray.H:260
BATOp(IntVect const &r) noexcept
Definition: AMReX_BoxArray.H:264
BATnull m_null
Definition: AMReX_BoxArray.H:270
BATbndryReg m_bndryReg
Definition: AMReX_BoxArray.H:274
BATindexType_coarsenRatio m_indexType_coarsenRatio
Definition: AMReX_BoxArray.H:273
BATOp(Orientation f, IndexType t, int in_rad, int out_rad, int extent_rad) noexcept
Definition: AMReX_BoxArray.H:268
BATcoarsenRatio m_coarsenRatio
Definition: AMReX_BoxArray.H:272
BATOp(IndexType t, IntVect const &r) noexcept
Definition: AMReX_BoxArray.H:266