Line data Source code
1 :
2 : #ifndef AMREX_BOX_H_
3 : #define AMREX_BOX_H_
4 : #include <AMReX_Config.H>
5 :
6 : #include <AMReX_Algorithm.H>
7 : #include <AMReX_ArrayLim.H>
8 : #include <AMReX_ccse-mpi.H>
9 : #include <AMReX_IntVect.H>
10 : #include <AMReX_IndexType.H>
11 : #include <AMReX_Orientation.H>
12 : #include <AMReX_SPACE.H>
13 : #include <AMReX_Array.H>
14 : #include <AMReX_Array4.H>
15 : #include <AMReX_Vector.H>
16 : #include <AMReX_GpuQualifiers.H>
17 : #include <AMReX_GpuControl.H>
18 : #include <AMReX_Math.H>
19 :
20 : #include <iosfwd>
21 :
22 : namespace amrex
23 : {
24 : class BoxCommHelper;
25 :
26 : /**
27 : * \brief A Rectangular Domain on an Integer Lattice
28 : *
29 : * A Box is an abstraction for defining discrete regions of
30 : * SPACEDIM indexing space. Boxes have an IndexType, which defines
31 : * IndexType::CELL or IndexType::NODE based points for each direction
32 : * and a low and high INTVECT which defines the lower and upper corners
33 : * of the Box. Boxes can exist in positive and negative indexing space.
34 : *
35 : * Box is a dimension dependent class, so SPACEDIM must be
36 : * defined as either 1, 2, or 3 when compiling.
37 : */
38 : class Box
39 : {
40 : friend MPI_Datatype ParallelDescriptor::Mpi_typemap<Box>::type();
41 : friend class BoxCommHelper;
42 :
43 : public:
44 : /*
45 : * \brief The default constructor. For safety, the constructed Box is
46 : * invalid and may be tested for validity with ok().
47 : * DO NOT CHANGE THIS BEHAVIOR!
48 : */
49 : AMREX_GPU_HOST_DEVICE
50 : constexpr Box () noexcept
51 : : smallend(1),
52 : bigend(0)
53 : {}
54 :
55 : //! Construct cell-centered type Box.
56 : AMREX_GPU_HOST_DEVICE
57 0 : constexpr Box (const IntVect& small, const IntVect& big) noexcept
58 0 : : smallend(small),
59 0 : bigend(big)
60 0 : {}
61 :
62 : //! Construct box with specified lengths.
63 : AMREX_GPU_HOST_DEVICE
64 : Box (const IntVect& small, const int* vec_len) noexcept
65 : : smallend(small),
66 : bigend(AMREX_D_DECL(small[0]+vec_len[0]-1,
67 : small[1]+vec_len[1]-1,
68 : small[2]+vec_len[2]-1))
69 : {}
70 :
71 : /**
72 : * \brief Construct Box with given type. small and big are expected
73 : * to be consistent with given type.
74 : */
75 : AMREX_GPU_HOST_DEVICE
76 1 : Box (const IntVect& small, const IntVect& big, const IntVect& typ) noexcept
77 1 : : smallend(small),
78 : bigend(big),
79 1 : btype(typ)
80 : {
81 : BL_ASSERT(typ.allGE(0) && typ.allLE(1));
82 1 : }
83 :
84 : //! Construct dimension specific Boxes.
85 : AMREX_GPU_HOST_DEVICE
86 0 : Box (const IntVect& small, const IntVect& big, IndexType t) noexcept
87 0 : : smallend(small),
88 : bigend(big),
89 0 : btype(t)
90 0 : {}
91 :
92 : template <typename T>
93 : AMREX_GPU_HOST_DEVICE
94 : explicit Box (Array4<T> const& a) noexcept
95 : : smallend(AMREX_D_DECL(a.begin.x,a.begin.y,a.begin.z)),
96 : bigend (AMREX_D_DECL(a.end.x-1,a.end.y-1,a.end.z-1))
97 : {}
98 :
99 : // dtor, copy-ctor, copy-op=, move-ctor, and move-op= are compiler generated.
100 :
101 : //! Get the smallend of the box.
102 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
103 0 : const IntVect& smallEnd () const& noexcept { return smallend; }
104 :
105 : //! Get the smallend of the box.
106 : [[nodiscard]] const IntVect& smallEnd () && = delete;
107 :
108 : //! Returns the coordinate of the low end in the given direction.
109 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
110 : int smallEnd (int dir) const& noexcept { return smallend[dir]; }
111 :
112 : //! Get the bigend.
113 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
114 0 : const IntVect& bigEnd () const& noexcept { return bigend; }
115 :
116 : //! Get the bigend.
117 : [[nodiscard]] const IntVect& bigEnd () && = delete;
118 :
119 : //! Returns the coordinate of the high end in the given direction.
120 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
121 : int bigEnd (int dir) const noexcept { return bigend[dir]; }
122 :
123 : //! Returns the indexing type.
124 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
125 : IndexType ixType () const noexcept { return btype; }
126 :
127 : //! Returns the indexing type.
128 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
129 : IntVect type () const noexcept { return btype.ixType(); }
130 :
131 : //! Returns the indexing type in the specified direction.
132 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
133 : IndexType::CellIndex type (int dir) const noexcept { return btype.ixType(dir); }
134 :
135 : //! Return the length of the Box.
136 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
137 : IntVect size () const noexcept
138 : {
139 : return IntVect(AMREX_D_DECL(bigend[0]-smallend[0] + 1,
140 : bigend[1]-smallend[1] + 1,
141 : bigend[2]-smallend[2] + 1));
142 : }
143 :
144 : //! Return the length of the Box.
145 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
146 : IntVect length () const noexcept
147 : {
148 : return IntVect(AMREX_D_DECL(bigend[0]-smallend[0] + 1,
149 : bigend[1]-smallend[1] + 1,
150 : bigend[2]-smallend[2] + 1));
151 : }
152 :
153 : //! Return the length of the Box in given direction.
154 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
155 1115620000 : int length (int dir) const noexcept { return bigend[dir] - smallend[dir] + 1; }
156 :
157 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
158 : GpuArray<int,3> length3d () const noexcept {
159 : #if (AMREX_SPACEDIM == 1)
160 : return {{bigend[0]-smallend[0]+1, 1, 1}};
161 : #elif (AMREX_SPACEDIM == 2)
162 : return {{bigend[0]-smallend[0]+1, bigend[1]-smallend[1]+1, 1}};
163 : #elif (AMREX_SPACEDIM == 3)
164 : return {{bigend[0]-smallend[0]+1, bigend[1]-smallend[1]+1, bigend[2]-smallend[2]+1}};
165 : #endif
166 : }
167 :
168 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
169 : GpuArray<int,3> loVect3d () const noexcept {
170 : #if (AMREX_SPACEDIM == 1)
171 : return {{smallend[0], 0, 0}};
172 : #elif (AMREX_SPACEDIM == 2)
173 : return {{smallend[0], smallend[1], 0}};
174 : #elif (AMREX_SPACEDIM == 3)
175 : return {{smallend[0], smallend[1], smallend[2]}};
176 : #endif
177 : }
178 :
179 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
180 : GpuArray<int,3> hiVect3d () const noexcept {
181 : #if (AMREX_SPACEDIM == 1)
182 : return {{bigend[0], 0, 0}};
183 : #elif (AMREX_SPACEDIM == 2)
184 : return {{bigend[0], bigend[1], 0}};
185 : #elif (AMREX_SPACEDIM == 3)
186 : return {{bigend[0], bigend[1], bigend[2]}};
187 : #endif
188 : }
189 :
190 : //! Returns a constant pointer the array of low end coordinates. Useful for calls to FORTRAN.
191 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
192 504824000 : const int* loVect () const& noexcept { return smallend.getVect(); }
193 : AMREX_GPU_HOST_DEVICE
194 : const int* loVect () && = delete;
195 : //! Returns a constant pointer the array of high end coordinates. Useful for calls to FORTRAN.
196 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
197 474399000 : const int* hiVect () const& noexcept { return bigend.getVect(); }
198 : AMREX_GPU_HOST_DEVICE
199 : const int* hiVect () && = delete;
200 :
201 : //! Returns the coordinate normal to given face.
202 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
203 : int operator[] (Orientation face) const noexcept {
204 : const int dir = face.coordDir();
205 : return face.isLow() ? smallend[dir] : bigend[dir];
206 : }
207 :
208 : //! Checks if it is an empty box.
209 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
210 : bool isEmpty () const noexcept { return !ok(); }
211 :
212 : //! Checks if it is a proper Box (including a valid type).
213 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
214 225123000 : bool ok () const noexcept { return bigend.allGE(smallend) && btype.ok(); }
215 :
216 : //! Returns true if argument is contained within Box.
217 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
218 : bool contains (const IntVect& p) const noexcept { return p.allGE(smallend) && p.allLE(bigend); }
219 :
220 : //! Returns true if argument is contained within Box.
221 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
222 : bool contains (const Dim3& p) const noexcept {
223 : return AMREX_D_TERM(p.x >= smallend[0] && p.x <= bigend[0],
224 : && p.y >= smallend[1] && p.y <= bigend[1],
225 : && p.z >= smallend[2] && p.z <= bigend[2]);
226 : }
227 :
228 : //! Returns true if argument is contained within Box.
229 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
230 : #if (AMREX_SPACEDIM == 1)
231 : bool contains (int i, int, int) const noexcept {
232 : #elif (AMREX_SPACEDIM == 2)
233 : bool contains (int i, int j, int) const noexcept {
234 : #else
235 : bool contains (int i, int j, int k) const noexcept {
236 : #endif
237 : return AMREX_D_TERM(i >= smallend[0] && i <= bigend[0],
238 : && j >= smallend[1] && j <= bigend[1],
239 : && k >= smallend[2] && k <= bigend[2]);
240 : }
241 :
242 : /** \brief Returns true if argument is contained within Box.
243 : * It is an error if the Boxes have different types.
244 : */
245 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
246 216 : bool contains (const Box& b) const noexcept
247 : {
248 : BL_ASSERT(sameType(b));
249 648 : return b.smallend.allGE(smallend) && b.bigend.allLE(bigend);
250 : }
251 :
252 : //! Returns true if argument is strictly contained within Box.
253 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
254 : bool strictly_contains (const IntVect& p) const noexcept { return p.allGT(smallend) && p.allLT(bigend); }
255 :
256 : /**
257 : * \brief Returns true if argument is strictly contained within Box.
258 : * It is an error if the Boxes have different types.
259 : */
260 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
261 : bool strictly_contains (const Box& b) const noexcept
262 : {
263 : BL_ASSERT(sameType(b));
264 : return b.smallend.allGT(smallend) && b.bigend.allLT(bigend);
265 : }
266 :
267 : //! Returns true if argument is strictly contained within Box.
268 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
269 : #if (AMREX_SPACEDIM == 1)
270 : bool strictly_contains (int i, int, int) const noexcept {
271 : #elif (AMREX_SPACEDIM == 2)
272 : bool strictly_contains (int i, int j, int) const noexcept {
273 : #else
274 : bool strictly_contains (int i, int j, int k) const noexcept {
275 : #endif
276 : return AMREX_D_TERM(i > smallend[0] && i < bigend[0],
277 : && j > smallend[1] && j < bigend[1],
278 : && k > smallend[2] && k < bigend[2]);
279 : }
280 :
281 : /**
282 : * \brief Returns true if Boxes have non-null intersections.
283 : * It is an error if the Boxes have different types.
284 : */
285 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
286 : bool intersects (const Box& b) const noexcept { Box isect(*this); isect &= b; return isect.ok(); }
287 :
288 : /**
289 : * \brief Returns true is Boxes same size, ie translates of each other,.
290 : * It is an error if they have different types.
291 : */
292 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
293 : bool sameSize (const Box& b) const noexcept {
294 : BL_ASSERT(sameType(b));
295 : return AMREX_D_TERM(length(0) == b.length(0),
296 : && length(1) == b.length(1),
297 : && length(2) == b.length(2));
298 : }
299 :
300 : //! Returns true if Boxes have same type.
301 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
302 : bool sameType (const Box &b) const noexcept { return btype == b.btype; }
303 :
304 : //! Returns true if Boxes are identical (including type).
305 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
306 464 : bool operator== (const Box& b) const noexcept { return smallend == b.smallend && bigend == b.bigend && b.btype == btype; }
307 :
308 : //! Returns true if Boxes differ (including type).
309 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
310 0 : bool operator!= (const Box& b) const noexcept { return !operator==(b); }
311 :
312 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
313 : bool operator< (const Box& rhs) const noexcept
314 : {
315 : return btype < rhs.btype ||
316 : ((btype == rhs.btype) &&
317 : ( (smallend < rhs.smallend) ||
318 : ((smallend == rhs.smallend) && (bigend < rhs.bigend)) ));
319 : }
320 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
321 : bool operator <= (const Box& rhs) const noexcept {
322 : return !(rhs < *this);
323 : }
324 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
325 : bool operator> (const Box& rhs) const noexcept {
326 : return rhs < *this;
327 : }
328 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
329 : bool operator>= (const Box& rhs) const noexcept {
330 : return !(*this < rhs);
331 : }
332 :
333 : //! Returns true if Box is cell-centered in all indexing directions.
334 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
335 : bool cellCentered () const noexcept { return !btype.any(); }
336 :
337 : /**
338 : * \brief Returns the number of points contained in the Box.
339 : */
340 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
341 74664200 : Long numPts () const noexcept {
342 74664200 : return ok() ? AMREX_D_TERM( static_cast<Long>(length(0)),
343 : *static_cast<Long>(length(1)),
344 : *static_cast<Long>(length(2)))
345 74664200 : : Long(0);
346 : }
347 :
348 : /**
349 : * \brief Returns the number of points contained in the Box.
350 : * This is intended for use only in diagnostic messages.
351 : */
352 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
353 : double d_numPts () const noexcept {
354 : return ok() ? AMREX_D_TERM( double(length(0)),
355 : *double(length(1)),
356 : *double(length(2)))
357 : : 0.0;
358 : }
359 :
360 : /**
361 : * \brief Return the volume, in indexing space, of region enclosed by
362 : * this Box. This is identical to numPts() for CELL centered
363 : * Box; otherwise, numPts() > volume().
364 : */
365 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
366 : Long volume () const noexcept {
367 : return ok() ? AMREX_D_TERM( static_cast<Long>(length(0)-btype[0]),
368 : *static_cast<Long>(length(1)-btype[1]),
369 : *static_cast<Long>(length(2)-btype[2]))
370 : : Long(0);
371 : }
372 :
373 : /**
374 : * \brief Returns length of longest side. dir is modified to give
375 : * direction with longest side: 0...SPACEDIM-1. Ignores type.
376 : */
377 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
378 : int longside (int& dir) const noexcept {
379 : int maxlen = length(0);
380 : dir = 0;
381 : for (int i = 1; i < AMREX_SPACEDIM; i++)
382 : {
383 : if (length(i) > maxlen)
384 : {
385 : maxlen = length(i);
386 : dir = i;
387 : }
388 : }
389 : return maxlen;
390 : }
391 :
392 : //! Returns length of longest side. Ignores type.
393 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
394 : int longside () const noexcept {
395 : int ignore = 0;
396 : return longside(ignore);
397 : }
398 :
399 : /**
400 : * \brief Returns length of shortest side. dir is modified to give
401 : * direction with shortest side: 0...SPACEDIM-1. Ignores type.
402 : */
403 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
404 : int shortside (int& dir) const noexcept {
405 : int minlen = length(0);
406 : dir = 0;
407 : for (int i = 1; i < AMREX_SPACEDIM; i++)
408 : {
409 : if (length(i) < minlen)
410 : {
411 : minlen = length(i);
412 : dir = i;
413 : }
414 : }
415 : return minlen;
416 : }
417 :
418 : //! Returns length of shortest side. Ignores type.
419 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
420 : int shortside () const noexcept {
421 : int ignore = 0;
422 : return shortside(ignore);
423 : }
424 :
425 : /**
426 : * \brief Returns offset of point from smallend; i.e.
427 : * index(smallend) -> 0, bigend would return numPts()-1.
428 : * Is used in accessing FArrayBox.
429 : */
430 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
431 : Long index (const IntVect& v) const noexcept;
432 :
433 : //! Given the offset, compute IntVect
434 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
435 : IntVect atOffset (Long offset) const noexcept;
436 :
437 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
438 : GpuArray<int,3> atOffset3d (Long offset) const noexcept;
439 :
440 : //! Redefine the small end of the Box.
441 : AMREX_GPU_HOST_DEVICE
442 : Box& setSmall (const IntVect& sm) noexcept { smallend = sm; return *this; }
443 :
444 : //! Redefine the small end of the Box.
445 : AMREX_GPU_HOST_DEVICE
446 : Box& setSmall (int dir, int sm_index) noexcept { smallend.setVal(dir,sm_index); return *this; }
447 :
448 : //! Redefine the big end of the Box.
449 : AMREX_GPU_HOST_DEVICE
450 : Box& setBig (const IntVect& bg) noexcept { bigend = bg; return *this; }
451 :
452 : //! Redefine the big end of the Box.
453 : AMREX_GPU_HOST_DEVICE
454 : Box& setBig (int dir, int bg_index) noexcept { bigend.setVal(dir,bg_index); return *this; }
455 :
456 : /**
457 : * \brief Set the entire range in a given direction, starting at
458 : * sm_index with length n_cells. NOTE: This will yield an
459 : * illegal Box if n_cells <= 0.
460 : */
461 : AMREX_GPU_HOST_DEVICE
462 : Box& setRange (int dir,
463 : int sm_index,
464 : int n_cells = 1) noexcept;
465 :
466 : //! Set indexing type
467 : AMREX_GPU_HOST_DEVICE
468 : Box& setType (const IndexType& t) noexcept { btype = t; return *this; }
469 :
470 : //! Shift this Box nzones indexing positions in coordinate direction dir.
471 : AMREX_GPU_HOST_DEVICE
472 : Box& shift (int dir, int nzones) noexcept { smallend.shift(dir,nzones); bigend.shift(dir,nzones); return *this; }
473 :
474 : //! Equivalent to b.shift(0,iv[0]).shift(1,iv[1]) ....
475 : AMREX_GPU_HOST_DEVICE
476 : Box& shift (const IntVect& iv) noexcept { smallend.shift(iv); bigend.shift(iv); return *this; }
477 :
478 : /**
479 : * \brief This member shifts the Box by "half" indices, thereby
480 : * converting the Box from type CELL to NODE and visa-versa.
481 : * b.shiftHalf(0,1) shifts b to the right by 1/2 cells.
482 : * b.shiftHalf(1,-3) shifts b in the -j direction by 3/2 cells.
483 : * NOTE: If num_halfs is EVEN the shift is num_halfs/2 full
484 : * zones and hence will not change the type.
485 : * This is: b.shifthalf(4) == b.shift(2).
486 : */
487 : AMREX_GPU_HOST_DEVICE
488 : Box& shiftHalf (int dir, int num_halfs) noexcept;
489 :
490 : //! Equivalent to b.shiftHalf(0,iv[0]).shiftHalf(1,iv[1]) ...
491 : AMREX_GPU_HOST_DEVICE
492 : Box& shiftHalf (const IntVect& iv) noexcept;
493 :
494 : /**
495 : * \brief Convert the Box from the current type into the
496 : * argument type. This may change the Box coordinates:
497 : * type CELL -> NODE : increase coordinate by one on high end
498 : * type NODE -> CELL : reduce coordinate by one on high end
499 : * other type mappings make no change.
500 : */
501 : AMREX_GPU_HOST_DEVICE
502 : Box& convert (IndexType typ) noexcept;
503 :
504 : /**
505 : * \brief Convert the Box from the current type into the
506 : * argument type. This may change the Box coordinates:
507 : * type CELL -> NODE : increase coordinate by one on high end
508 : * type NODE -> CELL : reduce coordinate by one on high end
509 : * other type mappings make no change.
510 : */
511 : AMREX_GPU_HOST_DEVICE
512 : Box& convert (const IntVect& typ) noexcept;
513 :
514 : //! Convert to NODE type in all directions.
515 : AMREX_GPU_HOST_DEVICE
516 : Box& surroundingNodes () noexcept;
517 :
518 : //! Convert to NODE type in given direction.
519 : AMREX_GPU_HOST_DEVICE
520 : Box& surroundingNodes (int dir) noexcept;
521 :
522 : AMREX_GPU_HOST_DEVICE
523 : Box& surroundingNodes (Direction d) noexcept { return surroundingNodes(static_cast<int>(d)); }
524 :
525 : //! Convert to CELL type in all directions.
526 : AMREX_GPU_HOST_DEVICE
527 : Box& enclosedCells () noexcept;
528 :
529 : //! Convert to CELL type in given direction.
530 : AMREX_GPU_HOST_DEVICE
531 : Box& enclosedCells (int dir) noexcept;
532 :
533 : AMREX_GPU_HOST_DEVICE
534 : Box& enclosedCells (Direction d) noexcept { return enclosedCells(static_cast<int>(d)); }
535 :
536 : /**
537 : * \brief Return Box that is intersection of this Box
538 : * and argument. The Boxes MUST be of same type.
539 : */
540 : AMREX_GPU_HOST_DEVICE
541 566450 : Box operator& (const Box& rhs) const noexcept { Box lhs(*this); lhs &= rhs; return lhs; }
542 :
543 : //! Intersect this Box with its argument. The Boxes MUST be of the same type.
544 : AMREX_GPU_HOST_DEVICE
545 566450 : Box& operator&= (const Box& rhs) noexcept
546 : {
547 : BL_ASSERT(sameType(rhs));
548 566450 : smallend.max(rhs.smallend);
549 566450 : bigend.min(rhs.bigend);
550 566450 : return *this;
551 : }
552 :
553 : /**
554 : * \brief Modify Box to that of the minimum Box containing both
555 : * the original Box and the argument.
556 : * Both Boxes must have identical type.
557 : */
558 : AMREX_GPU_HOST_DEVICE
559 : Box& minBox (const Box& b) noexcept {
560 : // BoxArray may call this with not ok boxes. BL_ASSERT(b.ok() && ok());
561 : BL_ASSERT(sameType(b));
562 : smallend.min(b.smallend);
563 : bigend.max(b.bigend);
564 : return *this;
565 : }
566 :
567 : //! Shift Box (relative) by given IntVect.
568 : AMREX_GPU_HOST_DEVICE
569 : Box& operator+= (const IntVect& v) noexcept { smallend += v; bigend += v; return *this; }
570 :
571 : //! Shift Box (relative) by given IntVect.
572 : AMREX_GPU_HOST_DEVICE
573 : Box operator+ (const IntVect& v) const noexcept { Box r(*this); r += v; return r; }
574 :
575 : //! Shift Box (relative) by given IntVect.
576 : AMREX_GPU_HOST_DEVICE
577 : Box& operator-= (const IntVect& v) noexcept { smallend -= v; bigend -= v; return *this; }
578 :
579 : //! Shift Box (relative) by given IntVect.
580 : AMREX_GPU_HOST_DEVICE
581 : Box operator- (const IntVect& v) const noexcept { Box r(*this); r -= v; return r; }
582 :
583 : /**
584 : * \brief Chop the Box at the chop_pnt in the dir direction
585 : * returns one Box, modifies the object Box.
586 : * The union of the two is the original Box.
587 : * The modified Box is the low end, the returned Box
588 : * is the high end. If type(dir) = CELL, the Boxes are disjoint
589 : * with the chop_pnt included in the high end (new Box).
590 : * It is an ERROR if chop_pnt is the low end of the orig Box.
591 : * If type(dir) = NODE, the chop_pnt is included in both Boxes
592 : * but is the only point in common. It is also an error if the
593 : * chop_pnt is an end node of the Box.
594 : */
595 : AMREX_GPU_HOST_DEVICE
596 : Box chop (int dir, int chop_pnt) noexcept;
597 :
598 : /*
599 : * \brief Grow Box in all directions by given amount.
600 : * NOTE: n_cell negative shrinks the Box by that number of cells.
601 : */
602 : AMREX_GPU_HOST_DEVICE
603 1087430 : Box& grow (int i) noexcept { smallend.diagShift(-i); bigend.diagShift(i); return *this; }
604 :
605 : //! Grow Box in each direction by specified amount.
606 : AMREX_GPU_HOST_DEVICE
607 846 : Box& grow (const IntVect& v) noexcept { smallend -= v; bigend += v; return *this;}
608 :
609 : /**
610 : * \brief Grow the Box on the low and high end by n_cell cells
611 : * in direction idir.
612 : */
613 : AMREX_GPU_HOST_DEVICE
614 0 : Box& grow (int idir, int n_cell) noexcept { smallend.shift(idir, -n_cell); bigend.shift(idir, n_cell); return *this; }
615 :
616 : AMREX_GPU_HOST_DEVICE
617 : Box& grow (Direction d, int n_cell) noexcept { return grow(static_cast<int>(d), n_cell); }
618 :
619 : /**
620 : * \brief Grow the Box on the low end by n_cell cells in direction idir.
621 : * NOTE: n_cell negative shrinks the Box by that number of cells.
622 : */
623 : AMREX_GPU_HOST_DEVICE
624 : Box& growLo (int idir, int n_cell = 1) noexcept { smallend.shift(idir, -n_cell); return *this; }
625 :
626 : AMREX_GPU_HOST_DEVICE
627 : Box& growLo (Direction d, int n_cell = 1) noexcept { return growLo(static_cast<int>(d), n_cell); }
628 :
629 : /**
630 : * \brief Grow the Box on the high end by n_cell cells in
631 : * direction idir. NOTE: n_cell negative shrinks the Box by that
632 : * number of cells.
633 : */
634 : AMREX_GPU_HOST_DEVICE
635 : Box& growHi (int idir, int n_cell = 1) noexcept { bigend.shift(idir,n_cell); return *this; }
636 :
637 : AMREX_GPU_HOST_DEVICE
638 : Box& growHi (Direction d, int n_cell = 1) noexcept { return growHi(static_cast<int>(d), n_cell); }
639 :
640 : //! Grow in the direction of the given face.
641 : AMREX_GPU_HOST_DEVICE
642 : Box& grow (Orientation face, int n_cell = 1) noexcept {
643 : int idir = face.coordDir();
644 : if (face.isLow()) {
645 : smallend.shift(idir, -n_cell);
646 : } else {
647 : bigend.shift(idir,n_cell);
648 : }
649 : return *this;
650 : }
651 :
652 : /**
653 : * \brief Refine Box by given (positive) refinement ratio.
654 : * NOTE: if type(dir) = CELL centered: lo <- lo*ratio and
655 : * hi <- (hi+1)*ratio - 1.
656 : * NOTE: if type(dir) = NODE centered: lo <- lo*ratio and
657 : * hi <- hi*ratio.
658 : */
659 : AMREX_GPU_HOST_DEVICE
660 : Box& refine (int ref_ratio) noexcept {
661 : return this->refine(IntVect(ref_ratio));
662 : }
663 :
664 : /*
665 : * \brief Refine Box by given (positive) refinement ratio.
666 : * NOTE: if type(dir) = CELL centered: lo <- lo*ratio and
667 : * hi <- (hi+1)*ratio - 1.
668 : * NOTE: if type(dir) = NODE centered: lo <- lo*ratio and
669 : * hi <- hi*ratio.
670 : */
671 : AMREX_GPU_HOST_DEVICE
672 : Box& refine (const IntVect& ref_ratio) noexcept;
673 :
674 : /**
675 : * \brief Coarsen Box by given (positive) refinement ratio.
676 : * NOTE: if type(dir) = CELL centered: lo <- lo/ratio and
677 : * hi <- hi/ratio.
678 : * NOTE: if type(dir) = NODE centered: lo <- lo/ratio and
679 : * hi <- hi/ratio + ((hi%ratio)==0 ? 0 : 1).
680 : * That is, refinement of coarsened Box must contain
681 : * the original Box.
682 : */
683 : AMREX_GPU_HOST_DEVICE
684 : Box& coarsen (int ref_ratio) noexcept {
685 : return this->coarsen(IntVect(ref_ratio));
686 : }
687 :
688 : /**
689 : * \brief Coarsen Box by given (positive) refinement ratio.
690 : * NOTE: if type(dir) = CELL centered: lo <- lo/ratio and
691 : * hi <- hi/ratio.
692 : * NOTE: if type(dir) = NODE centered: lo <- lo/ratio and
693 : * hi <- hi/ratio + ((hi%ratio)==0 ? 0 : 1).
694 : * That is, refinement of coarsened Box must contain
695 : * the original Box.
696 : */
697 : AMREX_GPU_HOST_DEVICE
698 : Box& coarsen (const IntVect& ref_ratio) noexcept;
699 :
700 : /**
701 : * \brief Step through the rectangle. It is a runtime error to give
702 : * a point not inside rectangle. Iteration may not be efficient.
703 : */
704 : AMREX_GPU_HOST_DEVICE
705 : void next (IntVect &) const noexcept;
706 :
707 : /**
708 : * \brief This static member function returns a constant reference to
709 : * an object of type Box representing the unit box in
710 : * AMREX_SPACEDIM-dimensional space.
711 : */
712 : AMREX_GPU_HOST_DEVICE
713 : static Box TheUnitBox () noexcept {
714 : return Box(IntVect::TheZeroVector(),IntVect::TheZeroVector());
715 : }
716 :
717 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
718 : bool isSquare() const noexcept;
719 :
720 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
721 : bool coarsenable(const IntVect& refrat, const IntVect& min_width) const noexcept
722 : {
723 : if (!size().allGE(refrat*min_width)) {
724 : return false;
725 : } else {
726 : Box testBox = *this;
727 : testBox.coarsen(refrat);
728 : testBox.refine (refrat);
729 : return (*this == testBox);
730 : }
731 : }
732 :
733 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
734 : bool coarsenable(int refrat, int min_width=1) const noexcept {
735 : return coarsenable(IntVect(refrat), IntVect(min_width));
736 : }
737 :
738 : [[nodiscard]] AMREX_GPU_HOST_DEVICE
739 : bool coarsenable(const IntVect& refrat, int min_width=1) const noexcept
740 : {
741 : return coarsenable(refrat, IntVect(min_width));
742 : }
743 :
744 : AMREX_GPU_HOST_DEVICE
745 : void normalize () noexcept
746 : {
747 : for (int idim=0; idim < AMREX_SPACEDIM; ++idim) {
748 : if (this->length(idim) == 0) {
749 : this->growHi(idim,1);
750 : }
751 : }
752 : }
753 :
754 : AMREX_GPU_HOST_DEVICE
755 : Box& makeSlab (int direction, int slab_index) noexcept
756 : {
757 : smallend[direction] = slab_index;
758 : bigend[direction] = slab_index;
759 : return *this;
760 : }
761 :
762 : AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 lbound (Box const& box) noexcept;
763 : AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 ubound (Box const& box) noexcept;
764 : AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 begin (Box const& box) noexcept;
765 : AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 end (Box const& box) noexcept;
766 : AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 length (Box const& box) noexcept;
767 : AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 max_lbound (Box const&, Box const&) noexcept;
768 : AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 max_lbound (Box const&, Dim3 const&) noexcept;
769 : AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 min_ubound (Box const&, Box const&) noexcept;
770 : AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 min_ubound (Box const&, Dim3 const&) noexcept;
771 : AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Box minBox (Box const&, Box const&, IndexType) noexcept;
772 :
773 : private:
774 : IntVect smallend;
775 : IntVect bigend;
776 : IndexType btype;
777 : };
778 :
779 : AMREX_GPU_HOST_DEVICE
780 : AMREX_FORCE_INLINE
781 : Box&
782 : Box::refine (const IntVect& ref_ratio) noexcept
783 : {
784 4770 : if (ref_ratio != 1) {
785 4770 : IntVect shft(1);
786 9540 : shft -= btype.ixType();
787 4770 : smallend *= ref_ratio;
788 4770 : bigend += shft;
789 4770 : bigend *= ref_ratio;
790 4770 : bigend -= shft;
791 : }
792 4770 : return *this;
793 : }
794 :
795 : AMREX_GPU_HOST_DEVICE
796 : AMREX_FORCE_INLINE
797 : Box&
798 : Box::coarsen (const IntVect& ref_ratio) noexcept
799 : {
800 497138 : if (ref_ratio != 1)
801 : {
802 497138 : smallend.coarsen(ref_ratio);
803 :
804 994276 : if (btype.any())
805 : {
806 4770 : IntVect off(0);
807 19080 : for (int dir = 0; dir < AMREX_SPACEDIM; dir++)
808 : {
809 28620 : if (btype[dir]) {
810 42930 : if (bigend[dir]%ref_ratio[dir]) {
811 : off.setVal(dir,1);
812 : }
813 : }
814 : }
815 4770 : bigend.coarsen(ref_ratio);
816 4770 : bigend += off;
817 : }
818 : else
819 : {
820 492368 : bigend.coarsen(ref_ratio);
821 : }
822 : }
823 :
824 497138 : return *this;
825 : }
826 :
827 : AMREX_GPU_HOST_DEVICE
828 : AMREX_FORCE_INLINE
829 : Box&
830 : Box::convert (const IntVect& typ) noexcept
831 : {
832 : BL_ASSERT(typ.allGE(0) && typ.allLE(1));
833 738490 : IntVect shft(typ - btype.ixType());
834 369245 : bigend += shft;
835 369245 : btype = IndexType(typ);
836 369245 : return *this;
837 : }
838 :
839 : AMREX_GPU_HOST_DEVICE
840 : AMREX_FORCE_INLINE
841 : Box&
842 : Box::convert (IndexType t) noexcept
843 : {
844 2611492 : for (int dir = 0; dir < AMREX_SPACEDIM; dir++)
845 : {
846 1958624 : const auto typ = t[dir];
847 1958624 : const auto bitval = btype[dir];
848 1958624 : const int off = typ - bitval;
849 1958624 : bigend.shift(dir,off);
850 1958624 : btype.setType(dir, (IndexType::CellIndex) typ);
851 : }
852 652874 : return *this;
853 : }
854 :
855 : AMREX_GPU_HOST_DEVICE
856 : AMREX_FORCE_INLINE
857 : Box&
858 : Box::surroundingNodes (int dir) noexcept
859 : {
860 : if (!(btype[dir]))
861 : {
862 : bigend.shift(dir,1);
863 : //
864 : // Set dir'th bit to 1 = IndexType::NODE.
865 : //
866 : btype.set(dir);
867 : }
868 : return *this;
869 : }
870 :
871 : AMREX_GPU_HOST_DEVICE
872 : AMREX_FORCE_INLINE
873 : Box&
874 : Box::surroundingNodes () noexcept
875 : {
876 : for (int i = 0; i < AMREX_SPACEDIM; ++i) {
877 : if ((btype[i] == 0)) {
878 : bigend.shift(i,1);
879 : }
880 : }
881 : btype.setall();
882 : return *this;
883 : }
884 :
885 : AMREX_GPU_HOST_DEVICE
886 : AMREX_FORCE_INLINE
887 : Box&
888 : Box::enclosedCells (int dir) noexcept
889 : {
890 : if (btype[dir])
891 : {
892 : bigend.shift(dir,-1);
893 : //
894 : // Set dir'th bit to 0 = IndexType::CELL.
895 : //
896 : btype.unset(dir);
897 : }
898 : return *this;
899 : }
900 :
901 : AMREX_GPU_HOST_DEVICE
902 : AMREX_FORCE_INLINE
903 : Box&
904 : Box::enclosedCells () noexcept
905 : {
906 0 : for (int i = 0 ; i < AMREX_SPACEDIM; ++i) {
907 0 : if (btype[i]) {
908 0 : bigend.shift(i,-1);
909 : }
910 : }
911 0 : btype.clear();
912 0 : return *this;
913 : }
914 :
915 : AMREX_GPU_HOST_DEVICE
916 : AMREX_FORCE_INLINE
917 : Long
918 : Box::index (const IntVect& v) const noexcept
919 : {
920 73925750 : Long result = v[0]-smallend[0];
921 : #if AMREX_SPACEDIM==2
922 : result += length(0)*Long(v[1]-smallend[1]);
923 : #elif AMREX_SPACEDIM==3
924 221776750 : result += length(0)*((v[1]-smallend[1])
925 147851500 : +Long(v[2]-smallend[2])*length(1));
926 : #endif
927 73925750 : return result;
928 : }
929 :
930 : AMREX_GPU_HOST_DEVICE
931 : AMREX_FORCE_INLINE
932 : IntVect
933 : Box::atOffset (Long offset) const noexcept
934 : {
935 : #if (AMREX_SPACEDIM == 1)
936 : return IntVect{static_cast<int>(offset+smallend[0])};
937 : #elif (AMREX_SPACEDIM == 2)
938 : int xlen = bigend[0]-smallend[0]+1;
939 : Long j = offset / xlen;
940 : Long i = offset - j*xlen;
941 : return IntVect{static_cast<int>(i+smallend[0]),
942 : static_cast<int>(j+smallend[1])};
943 : #elif (AMREX_SPACEDIM == 3)
944 : int xlen = bigend[0]-smallend[0]+1;
945 : int ylen = bigend[1]-smallend[1]+1;
946 : Long k = offset / (xlen*ylen);
947 : Long j = (offset - k*(xlen*ylen)) / xlen;
948 : Long i = (offset - k*(xlen*ylen)) - j*xlen;
949 : return IntVect{static_cast<int>(i+smallend[0]),
950 : static_cast<int>(j+smallend[1]),
951 : static_cast<int>(k+smallend[2])};
952 : #endif
953 : }
954 :
955 : AMREX_GPU_HOST_DEVICE
956 : AMREX_FORCE_INLINE
957 : GpuArray<int,3>
958 : Box::atOffset3d (Long offset) const noexcept
959 : {
960 : #if (AMREX_SPACEDIM == 1)
961 : return {{static_cast<int>(offset+smallend[0]),
962 : static_cast<int>(0),
963 : static_cast<int>(0)}};
964 : #elif (AMREX_SPACEDIM == 2)
965 : int xlen = bigend[0]-smallend[0]+1;
966 : Long j = offset / xlen;
967 : Long i = offset - j*xlen;
968 : return {{static_cast<int>(i+smallend[0]),
969 : static_cast<int>(j+smallend[1]),
970 : static_cast<int>(0)}};
971 : #elif (AMREX_SPACEDIM == 3)
972 : int xlen = bigend[0]-smallend[0]+1;
973 : int ylen = bigend[1]-smallend[1]+1;
974 : Long k = offset / (xlen*ylen);
975 : Long j = (offset - k*(xlen*ylen)) / xlen;
976 : Long i = (offset - k*(xlen*ylen)) - j*xlen;
977 : return {{static_cast<int>(i+smallend[0]),
978 : static_cast<int>(j+smallend[1]),
979 : static_cast<int>(k+smallend[2])}};
980 : #endif
981 : }
982 :
983 : AMREX_GPU_HOST_DEVICE
984 : AMREX_FORCE_INLINE
985 : Box&
986 : Box::setRange (int dir,
987 : int sm_index,
988 : int n_cells) noexcept
989 : {
990 : smallend.setVal(dir,sm_index);
991 : bigend.setVal(dir,sm_index+n_cells-1);
992 : return *this;
993 : }
994 :
995 : AMREX_GPU_HOST_DEVICE
996 : AMREX_FORCE_INLINE
997 : void
998 : Box::next (IntVect& p) const noexcept // NOLINT(readability-convert-member-functions-to-static)
999 : {
1000 : BL_ASSERT(contains(p));
1001 :
1002 : ++p[0];
1003 :
1004 : #if (AMREX_SPACEDIM >= 2)
1005 : if (p[0] > bigend[0])
1006 : {
1007 : p[0] = smallend[0];
1008 : ++p[1];
1009 : #if (AMREX_SPACEDIM == 3)
1010 : if (p[1] > bigend[1])
1011 : {
1012 : p[1] = smallend[1];
1013 : ++p[2];
1014 : }
1015 : #endif
1016 : }
1017 : #endif
1018 : }
1019 :
1020 : AMREX_GPU_HOST_DEVICE
1021 : AMREX_FORCE_INLINE
1022 : bool
1023 : Box::isSquare () const noexcept // NOLINT(readability-convert-member-functions-to-static)
1024 : {
1025 : #if AMREX_SPACEDIM==1
1026 : return false; // can't build a square in 1-D
1027 : #elif AMREX_SPACEDIM==2
1028 : const IntVect& sz = this->size();
1029 : return (sz[0] == sz[1]);
1030 : #elif AMREX_SPACEDIM==3
1031 : const IntVect& sz = this->size();
1032 : return (sz[0] == sz[1] && (sz[1] == sz[2]));
1033 : #endif
1034 : }
1035 :
1036 : //
1037 : // Modified Box is low end, returned Box is high end.
1038 : // If CELL: chop_pnt included in high end.
1039 : // If NODE: chop_pnt included in both Boxes.
1040 : //
1041 :
1042 : AMREX_GPU_HOST_DEVICE
1043 : inline
1044 : Box
1045 : Box::chop (int dir, int chop_pnt) noexcept
1046 : {
1047 : //
1048 : // Define new high end Box including chop_pnt.
1049 : //
1050 : IntVect sm(smallend);
1051 : IntVect bg(bigend);
1052 : sm.setVal(dir,chop_pnt);
1053 : if (btype[dir])
1054 : {
1055 : //
1056 : // NODE centered Box.
1057 : //
1058 : BL_ASSERT(chop_pnt > smallend[dir] && chop_pnt < bigend[dir]);
1059 : //
1060 : // Shrink original Box to just contain chop_pnt.
1061 : //
1062 : bigend.setVal(dir,chop_pnt);
1063 : }
1064 : else
1065 : {
1066 : //
1067 : // CELL centered Box.
1068 : //
1069 : BL_ASSERT(chop_pnt > smallend[dir] && chop_pnt <= bigend[dir]);
1070 : //
1071 : // Shrink original Box to one below chop_pnt.
1072 : //
1073 : bigend.setVal(dir,chop_pnt-1);
1074 : }
1075 : return Box(sm,bg,btype);
1076 : }
1077 :
1078 : AMREX_GPU_HOST_DEVICE
1079 : inline
1080 : Box&
1081 : Box::shiftHalf (int dir, int num_halfs) noexcept
1082 : {
1083 : const int nbit = (num_halfs<0 ? -num_halfs : num_halfs)%2;
1084 : int nshift = num_halfs/2;
1085 : //
1086 : // Toggle btyp bit if num_halfs is odd.
1087 : //
1088 : const unsigned int bit_dir = btype[dir];
1089 : if (nbit) {
1090 : btype.flip(dir);
1091 : }
1092 : if (num_halfs < 0) {
1093 : nshift -= (bit_dir ? nbit : 0);
1094 : } else {
1095 : nshift += (bit_dir ? 0 : nbit);
1096 : }
1097 : smallend.shift(dir,nshift);
1098 : bigend.shift(dir,nshift);
1099 : return *this;
1100 : }
1101 :
1102 : AMREX_GPU_HOST_DEVICE
1103 : inline
1104 : Box&
1105 : Box::shiftHalf (const IntVect& iv) noexcept
1106 : {
1107 : for (int i = 0; i < AMREX_SPACEDIM; i++) {
1108 : shiftHalf(i,iv[i]);
1109 : }
1110 : return *this;
1111 : }
1112 :
1113 : class BoxCommHelper
1114 : {
1115 : public:
1116 :
1117 : explicit BoxCommHelper (const Box& bx, int* p_ = nullptr);
1118 :
1119 : [[nodiscard]] int* data () const& noexcept { return p; }
1120 : int* data () && = delete;
1121 :
1122 : [[nodiscard]] Box make_box () const noexcept {
1123 : return Box(IntVect(p), IntVect(p+AMREX_SPACEDIM), IntVect(p+2*AMREX_SPACEDIM));
1124 : }
1125 :
1126 : [[nodiscard]] static int size () noexcept { return 3*AMREX_SPACEDIM; }
1127 :
1128 : private:
1129 : int* p;
1130 : std::vector<int> v;
1131 : };
1132 :
1133 : class BoxConverter { // NOLINT
1134 : public:
1135 : virtual Box doit (const Box& fine) const = 0; // NOLINT
1136 : virtual BoxConverter* clone () const = 0; // NOLINT
1137 74 : virtual ~BoxConverter () = default;
1138 : };
1139 :
1140 : void AllGatherBoxes (Vector<Box>& bxs, int n_extra_reserve=0);
1141 :
1142 : /**
1143 : * \brief Grow Box in all directions by given amount.
1144 :
1145 : * NOTE: n_cell negative shrinks the Box by that number of cells.
1146 : */
1147 : [[nodiscard]]
1148 : AMREX_GPU_HOST_DEVICE
1149 : AMREX_FORCE_INLINE
1150 : Box grow (const Box& b, int i) noexcept
1151 : {
1152 : Box result = b;
1153 : result.grow(i);
1154 : return result;
1155 : }
1156 :
1157 : //! Grow Box in each direction by specified amount.
1158 : [[nodiscard]]
1159 : AMREX_GPU_HOST_DEVICE
1160 : AMREX_FORCE_INLINE
1161 : Box grow (const Box& b, const IntVect& v) noexcept
1162 : {
1163 282 : Box result = b;
1164 282 : result.grow(v);
1165 282 : return result;
1166 : }
1167 :
1168 : //! Grow Box in direction idir be n_cell cells
1169 : [[nodiscard]]
1170 : AMREX_GPU_HOST_DEVICE
1171 : AMREX_FORCE_INLINE
1172 : Box grow (const Box& b, int idir, int n_cell) noexcept
1173 : {
1174 : Box result = b;
1175 : result.grow(idir, n_cell);
1176 : return result;
1177 : }
1178 :
1179 : [[nodiscard]]
1180 : AMREX_GPU_HOST_DEVICE
1181 : AMREX_FORCE_INLINE
1182 : Box grow (const Box& b, Direction d, int n_cell) noexcept
1183 : {
1184 : return grow(b, static_cast<int>(d), n_cell);
1185 : }
1186 :
1187 : [[nodiscard]]
1188 : AMREX_GPU_HOST_DEVICE
1189 : AMREX_FORCE_INLINE
1190 : Box growLo (const Box& b, int idir, int n_cell) noexcept
1191 : {
1192 : Box result = b;
1193 : result.growLo(idir, n_cell);
1194 : return result;
1195 : }
1196 :
1197 : [[nodiscard]]
1198 : AMREX_GPU_HOST_DEVICE
1199 : AMREX_FORCE_INLINE
1200 : Box growLo (const Box& b, Direction d, int n_cell) noexcept
1201 : {
1202 : return growLo(b, static_cast<int>(d), n_cell);
1203 : }
1204 :
1205 : [[nodiscard]]
1206 : AMREX_GPU_HOST_DEVICE
1207 : AMREX_FORCE_INLINE
1208 : Box growHi (const Box& b, int idir, int n_cell) noexcept
1209 : {
1210 : Box result = b;
1211 : result.growHi(idir, n_cell);
1212 : return result;
1213 : }
1214 :
1215 : [[nodiscard]]
1216 : AMREX_GPU_HOST_DEVICE
1217 : AMREX_FORCE_INLINE
1218 : Box growHi (const Box& b, Direction d, int n_cell) noexcept
1219 : {
1220 : return growHi(b, static_cast<int>(d), n_cell);
1221 : }
1222 :
1223 : /**
1224 : * \brief Coarsen Box by given (positive) refinement ratio.
1225 : * NOTE: if type(dir) = CELL centered: lo <- lo/ratio and
1226 : * hi <- hi/ratio.
1227 : * NOTE: if type(dir) = NODE centered: lo <- lo/ratio and
1228 : * hi <- hi/ratio + ((hi%ratio)==0 ? 0 : 1).
1229 : * That is, refinement of coarsened Box must contain
1230 : * the original Box.
1231 : */
1232 : [[nodiscard]]
1233 : AMREX_GPU_HOST_DEVICE
1234 : AMREX_FORCE_INLINE
1235 : Box coarsen (const Box& b, int ref_ratio) noexcept
1236 : {
1237 4770 : Box result = b;
1238 4770 : result.coarsen(IntVect(ref_ratio));
1239 4770 : return result;
1240 : }
1241 :
1242 : /**
1243 : * \brief Coarsen Box by given (positive) refinement ratio.
1244 : * NOTE: if type(dir) = CELL centered: lo <- lo/ratio and
1245 : * hi <- hi/ratio.
1246 : * NOTE: if type(dir) = NODE centered: lo <- lo/ratio and
1247 : * hi <- hi/ratio + ((hi%ratio)==0 ? 0 : 1).
1248 : * That is, refinement of coarsened Box must contain
1249 : * the original Box.
1250 : */
1251 : [[nodiscard]]
1252 : AMREX_GPU_HOST_DEVICE
1253 : AMREX_FORCE_INLINE
1254 : Box coarsen (const Box& b, const IntVect& ref_ratio) noexcept
1255 : {
1256 492368 : Box result = b;
1257 : result.coarsen(ref_ratio);
1258 492368 : return result;
1259 : }
1260 :
1261 : /**
1262 : * Refine Box by given (positive) refinement ratio.
1263 : * NOTE: if type(dir) = CELL centered: lo <- lo*ratio and
1264 : * hi <- (hi+1)*ratio - 1.
1265 : * NOTE: if type(dir) = NODE centered: lo <- lo*ratio and
1266 : * hi <- hi*ratio.
1267 : */
1268 : [[nodiscard]]
1269 : AMREX_GPU_HOST_DEVICE
1270 : AMREX_FORCE_INLINE
1271 : Box refine (const Box& b, int ref_ratio) noexcept
1272 : {
1273 4770 : Box result = b;
1274 4770 : result.refine(IntVect(ref_ratio));
1275 4770 : return result;
1276 : }
1277 :
1278 : /**
1279 : * \brief Refine Box by given (positive) refinement ratio.
1280 : * NOTE: if type(dir) = CELL centered: lo <- lo*ratio and
1281 : * hi <- (hi+1)*ratio - 1.
1282 : * NOTE: if type(dir) = NODE centered: lo <- lo*ratio and
1283 : * hi <- hi*ratio.
1284 : */
1285 : [[nodiscard]]
1286 : AMREX_GPU_HOST_DEVICE
1287 : AMREX_FORCE_INLINE
1288 : Box refine (const Box& b, const IntVect& ref_ratio) noexcept
1289 : {
1290 : Box result = b;
1291 : result.refine(ref_ratio);
1292 : return result;
1293 : }
1294 :
1295 : //! Return a Box with indices shifted by nzones in dir direction.
1296 : [[nodiscard]]
1297 : AMREX_GPU_HOST_DEVICE
1298 : AMREX_FORCE_INLINE
1299 : Box shift (const Box& b, int dir, int nzones) noexcept
1300 : {
1301 : Box result = b;
1302 : result.shift(dir, nzones);
1303 : return result;
1304 : }
1305 :
1306 : [[nodiscard]]
1307 : AMREX_GPU_HOST_DEVICE
1308 : AMREX_FORCE_INLINE
1309 : Box shift (const Box& b, const IntVect& nzones) noexcept
1310 : {
1311 : Box result = b;
1312 : result.shift(nzones);
1313 : return result;
1314 : }
1315 :
1316 : /**
1317 : * \brief Returns a Box with NODE based coordinates in direction dir
1318 : * that encloses Box b. NOTE: equivalent to b.convert(dir,NODE)
1319 : * NOTE: error if b.type(dir) == NODE.
1320 : */
1321 : [[nodiscard]]
1322 : AMREX_GPU_HOST_DEVICE
1323 : AMREX_FORCE_INLINE
1324 : Box surroundingNodes (const Box& b, int dir) noexcept
1325 : {
1326 : Box bx(b);
1327 : bx.surroundingNodes(dir);
1328 : return bx;
1329 : }
1330 :
1331 : [[nodiscard]]
1332 : AMREX_GPU_HOST_DEVICE
1333 : AMREX_FORCE_INLINE
1334 : Box surroundingNodes (const Box& b, Direction d) noexcept
1335 : {
1336 : return surroundingNodes(b, static_cast<int>(d));
1337 : }
1338 :
1339 : /**
1340 : * \brief Returns a Box with NODE based coordinates in all
1341 : * directions that encloses Box b.
1342 : */
1343 : [[nodiscard]]
1344 : AMREX_GPU_HOST_DEVICE
1345 : AMREX_FORCE_INLINE
1346 : Box surroundingNodes (const Box& b) noexcept
1347 : {
1348 : Box bx(b);
1349 : bx.surroundingNodes();
1350 : return bx;
1351 : }
1352 :
1353 : //! Returns a Box with different type
1354 : [[nodiscard]]
1355 : AMREX_GPU_HOST_DEVICE
1356 : AMREX_FORCE_INLINE
1357 : Box convert (const Box& b, const IntVect& typ) noexcept
1358 : {
1359 : Box bx(b);
1360 : bx.convert(typ);
1361 : return bx;
1362 : }
1363 :
1364 : [[nodiscard]]
1365 : AMREX_GPU_HOST_DEVICE
1366 : AMREX_FORCE_INLINE
1367 : Box convert (const Box& b, const IndexType& typ) noexcept
1368 : {
1369 652874 : Box bx(b);
1370 652874 : bx.convert(typ);
1371 652874 : return bx;
1372 : }
1373 :
1374 : /**
1375 : * \brief Returns a Box with CELL based coordinates in
1376 : * direction dir that is enclosed by b.
1377 : * NOTE: equivalent to b.convert(dir,CELL)
1378 : * NOTE: error if b.type(dir) == CELL.
1379 : */
1380 : [[nodiscard]]
1381 : AMREX_GPU_HOST_DEVICE
1382 : AMREX_FORCE_INLINE
1383 : Box enclosedCells (const Box& b, int dir) noexcept
1384 : {
1385 : Box bx(b);
1386 : bx.enclosedCells(dir);
1387 : return bx;
1388 : }
1389 :
1390 : [[nodiscard]]
1391 : AMREX_GPU_HOST_DEVICE
1392 : AMREX_FORCE_INLINE
1393 : Box enclosedCells (const Box& b, Direction d) noexcept
1394 : {
1395 : return enclosedCells(b, static_cast<int>(d));
1396 : }
1397 :
1398 : /**
1399 : * \brief Returns a Box with CELL based coordinates in all
1400 : * directions that is enclosed by b.
1401 : */
1402 : [[nodiscard]]
1403 : AMREX_GPU_HOST_DEVICE
1404 : AMREX_FORCE_INLINE
1405 : Box enclosedCells (const Box& b) noexcept
1406 : {
1407 0 : Box bx(b);
1408 : bx.enclosedCells();
1409 0 : return bx;
1410 : }
1411 :
1412 : /**
1413 : * \brief Returns the edge-centered Box (in direction dir) defining
1414 : * the low side of Box b.
1415 : */
1416 : [[nodiscard]]
1417 : AMREX_GPU_HOST_DEVICE
1418 : AMREX_FORCE_INLINE
1419 : Box bdryLo (const Box& b, int dir, int len=1) noexcept
1420 : {
1421 : IntVect low(b.smallEnd());
1422 : IntVect hi(b.bigEnd());
1423 : int sm = low[dir];
1424 : low.setVal(dir,sm-len+1);
1425 : hi.setVal(dir,sm);
1426 : //
1427 : // set dir'th bit to 1 = IndexType::NODE.
1428 : //
1429 : IndexType typ(b.ixType());
1430 : typ.set(dir);
1431 : return Box(low,hi,typ);
1432 : }
1433 :
1434 : /**
1435 : * \brief Returns the edge-centered Box (in direction dir) defining
1436 : * the high side of Box b.
1437 : */
1438 : [[nodiscard]]
1439 : AMREX_GPU_HOST_DEVICE
1440 : AMREX_FORCE_INLINE
1441 : Box bdryHi (const Box& b, int dir, int len=1) noexcept
1442 : {
1443 : IntVect low(b.smallEnd());
1444 : IntVect hi(b.bigEnd());
1445 : auto const bitval = b.type()[dir];
1446 : int bg = hi[dir] + 1 - bitval%2;
1447 : low.setVal(dir,bg);
1448 : hi.setVal(dir,bg+len-1);
1449 : //
1450 : // Set dir'th bit to 1 = IndexType::NODE.
1451 : //
1452 : IndexType typ(b.ixType());
1453 : typ.set(dir);
1454 : return Box(low,hi,typ);
1455 : }
1456 :
1457 : /**
1458 : * \brief Similar to bdryLo and bdryHi except that it operates on the
1459 : * given face of box b.
1460 : */
1461 : [[nodiscard]]
1462 : AMREX_GPU_HOST_DEVICE
1463 : AMREX_FORCE_INLINE
1464 : Box bdryNode (const Box& b, Orientation face, int len=1) noexcept
1465 : {
1466 : int dir = face.coordDir();
1467 : IntVect low(b.smallEnd());
1468 : IntVect hi(b.bigEnd());
1469 : if (face.isLow())
1470 : {
1471 : int sm = low[dir];
1472 : low.setVal(dir,sm-len+1);
1473 : hi.setVal(dir,sm);
1474 : }
1475 : else
1476 : {
1477 : int bitval = b.type()[dir];
1478 : int bg = hi[dir] + 1 - bitval%2;
1479 : low.setVal(dir,bg);
1480 : hi.setVal(dir,bg+len-1);
1481 : }
1482 : //
1483 : // Set dir'th bit to 1 = IndexType::NODE.
1484 : //
1485 : IndexType typ(b.ixType());
1486 : typ.set(dir);
1487 : return Box(low,hi,typ);
1488 : }
1489 :
1490 : /**
1491 : * \brief Returns the cell centered Box of length len adjacent
1492 : * to b on the low end along the coordinate direction dir.
1493 : * The return Box is identical to b in the other directions.
1494 : * The return Box and b have an empty intersection.
1495 : * NOTE: len >= 1
1496 : * NOTE: Box retval = b.adjCellLo(b,dir,len)
1497 : * is equivalent to the following set of operations:
1498 : * Box retval(b);
1499 : * retval.convert(dir,Box::CELL);
1500 : * retval.setrange(dir,retval.smallEnd(dir)-len,len);
1501 : */
1502 : [[nodiscard]]
1503 : AMREX_GPU_HOST_DEVICE
1504 : AMREX_FORCE_INLINE
1505 : Box adjCellLo (const Box& b, int dir, int len=1) noexcept
1506 : {
1507 : BL_ASSERT(len > 0);
1508 : IntVect low(b.smallEnd());
1509 : IntVect hi(b.bigEnd());
1510 : int sm = low[dir];
1511 : low.setVal(dir,sm - len);
1512 : hi.setVal(dir,sm - 1);
1513 : //
1514 : // Set dir'th bit to 0 = IndexType::CELL.
1515 : //
1516 : IndexType typ(b.ixType());
1517 : typ.unset(dir);
1518 : return Box(low,hi,typ);
1519 : }
1520 :
1521 : //! Similar to adjCellLo but builds an adjacent Box on the high end.
1522 : [[nodiscard]]
1523 : AMREX_GPU_HOST_DEVICE
1524 : AMREX_FORCE_INLINE
1525 : Box adjCellHi (const Box& b, int dir, int len=1) noexcept
1526 : {
1527 : BL_ASSERT(len > 0);
1528 : IntVect low(b.smallEnd());
1529 : IntVect hi(b.bigEnd());
1530 : int bitval = b.type()[dir];
1531 : int bg = hi[dir] + 1 - bitval%2;
1532 : low.setVal(dir,bg);
1533 : hi.setVal(dir,bg + len - 1);
1534 : //
1535 : // Set dir'th bit to 0 = IndexType::CELL.
1536 : //
1537 : IndexType typ(b.ixType());
1538 : typ.unset(dir);
1539 : return Box(low,hi,typ);
1540 : }
1541 :
1542 : //! Similar to adjCellLo and adjCellHi; operates on given face.
1543 : [[nodiscard]]
1544 : AMREX_GPU_HOST_DEVICE
1545 : AMREX_FORCE_INLINE
1546 : Box adjCell (const Box& b, Orientation face, int len=1) noexcept
1547 : {
1548 : BL_ASSERT(len > 0);
1549 : IntVect low(b.smallEnd());
1550 : IntVect hi(b.bigEnd());
1551 : int dir = face.coordDir();
1552 : if (face.isLow())
1553 : {
1554 : int sm = low[dir];
1555 : low.setVal(dir,sm - len);
1556 : hi.setVal(dir,sm - 1);
1557 : }
1558 : else
1559 : {
1560 : int bitval = b.type()[dir];
1561 : int bg = hi[dir] + 1 - bitval%2;
1562 : low.setVal(dir,bg);
1563 : hi.setVal(dir,bg + len - 1);
1564 : }
1565 : //
1566 : // Set dir'th bit to 0 = IndexType::CELL.
1567 : //
1568 : IndexType typ(b.ixType());
1569 : typ.unset(dir);
1570 : return Box(low,hi,typ);
1571 : }
1572 :
1573 : /**
1574 : * \brief Modify Box to that of the minimum Box containing both
1575 : * the original Box and the argument.
1576 : * Both Boxes must have identical type.
1577 : */
1578 : [[nodiscard]]
1579 : AMREX_GPU_HOST_DEVICE
1580 : AMREX_FORCE_INLINE
1581 : Box minBox (const Box& b1, const Box& b2) noexcept
1582 : {
1583 : Box result = b1;
1584 : result.minBox(b2);
1585 : return result;
1586 : }
1587 :
1588 : //! Write an ASCII representation to the ostream.
1589 : std::ostream& operator<< (std::ostream& os, const Box& bx);
1590 :
1591 : //! Read from istream.
1592 : std::istream& operator>> (std::istream& is, Box& bx);
1593 :
1594 : [[nodiscard]]
1595 : AMREX_GPU_HOST_DEVICE
1596 : AMREX_FORCE_INLINE
1597 : Dim3 lbound (Box const& box) noexcept
1598 : {
1599 : #if (AMREX_SPACEDIM == 1)
1600 : return {box.smallend[0], 0, 0};
1601 : #elif (AMREX_SPACEDIM == 2)
1602 : return {box.smallend[0], box.smallend[1], 0};
1603 : #elif (AMREX_SPACEDIM == 3)
1604 88504178 : return {box.smallend[0], box.smallend[1], box.smallend[2]};
1605 : #endif
1606 : }
1607 :
1608 : [[nodiscard]]
1609 : AMREX_GPU_HOST_DEVICE
1610 : AMREX_FORCE_INLINE
1611 : Dim3 ubound (Box const& box) noexcept
1612 : {
1613 : #if (AMREX_SPACEDIM == 1)
1614 : return {box.bigend[0], 0, 0};
1615 : #elif (AMREX_SPACEDIM == 2)
1616 : return {box.bigend[0], box.bigend[1], 0};
1617 : #elif (AMREX_SPACEDIM == 3)
1618 88359878 : return {box.bigend[0], box.bigend[1], box.bigend[2]};
1619 : #endif
1620 : }
1621 :
1622 : [[nodiscard]]
1623 : AMREX_GPU_HOST_DEVICE
1624 : AMREX_FORCE_INLINE
1625 : Dim3 begin (Box const& box) noexcept
1626 : {
1627 : #if (AMREX_SPACEDIM == 1)
1628 : return {box.smallend[0], 0, 0};
1629 : #elif (AMREX_SPACEDIM == 2)
1630 : return {box.smallend[0], box.smallend[1], 0};
1631 : #elif (AMREX_SPACEDIM == 3)
1632 18576058 : return {box.smallend[0], box.smallend[1], box.smallend[2]};
1633 : #endif
1634 : }
1635 :
1636 : [[nodiscard]]
1637 : AMREX_GPU_HOST_DEVICE
1638 : AMREX_FORCE_INLINE
1639 : Dim3 end (Box const& box) noexcept
1640 : {
1641 : #if (AMREX_SPACEDIM == 1)
1642 : return {box.bigend[0]+1, 1, 1};
1643 : #elif (AMREX_SPACEDIM == 2)
1644 : return {box.bigend[0]+1, box.bigend[1]+1, 1};
1645 : #elif (AMREX_SPACEDIM == 3)
1646 18576058 : return {box.bigend[0]+1, box.bigend[1]+1, box.bigend[2]+1};
1647 : #endif
1648 : }
1649 :
1650 : [[nodiscard]]
1651 : AMREX_GPU_HOST_DEVICE
1652 : AMREX_FORCE_INLINE
1653 : Dim3 length (Box const& box) noexcept
1654 : {
1655 : #if (AMREX_SPACEDIM == 1)
1656 : return {box.bigend[0]-box.smallend[0]+1, 1, 1};
1657 : #elif (AMREX_SPACEDIM == 2)
1658 : return {box.bigend[0]-box.smallend[0]+1,
1659 : box.bigend[1]-box.smallend[1]+1, 1};
1660 : #elif (AMREX_SPACEDIM == 3)
1661 : return {box.bigend[0]-box.smallend[0]+1,
1662 : box.bigend[1]-box.smallend[1]+1,
1663 : box.bigend[2]-box.smallend[2]+1};
1664 : #endif
1665 : }
1666 :
1667 : // Max of lower bound
1668 : [[nodiscard]]
1669 : AMREX_GPU_HOST_DEVICE
1670 : AMREX_FORCE_INLINE
1671 : Dim3 max_lbound (Box const& b1, Box const& b2) noexcept
1672 : {
1673 : #if (AMREX_SPACEDIM == 1)
1674 : return {amrex::max(b1.smallend[0], b2.smallend[0]),
1675 : 0, 0};
1676 : #elif (AMREX_SPACEDIM == 2)
1677 : return {amrex::max(b1.smallend[0], b2.smallend[0]),
1678 : amrex::max(b1.smallend[1], b2.smallend[1]),
1679 : 0};
1680 : #elif (AMREX_SPACEDIM == 3)
1681 : return {amrex::max(b1.smallend[0], b2.smallend[0]),
1682 : amrex::max(b1.smallend[1], b2.smallend[1]),
1683 : amrex::max(b1.smallend[2], b2.smallend[2])};
1684 : #endif
1685 : }
1686 :
1687 : [[nodiscard]]
1688 : AMREX_GPU_HOST_DEVICE
1689 : AMREX_FORCE_INLINE
1690 : Dim3 max_lbound (Box const& b1, Dim3 const& lo) noexcept
1691 : {
1692 : #if (AMREX_SPACEDIM == 1)
1693 : return {amrex::max(b1.smallend[0], lo.x),
1694 : 0, 0};
1695 : #elif (AMREX_SPACEDIM == 2)
1696 : return {amrex::max(b1.smallend[0], lo.x),
1697 : amrex::max(b1.smallend[1], lo.y),
1698 : 0};
1699 : #elif (AMREX_SPACEDIM == 3)
1700 : return {amrex::max(b1.smallend[0], lo.x),
1701 : amrex::max(b1.smallend[1], lo.y),
1702 : amrex::max(b1.smallend[2], lo.z)};
1703 : #endif
1704 : }
1705 :
1706 : // Min of upper bound
1707 : [[nodiscard]]
1708 : AMREX_GPU_HOST_DEVICE
1709 : AMREX_FORCE_INLINE
1710 : Dim3 min_ubound (Box const& b1, Box const& b2) noexcept
1711 : {
1712 : #if (AMREX_SPACEDIM == 1)
1713 : return {amrex::min(b1.bigend[0], b2.bigend[0]),
1714 : 0, 0};
1715 : #elif (AMREX_SPACEDIM == 2)
1716 : return {amrex::min(b1.bigend[0], b2.bigend[0]),
1717 : amrex::min(b1.bigend[1], b2.bigend[1]),
1718 : 0};
1719 : #elif (AMREX_SPACEDIM == 3)
1720 : return {amrex::min(b1.bigend[0], b2.bigend[0]),
1721 : amrex::min(b1.bigend[1], b2.bigend[1]),
1722 : amrex::min(b1.bigend[2], b2.bigend[2])};
1723 : #endif
1724 : }
1725 :
1726 : [[nodiscard]]
1727 : AMREX_GPU_HOST_DEVICE
1728 : AMREX_FORCE_INLINE
1729 : Dim3 min_ubound (Box const& b1, Dim3 const& hi) noexcept
1730 : {
1731 : #if (AMREX_SPACEDIM == 1)
1732 : return {amrex::min(b1.bigend[0], hi.x),
1733 : 0, 0};
1734 : #elif (AMREX_SPACEDIM == 2)
1735 : return {amrex::min(b1.bigend[0], hi.x),
1736 : amrex::min(b1.bigend[1], hi.y),
1737 : 0};
1738 : #elif (AMREX_SPACEDIM == 3)
1739 : return {amrex::min(b1.bigend[0], hi.x),
1740 : amrex::min(b1.bigend[1], hi.y),
1741 : amrex::min(b1.bigend[2], hi.z)};
1742 : #endif
1743 : }
1744 :
1745 : [[nodiscard]]
1746 : AMREX_GPU_HOST_DEVICE
1747 : AMREX_FORCE_INLINE
1748 : Box minBox (Box const& b1, Box const& b2, IndexType typ) noexcept
1749 : {
1750 : #if (AMREX_SPACEDIM == 1)
1751 : return Box(IntVect(amrex::max(b1.smallend[0], b2.smallend[0])),
1752 : IntVect(amrex::min(b1.bigend [0], b2.bigend [0])),
1753 : typ);
1754 : #elif (AMREX_SPACEDIM == 2)
1755 : return Box(IntVect(amrex::max(b1.smallend[0], b2.smallend[0]),
1756 : amrex::max(b1.smallend[1], b2.smallend[1])),
1757 : IntVect(amrex::min(b1.bigend [0], b2.bigend [0]),
1758 : amrex::min(b1.bigend [1], b2.bigend [1])),
1759 : typ);
1760 : #elif (AMREX_SPACEDIM == 3)
1761 : return Box(IntVect(amrex::max(b1.smallend[0], b2.smallend[0]),
1762 : amrex::max(b1.smallend[1], b2.smallend[1]),
1763 : amrex::max(b1.smallend[2], b2.smallend[2])),
1764 : IntVect(amrex::min(b1.bigend [0], b2.bigend [0]),
1765 : amrex::min(b1.bigend [1], b2.bigend [1]),
1766 : amrex::min(b1.bigend [2], b2.bigend [2])),
1767 : typ);
1768 : #endif
1769 : }
1770 :
1771 : // Returns a Box that covers all the argument Boxes in index
1772 : // space. The types are ignored. Thus, the arguments can have
1773 : // different index types, and the returned Box's index type has no
1774 : // meaning.
1775 : [[nodiscard]]
1776 : AMREX_FORCE_INLINE
1777 : Box getIndexBounds (Box const& b1) noexcept
1778 : {
1779 : return b1;
1780 : }
1781 :
1782 : [[nodiscard]]
1783 : AMREX_FORCE_INLINE
1784 : Box getIndexBounds (Box const& b1, Box const& b2) noexcept
1785 : {
1786 : Box b = b1;
1787 : b.setType(b2.ixType());
1788 : b.minBox(b2);
1789 : return b;
1790 : }
1791 :
1792 : template <class T, class ... Ts>
1793 : [[nodiscard]]
1794 : AMREX_FORCE_INLINE
1795 : Box getIndexBounds (T const& b1, T const& b2, Ts const& ... b3) noexcept
1796 : {
1797 : return getIndexBounds(getIndexBounds(b1,b2),b3...);
1798 : }
1799 :
1800 :
1801 : [[nodiscard]]
1802 : AMREX_GPU_HOST_DEVICE
1803 : AMREX_FORCE_INLINE
1804 : IntVect getCell (Box const* boxes, int nboxes, Long icell) noexcept
1805 : {
1806 : int ibox;
1807 : Long ncells_subtotal = 0;
1808 : Long offset = 0;
1809 : for (ibox = 0; ibox < nboxes; ++ibox) {
1810 : const Long n = boxes[ibox].numPts();
1811 : ncells_subtotal += n;
1812 : if (icell < ncells_subtotal) {
1813 : offset = icell - (ncells_subtotal - n);
1814 : break;
1815 : }
1816 : }
1817 : return boxes[ibox].atOffset(offset);
1818 : }
1819 :
1820 : [[nodiscard]]
1821 : AMREX_GPU_HOST_DEVICE
1822 : AMREX_FORCE_INLINE
1823 : Box makeSlab (Box const& b, int direction, int slab_index) noexcept
1824 : {
1825 : Box r = b;
1826 : r.makeSlab(direction,slab_index);
1827 : return r;
1828 : }
1829 :
1830 : [[nodiscard]]
1831 : AMREX_GPU_HOST_DEVICE
1832 : AMREX_FORCE_INLINE
1833 : Box makeSingleCellBox (int i, int j, int k, IndexType typ = IndexType::TheCellType())
1834 : {
1835 : #if (AMREX_SPACEDIM == 1)
1836 : amrex::ignore_unused(j,k);
1837 : #elif (AMREX_SPACEDIM == 2)
1838 : amrex::ignore_unused(k);
1839 : #endif
1840 : return Box(IntVect(AMREX_D_DECL(i,j,k)),IntVect(AMREX_D_DECL(i,j,k)),typ);
1841 : }
1842 :
1843 : struct BoxIndexer
1844 : {
1845 : std::uint64_t npts;
1846 :
1847 : #if (AMREX_SPACEDIM == 3)
1848 : Math::FastDivmodU64 fdxy;
1849 : Math::FastDivmodU64 fdx;
1850 : IntVect lo;
1851 :
1852 : BoxIndexer (Box const& box)
1853 : : npts(box.numPts()),
1854 : fdxy(std::uint64_t(box.length(0))*std::uint64_t(box.length(1))),
1855 : fdx (std::uint64_t(box.length(0))),
1856 : lo (box.smallEnd())
1857 : {}
1858 :
1859 : [[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
1860 : Dim3 operator() (std::uint64_t icell) const
1861 : {
1862 : std::uint64_t x, y, z, rem;
1863 : fdxy(z, rem, icell);
1864 : fdx(y, x, rem);
1865 : return {int(x)+lo[0], int(y)+lo[1], int(z)+lo[2]};
1866 : }
1867 :
1868 : [[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
1869 : IntVect intVect (std::uint64_t icell) const
1870 : {
1871 : std::uint64_t x, y, z, rem;
1872 : fdxy(z, rem, icell);
1873 : fdx(y, x, rem);
1874 : return {int(x)+lo[0], int(y)+lo[1], int(z)+lo[2]};
1875 : }
1876 :
1877 : #elif (AMREX_SPACEDIM == 2)
1878 :
1879 : Math::FastDivmodU64 fdx;
1880 : IntVect lo;
1881 :
1882 : BoxIndexer (Box const& box)
1883 : : npts(box.numPts()),
1884 : fdx (std::uint64_t(box.length(0))),
1885 : lo (box.smallEnd())
1886 : {}
1887 :
1888 : [[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
1889 : Dim3 operator() (std::uint64_t icell) const
1890 : {
1891 : std::uint64_t x, y;
1892 : fdx(y, x, icell);
1893 : return {int(x)+lo[0], int(y)+lo[1], 0};
1894 : }
1895 :
1896 : [[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
1897 : IntVect intVect (std::uint64_t icell) const
1898 : {
1899 : std::uint64_t x, y;
1900 : fdx(y, x, icell);
1901 : return {int(x)+lo[0], int(y)+lo[1]};
1902 : }
1903 :
1904 : #elif (AMREX_SPACEDIM == 1)
1905 :
1906 : int lo;
1907 :
1908 : BoxIndexer (Box const& box)
1909 : : npts(box.numPts()),
1910 : lo(box.smallEnd(0))
1911 : {}
1912 :
1913 : [[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
1914 : Dim3 operator() (std::uint64_t icell) const
1915 : {
1916 : return {int(icell)+lo, 0, 0};
1917 : }
1918 :
1919 : [[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
1920 : IntVect intVect (std::uint64_t icell) const
1921 : {
1922 : return IntVect{int(icell)+lo};
1923 : }
1924 :
1925 : #endif
1926 :
1927 : [[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
1928 : std::uint64_t numPts () const { return npts; }
1929 : };
1930 :
1931 : }
1932 :
1933 : #endif /*AMREX_BOX_H*/
|