Line data Source code
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 :
21 : //! Make a BoxArray from the the complement of b2 in b1in.
22 : [[nodiscard]] BoxArray boxComplement (const Box& b1in, const Box& b2);
23 :
24 : //! Make a BoxArray from the complement of BoxArray ba in Box b.
25 : [[nodiscard]] BoxArray complementIn (const Box& b, const BoxArray& ba);
26 :
27 : //! Make a BoxArray from the intersection of Box b and BoxArray(+ghostcells).
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 :
32 : //! Make a BoxArray from the intersection of two BoxArrays.
33 : [[nodiscard]] BoxArray intersect (const BoxArray& lhs, const BoxArray& rhs);
34 :
35 : //! Make a BoxList from the intersection of BoxArray and BoxList.
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 :
47 : //! Find the ghost cells of a given BoxArray.
48 : [[nodiscard]] BoxList GetBndryCells (const BoxArray& ba, int ngrow);
49 :
50 : //! Read a BoxArray from a stream. If b is true, read in a special way
51 : void readBoxArray (BoxArray& ba, std::istream& s, bool b = false);
52 :
53 : //! Note that two BoxArrays that match are not necessarily equal.
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);
75 : //!
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 : //
92 : //! The data.
93 : Vector<Box> m_abox;
94 : //
95 : //! Box hash stuff.
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 6842 : [[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 1836 : [[nodiscard]] static constexpr IndexType index_type () { return IndexType(); }
126 : [[nodiscard]] static constexpr IntVect coarsen_ratio () { return IntVect::TheUnitVector(); }
127 : };
128 :
129 : struct BATindexType
130 : {
131 : explicit BATindexType (IndexType a_typ) : m_typ(a_typ) {}
132 329028 : [[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 1500 : [[nodiscard]] IndexType index_type () const noexcept { return m_typ; }
137 : [[nodiscard]] static constexpr IntVect coarsen_ratio () { return IntVect::TheUnitVector(); }
138 : IndexType m_typ;
139 : };
140 :
141 : struct BATcoarsenRatio
142 : {
143 : explicit BATcoarsenRatio (IntVect const& a_crse_ratio) : m_crse_ratio(a_crse_ratio) {}
144 8412 : [[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 1234 : [[nodiscard]] static constexpr IndexType index_type () { return IndexType(); }
149 : [[nodiscard]] IntVect coarsen_ratio () const noexcept { return m_crse_ratio; }
150 : IntVect m_crse_ratio;
151 : };
152 :
153 : struct BATindexType_coarsenRatio
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 488162 : [[nodiscard]] Box operator() (const Box& bx) const noexcept {
159 976324 : return amrex::convert(amrex::coarsen(bx,m_crse_ratio),m_typ);
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 2400 : [[nodiscard]] IndexType index_type () const noexcept { return m_typ; }
168 : [[nodiscard]] IntVect coarsen_ratio () const noexcept { return m_crse_ratio; }
169 :
170 : IndexType m_typ;
171 : IntVect m_crse_ratio;
172 : };
173 :
174 : struct BATbndryReg
175 : {
176 : BATbndryReg (Orientation a_face, IndexType a_typ,
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 0 : [[nodiscard]] Box operator() (const Box& a_bx) const noexcept {
219 0 : IntVect lo = amrex::coarsen(a_bx.smallEnd(), m_crse_ratio);
220 0 : IntVect hi = amrex::coarsen(a_bx.bigEnd(), m_crse_ratio);
221 0 : const int d = m_face.coordDir();
222 0 : if (m_face.isLow()) {
223 0 : hi[d] = lo[d];
224 : } else {
225 0 : lo[d] = hi[d];
226 : }
227 0 : lo += m_loshft;
228 0 : hi += m_hishft;
229 0 : 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 0 : [[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 :
246 : Orientation m_face;
247 : IndexType m_typ;
248 : IntVect m_crse_ratio;
249 : IntVect m_loshft;
250 : IntVect m_hishft;
251 : IntVect m_doilo;
252 : IntVect m_doihi;
253 : };
254 :
255 : struct BATransformer
256 : {
257 : enum struct BATType { null, indexType, coarsenRatio, indexType_coarsenRatio, bndryReg};
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) {}
270 : BATnull m_null;
271 : BATindexType m_indexType;
272 : BATcoarsenRatio m_coarsenRatio;
273 : BATindexType_coarsenRatio m_indexType_coarsenRatio;
274 : BATbndryReg m_bndryReg;
275 : };
276 :
277 : BATransformer () = default;
278 :
279 : BATransformer (IndexType t)
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 663724 : [[nodiscard]] Box operator() (Box const& ab) const noexcept {
288 663724 : switch (m_bat_type)
289 : {
290 6842 : case BATType::null:
291 6842 : return m_op.m_null(ab);
292 164514 : case BATType::indexType:
293 164514 : return m_op.m_indexType(ab);
294 4206 : case BATType::coarsenRatio:
295 4206 : return m_op.m_coarsenRatio(ab);
296 488162 : case BATType::indexType_coarsenRatio:
297 488162 : return m_op.m_indexType_coarsenRatio(ab);
298 0 : default:
299 0 : 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);
310 : case BATType::coarsenRatio:
311 : return m_op.m_coarsenRatio.coarsen(a_box);
312 : case BATType::indexType_coarsenRatio:
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:
325 : return amrex::BATindexType::doiLo();
326 : case BATType::coarsenRatio:
327 : return amrex::BATcoarsenRatio::doiLo();
328 : case BATType::indexType_coarsenRatio:
329 : return amrex::BATindexType_coarsenRatio::doiLo();
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();
342 : case BATType::coarsenRatio:
343 : return amrex::BATcoarsenRatio::doiHi();
344 : case BATType::indexType_coarsenRatio:
345 : return m_op.m_indexType_coarsenRatio.doiHi();
346 : default:
347 : return m_op.m_bndryReg.doiHi();
348 : }
349 : }
350 :
351 6970 : [[nodiscard]] IndexType index_type () const noexcept {
352 6970 : switch (m_bat_type)
353 : {
354 1836 : case BATType::null:
355 1836 : return amrex::BATnull::index_type();
356 1500 : case BATType::indexType:
357 1500 : return m_op.m_indexType.index_type();
358 1234 : case BATType::coarsenRatio:
359 1234 : return amrex::BATcoarsenRatio::index_type();
360 2400 : case BATType::indexType_coarsenRatio:
361 2400 : return m_op.m_indexType_coarsenRatio.index_type();
362 0 : default:
363 0 : 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:
371 : return amrex::BATnull::coarsen_ratio();
372 : case BATType::indexType:
373 : return amrex::BATindexType::coarsen_ratio();
374 : case BATType::coarsenRatio:
375 : return m_op.m_coarsenRatio.coarsen_ratio();
376 : case BATType::indexType_coarsenRatio:
377 : return m_op.m_indexType_coarsenRatio.coarsen_ratio();
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 {
399 : m_bat_type = BATType::coarsenRatio;
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 {
409 : m_bat_type = BATType::indexType_coarsenRatio;
410 : auto t = m_op.m_indexType.m_typ;
411 : m_op.m_indexType_coarsenRatio.m_typ = t;
412 : m_op.m_indexType_coarsenRatio.m_crse_ratio = a_ratio;
413 : return;
414 : }
415 : }
416 : case BATType::coarsenRatio:
417 : {
418 : if (a_ratio == IntVect::TheUnitVector()) {
419 : m_bat_type = BATType::null;
420 : return;
421 : } else {
422 : m_op.m_coarsenRatio.m_crse_ratio = a_ratio;
423 : return;
424 : }
425 : }
426 : case BATType::indexType_coarsenRatio:
427 : {
428 : if (a_ratio == IntVect::TheUnitVector()) {
429 : m_bat_type = BATType::indexType;
430 : auto t = m_op.m_indexType_coarsenRatio.m_typ;
431 : m_op.m_indexType.m_typ = t;
432 : return;
433 : } else {
434 : m_op.m_indexType_coarsenRatio.m_crse_ratio = a_ratio;
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 {
454 : m_bat_type = BATType::indexType;
455 : m_op.m_indexType.m_typ = typ;
456 : return;
457 : }
458 : }
459 : case BATType::indexType:
460 : {
461 : if (typ.cellCentered()) {
462 : m_bat_type = BATType::null;
463 : return;
464 : } else {
465 : m_op.m_indexType.m_typ = typ;
466 : return;
467 : }
468 : }
469 : case BATType::coarsenRatio:
470 : {
471 : if (typ.cellCentered()) {
472 : return;
473 : } else {
474 : m_bat_type = BATType::indexType_coarsenRatio;
475 : auto r = m_op.m_coarsenRatio.m_crse_ratio;
476 : m_op.m_indexType_coarsenRatio.m_typ = typ;
477 : m_op.m_indexType_coarsenRatio.m_crse_ratio = r;
478 : return;
479 : }
480 : }
481 : case BATType::indexType_coarsenRatio:
482 : {
483 : if (typ.cellCentered()) {
484 : m_bat_type = BATType::coarsenRatio;
485 : auto r = m_op.m_indexType_coarsenRatio.m_crse_ratio;
486 : m_op.m_coarsenRatio.m_crse_ratio = r;
487 : return;
488 : } else {
489 : m_op.m_indexType_coarsenRatio.m_typ = typ;
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 :
512 : BATType m_bat_type{BATType::null};
513 : BATOp m_op;
514 : };
515 :
516 : // for backward compatibility
517 : using BndryBATransformer = BATransformer;
518 :
519 : class MFIter;
520 : class AmrMesh;
521 : class FabArrayBase;
522 :
523 : /**
524 : * \brief A collection of Boxes stored in an Array.
525 : *
526 : * It is a reference-counted concrete class, not a polymorphic one; i.e. you
527 : * cannot use any of the List member functions with a BoxList.
528 : */
529 11120 : class BoxArray
530 : {
531 : public:
532 :
533 : BoxArray () noexcept;
534 6366 : BoxArray (const BoxArray& rhs) = default;
535 0 : BoxArray (BoxArray&& rhs) noexcept = default;
536 : BoxArray& operator= (BoxArray const& rhs) = default;
537 : BoxArray& operator= (BoxArray&& rhs) noexcept = default;
538 515104 : ~BoxArray() noexcept = default;
539 :
540 : //! Make a boxarray out of a single box
541 : explicit BoxArray (const Box& bx);
542 :
543 : //! Construct a BoxArray of the specified size.
544 : explicit BoxArray (size_t n);
545 :
546 : //! Construct a BoxArray from an array of Boxes of size nbox.
547 : BoxArray (const Box* bxvec,
548 : int nbox);
549 :
550 : //! Construct a BoxArray from a BoxList.
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 :
558 : /**
559 : * \brief Initialize the BoxArray from a single box.
560 : * It is an error if the BoxArray has already been initialized.
561 : */
562 : void define (const Box& bx);
563 : /**
564 : * \brief Initialize the BoxArray from the supplied BoxList.
565 : * It is an error if the BoxArray has already been initialized.
566 : */
567 : void define (const BoxList& bl);
568 : void define (BoxList&& bl) noexcept;
569 :
570 : //! Remove all Boxes from the BoxArray.
571 : void clear ();
572 :
573 : //! Resize the BoxArray. See Vector<T>::resize() for the gory details.
574 : void resize (Long len);
575 :
576 : //! Return the number of boxes in the BoxArray.
577 300 : [[nodiscard]] Long size () const noexcept { return m_ref->m_abox.size(); }
578 :
579 : //! Return the number of boxes that can be held in the current allocated storage
580 : [[nodiscard]] Long capacity () const noexcept { return static_cast<Long>(m_ref->m_abox.capacity()); }
581 :
582 : //! Return whether the BoxArray is empty
583 13018 : [[nodiscard]] bool empty () const noexcept { return m_ref->m_abox.empty(); }
584 :
585 : //! Returns the total number of cells contained in all boxes in the BoxArray.
586 : [[nodiscard]] Long numPts() const noexcept;
587 :
588 : //! Returns the total number of cells (in double type) contained in all boxes in the BoxArray.
589 : [[nodiscard]] double d_numPts () const noexcept;
590 : /**
591 : * \brief Initialize the BoxArray from the supplied istream.
592 : * It is an error if the BoxArray has already been initialized.
593 : * Note that the BoxArray in the istream must have been written
594 : * using writeOn().
595 : */
596 : int readFrom (std::istream& is);
597 :
598 : //! Output this BoxArray to a checkpoint file.
599 : std::ostream& writeOn (std::ostream&) const;
600 :
601 : //! Are the BoxArrays equal?
602 : [[nodiscard]] bool operator== (const BoxArray& rhs) const noexcept;
603 :
604 : //! Are the BoxArrays not equal?
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 :
610 : //! Are the BoxArrays equal after conversion to cell-centered
611 : [[nodiscard]] bool CellEqual (const BoxArray& rhs) const noexcept;
612 :
613 : //! Forces each Box in BoxArray to have sides <= block_size.
614 : BoxArray& maxSize (int block_size);
615 :
616 : BoxArray& maxSize (const IntVect& block_size);
617 :
618 : //! Refine each Box in the BoxArray to the specified ratio.
619 : BoxArray& refine (int refinement_ratio);
620 :
621 : //! Refine each Box in the BoxArray to the specified ratio.
622 : BoxArray& refine (const IntVect& iv);
623 :
624 : //! Coarsen each Box in the BoxArray to the specified ratio.
625 : BoxArray& coarsen (int refinement_ratio);
626 :
627 : //! Coarsen each Box in the BoxArray to the specified ratio.
628 : [[nodiscard]] bool coarsenable (int refinement_ratio, int min_width=1) const;
629 : [[nodiscard]] bool coarsenable (const IntVect& refinement_ratio, int min_width=1) const;
630 : [[nodiscard]] bool coarsenable (const IntVect& refinement_ratio, const IntVect& min_width) const;
631 :
632 : //! Coarsen each Box in the BoxArray to the specified ratio.
633 : BoxArray& coarsen (const IntVect& iv);
634 :
635 : //! Grow and then coarsen each Box in the BoxArray.
636 : BoxArray& growcoarsen (int n, const IntVect& iv);
637 : BoxArray& growcoarsen (IntVect const& ngrow, const IntVect& iv);
638 :
639 : //! Grow each Box in the BoxArray by the specified amount.
640 : BoxArray& grow (int n);
641 :
642 : //! Grow each Box in the BoxArray by the specified amount.
643 : BoxArray& grow (const IntVect& iv);
644 : /**
645 : * \brief Grow each Box in the BoxArray on the low and high ends
646 : * by n_cell cells in the idir direction.
647 : */
648 : BoxArray& grow (int idir, int n_cell);
649 : /**
650 : * \brief Grow each Box in the BoxArray on the low end
651 : * by n_cell cells in the idir direction.
652 : */
653 : BoxArray& growLo (int idir, int n_cell);
654 : /**
655 : * \brief Grow each Box in the BoxArray on the high end
656 : * by n_cell cells in the idir direction.
657 : */
658 : BoxArray& growHi (int idir, int n_cell);
659 : /**
660 : * \brief Apply surroundingNodes(Box) to each Box in BoxArray.
661 : * See the documentation of Box for details.
662 : */
663 : BoxArray& surroundingNodes ();
664 : /**
665 : * \brief Apply surroundingNodes(Box,int) to each Box in
666 : * BoxArray. See the documentation of Box for details.
667 : */
668 : BoxArray& surroundingNodes (int dir);
669 :
670 : //! Apply Box::enclosedCells() to each Box in the BoxArray.
671 : BoxArray& enclosedCells ();
672 :
673 : //! Apply Box::enclosedCells(int) to each Box in the BoxArray.
674 : BoxArray& enclosedCells (int dir);
675 :
676 : //! Apply Box::convert(IndexType) to each Box in the BoxArray.
677 : BoxArray& convert (IndexType typ);
678 :
679 : BoxArray& convert (const IntVect& iv);
680 :
681 : //! Apply function (*fp)(Box) to each Box in the BoxArray.
682 : BoxArray& convert (Box (*fp)(const Box&));
683 :
684 : //! Apply Box::shift(int,int) to each Box in the BoxArray.
685 : BoxArray& shift (int dir, int nzones);
686 :
687 : //! Apply Box::shift(const IntVect &iv) to each Box in the BoxArray.
688 : BoxArray& shift (const IntVect &iv);
689 :
690 : //! Set element i in this BoxArray to Box ibox.
691 : void set (int i, const Box& ibox);
692 :
693 : //! Return element index of this BoxArray.
694 246523 : [[nodiscard]] Box operator[] (int index) const noexcept {
695 246523 : return m_bat(m_ref->m_abox[index]);
696 : }
697 :
698 : //! Return element index of this BoxArray.
699 : [[nodiscard]] Box operator[] (const MFIter& mfi) const noexcept;
700 :
701 : //! Return element index of this BoxArray.
702 : [[nodiscard]] Box get (int index) const noexcept { return operator[](index); }
703 :
704 : //! Return cell-centered box at element index of this BoxArray.
705 : [[nodiscard]] Box getCellCenteredBox (int index) const noexcept {
706 : return m_bat.coarsen(m_ref->m_abox[index]);
707 : }
708 :
709 : /**
710 : * \brief Return true if Box is valid and they all have the same
711 : * IndexType. Is true by default if the BoxArray is empty.
712 : */
713 : [[nodiscard]] bool ok () const;
714 :
715 : //! Return true if set of intersecting Boxes in BoxArray is null.
716 : [[nodiscard]] bool isDisjoint () const;
717 :
718 : //! Create a BoxList from this BoxArray.
719 : [[nodiscard]] BoxList boxList () const;
720 :
721 : //! True if the IntVect is within any of the Boxes in this BoxArray.
722 : [[nodiscard]] bool contains (const IntVect& v) const;
723 :
724 : /**
725 : * \brief True if the Box is contained in this BoxArray(+ng).
726 : * The Box must also have the same IndexType as those in this BoxArray.
727 : */
728 : [[nodiscard]]
729 : bool contains (const Box& b, bool assume_disjoint_ba = false,
730 : const IntVect& ng = IntVect(0)) const;
731 :
732 : /**
733 : * \brief True if all Boxes in ba are contained in this BoxArray(+ng).
734 : */
735 : [[nodiscard]]
736 : bool contains (const BoxArray& ba, bool assume_disjoint_ba = false,
737 : const IntVect& ng = IntVect(0)) const;
738 :
739 : /**
740 : * \brief True if all cells in ba are periodically contained in this
741 : * BoxArray.
742 : *
743 : * If a cell after being periodically shifted is contained in this
744 : * BoxArray, it's considered being periodically contained.
745 : */
746 : [[nodiscard]]
747 : bool contains (const BoxArray& ba, Periodicity const& period) const;
748 :
749 : //! Return smallest Box that contains all Boxes in this BoxArray.
750 : [[nodiscard]] Box minimalBox () const;
751 : [[nodiscard]] Box minimalBox (Long& npts_avg_box) const;
752 :
753 : /**
754 : * \brief True if the Box intersects with this BoxArray(+ghostcells).
755 : * The Box must have the same IndexType as those in this BoxArray.
756 : */
757 : [[nodiscard]]
758 : bool intersects (const Box& b, int ng = 0) const;
759 :
760 : [[nodiscard]]
761 : bool intersects (const Box& b, const IntVect& ng) const;
762 :
763 : //! Return intersections of Box and BoxArray
764 : [[nodiscard]]
765 : std::vector< std::pair<int,Box> > intersections (const Box& bx) const;
766 :
767 : //! Return intersections of Box and BoxArray(+ghostcells).
768 : [[nodiscard]]
769 : std::vector< std::pair<int,Box> > intersections (const Box& bx, bool first_only, int ng) const;
770 :
771 : [[nodiscard]]
772 : std::vector< std::pair<int,Box> > intersections (const Box& bx, bool first_only, const IntVect& ng) const;
773 :
774 : //! intersect Box and BoxArray, then store the result in isects
775 : void intersections (const Box& bx, std::vector< std::pair<int,Box> >& isects) const;
776 :
777 : //! intersect Box and BoxArray(+ghostcells), then store the result in isects
778 : void intersections (const Box& bx, std::vector< std::pair<int,Box> >& isects,
779 : bool first_only, int ng) const;
780 :
781 : void intersections (const Box& bx, std::vector< std::pair<int,Box> >& isects,
782 : bool first_only, const IntVect& ng) const;
783 :
784 : //! Return box - boxarray
785 : [[nodiscard]] BoxList complementIn (const Box& b) const;
786 : void complementIn (BoxList& bl, const Box& b) const;
787 :
788 : //! Clear out the internal hash table used by intersections.
789 : void clear_hash_bin () const;
790 :
791 : //! Change the BoxArray to one with no overlap and then simplify it (see the simplify function in BoxList).
792 : void removeOverlap (bool simplify=true);
793 :
794 : //! whether two BoxArrays share the same data
795 9540 : [[nodiscard]] static bool SameRefs (const BoxArray& lhs, const BoxArray& rhs) { return lhs.m_ref == rhs.m_ref; }
796 :
797 : struct RefID {
798 403505 : RefID () noexcept = default;
799 0 : explicit RefID (BARef* data_) noexcept : data(data_) {}
800 : bool operator< (const RefID& rhs) const noexcept { return std::less<>()(data,rhs.data); }
801 : bool operator== (const RefID& rhs) const noexcept { return data == rhs.data; }
802 0 : bool operator!= (const RefID& rhs) const noexcept { return data != rhs.data; }
803 : friend std::ostream& operator<< (std::ostream& os, const RefID& id);
804 : private:
805 : BARef* data{nullptr};
806 : };
807 :
808 : //! Return a unique ID of the reference
809 0 : [[nodiscard]] RefID getRefID () const noexcept { return RefID { m_ref.get() }; }
810 :
811 : //! Return index type of this BoxArray
812 6970 : [[nodiscard]] IndexType ixType () const noexcept { return m_bat.index_type(); }
813 :
814 : //! Return crse ratio of this BoxArray
815 : [[nodiscard]] IntVect crseRatio () const noexcept { return m_bat.coarsen_ratio(); }
816 :
817 : static void Initialize ();
818 : static void Finalize ();
819 : static bool initialized;
820 :
821 : //! Make ourselves unique.
822 : void uniqify ();
823 :
824 : [[nodiscard]] BoxList const& simplified_list () const; // For regular AMR grids only!!!
825 : [[nodiscard]] BoxArray simplified () const; // For regular AMR grids only!!!
826 :
827 : [[nodiscard]] BATransformer const& transformer () const;
828 :
829 : friend class AmrMesh;
830 : friend class FabArrayBase;
831 :
832 : private:
833 : //! Update BoxArray index type according the box type, and then convert boxes to cell-centered.
834 : void type_update ();
835 :
836 : [[nodiscard]] BARef::HashType& getHashMap () const;
837 :
838 : [[nodiscard]] IntVect getDoiLo () const noexcept;
839 : [[nodiscard]] IntVect getDoiHi () const noexcept;
840 :
841 : BATransformer m_bat;
842 : //! The data -- a reference-counted pointer to a Ref.
843 : std::shared_ptr<BARef> m_ref;
844 : mutable std::shared_ptr<BoxList> m_simplified_list;
845 : };
846 :
847 : //! Write a BoxArray to an ostream in ASCII format.
848 : std::ostream& operator<< (std::ostream& os, const BoxArray& ba);
849 :
850 : std::ostream& operator<< (std::ostream& os, const BoxArray::RefID& id);
851 :
852 : }
853 :
854 : #endif /*BL_BOXARRAY_H*/
|