Block-Structured AMR Software Framework
 
Loading...
Searching...
No Matches
AMReX_Box.H
Go to the documentation of this file.
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_BoxIterator.H>
9#include <AMReX_ccse-mpi.H>
10#include <AMReX_IntVect.H>
11#include <AMReX_IndexType.H>
12#include <AMReX_Orientation.H>
13#include <AMReX_SPACE.H>
14#include <AMReX_Array.H>
15#include <AMReX_Array4.H>
16#include <AMReX_Vector.H>
17#include <AMReX_GpuQualifiers.H>
18#include <AMReX_GpuControl.H>
19#include <AMReX_Math.H>
20
21#include <iosfwd>
22
23namespace amrex
24{
25template<int dim>
26class BoxND;
28using Box = BoxND<AMREX_SPACEDIM>;
29class BoxCommHelper;
30template<int dim>
31class BoxIteratorND;
32
47template<int dim>
48class BoxND
49{
51 friend class BoxCommHelper;
52
53public:
54 /*
55 * \brief The default constructor. For safety, the constructed BoxND is
56 * invalid and may be tested for validity with ok().
57 * DO NOT CHANGE THIS BEHAVIOR!
58 */
60 constexpr BoxND () noexcept
61 : smallend(1),
62 bigend(0)
63 {}
64
67 constexpr BoxND (const IntVectND<dim>& small, const IntVectND<dim>& big) noexcept
68 : smallend(small),
69 bigend(big)
70 {}
71
74 BoxND (const IntVectND<dim>& small, const int* vec_len) noexcept
75 : smallend(small),
76 bigend(small + IntVectND<dim>(vec_len) - 1)
77 {}
78
84 BoxND (const IntVectND<dim>& small, const IntVectND<dim>& big, const IntVectND<dim>& typ) noexcept
85 : smallend(small),
86 bigend(big),
87 btype(typ)
88 {
89 BL_ASSERT(typ.allGE(0) && typ.allLE(1));
90 }
91
94 BoxND (const IntVectND<dim>& small, const IntVectND<dim>& big, IndexTypeND<dim> t) noexcept
95 : smallend(small),
96 bigend(big),
97 btype(t)
98 {}
99
100 template <typename T, int Tdim=dim, std::enable_if_t<( 1<=Tdim && Tdim<=3 ), int> = 0>
101 AMREX_GPU_HOST_DEVICE
102 explicit BoxND (Array4<T> const& a) noexcept
103 : smallend(a.begin),
104 bigend(IntVectND<dim>(a.end) - 1)
105 {}
106
107 // dtor, copy-ctor, copy-op=, move-ctor, and move-op= are compiler generated.
108
110 [[nodiscard]] AMREX_GPU_HOST_DEVICE
111 const IntVectND<dim>& smallEnd () const& noexcept { return smallend; }
112
114 [[nodiscard]] const IntVectND<dim>& smallEnd () && = delete;
116
118 [[nodiscard]] AMREX_GPU_HOST_DEVICE
119 int smallEnd (int dir) const& noexcept { return smallend[dir]; }
120
122 [[nodiscard]] AMREX_GPU_HOST_DEVICE
123 const IntVectND<dim>& bigEnd () const& noexcept { return bigend; }
124
126 [[nodiscard]] const IntVectND<dim>& bigEnd () && = delete;
128
130 [[nodiscard]] AMREX_GPU_HOST_DEVICE
131 int bigEnd (int dir) const noexcept { return bigend[dir]; }
132
134 [[nodiscard]] AMREX_GPU_HOST_DEVICE
135 IndexTypeND<dim> ixType () const noexcept { return btype; }
136
138 [[nodiscard]] AMREX_GPU_HOST_DEVICE
139 IntVectND<dim> type () const noexcept { return btype.ixType(); }
140
142 [[nodiscard]] AMREX_GPU_HOST_DEVICE
143 IndexType::CellIndex type (int dir) const noexcept { return btype.ixType(dir); }
144
146 [[nodiscard]] AMREX_GPU_HOST_DEVICE
147 IntVectND<dim> size () const noexcept
148 {
149 return bigend - smallend + 1;
150 }
151
153 [[nodiscard]] AMREX_GPU_HOST_DEVICE
154 IntVectND<dim> length () const noexcept
155 {
156 return bigend - smallend + 1;
157 }
158
160 [[nodiscard]] AMREX_GPU_HOST_DEVICE
161 int length (int dir) const noexcept { return bigend[dir] - smallend[dir] + 1; }
162
163 template <int N=dim, std::enable_if_t<( 1<=N && N<=3 ), int> = 0>
164 [[nodiscard]] AMREX_GPU_HOST_DEVICE
165 GpuArray<int,3> length3d () const noexcept {
166 Dim3 len3d = length().dim3(1);
167 return {{len3d.x, len3d.y, len3d.z}};
168 }
169
170 template <int N=dim, std::enable_if_t<( 1<=N && N<=3 ), int> = 0>
171 [[nodiscard]] AMREX_GPU_HOST_DEVICE
172 GpuArray<int,3> loVect3d () const noexcept {
173 Dim3 lo3d = smallend.dim3(0);
174 return {{lo3d.x, lo3d.y, lo3d.z}};
175 }
176
177 template <int N=dim, std::enable_if_t<( 1<=N && N<=3 ), int> = 0>
178 [[nodiscard]] AMREX_GPU_HOST_DEVICE
179 GpuArray<int,3> hiVect3d () const noexcept {
180 Dim3 hi3d = bigend.dim3(0);
181 return {{hi3d.x, hi3d.y, hi3d.z}};
182 }
183
185 [[nodiscard]] AMREX_GPU_HOST_DEVICE
186 const int* loVect () const& noexcept { return smallend.getVect(); }
187 AMREX_GPU_HOST_DEVICE
188 const int* loVect () && = delete;
190 [[nodiscard]] AMREX_GPU_HOST_DEVICE
191 const int* hiVect () const& noexcept { return bigend.getVect(); }
192 AMREX_GPU_HOST_DEVICE
193 const int* hiVect () && = delete;
194
196 [[nodiscard]] AMREX_GPU_HOST_DEVICE
197 int operator[] (Orientation face) const noexcept {
198 const int dir = face.coordDir();
199 return face.isLow() ? smallend[dir] : bigend[dir];
200 }
201
203 [[nodiscard]] AMREX_GPU_HOST_DEVICE
204 bool isEmpty () const noexcept { return !ok(); }
205
207 [[nodiscard]] AMREX_GPU_HOST_DEVICE
208 bool ok () const noexcept { return bigend.allGE(smallend) && btype.ok(); }
209
211 [[nodiscard]] AMREX_GPU_HOST_DEVICE
212 bool contains (const IntVectND<dim>& p) const noexcept {
213 return p.allGE(smallend) && p.allLE(bigend);
214 }
215
217 template <int N=dim, std::enable_if_t<( 1<=N && N<=3 ), int> = 0>
218 [[nodiscard]] AMREX_GPU_HOST_DEVICE
219 bool contains (const Dim3& p) const noexcept {
220 IntVectND<dim> piv{p};
221 return contains(piv);
222 }
223
225 template <int N=dim, std::enable_if_t<( 1<=N && N<=3 ), int> = 0>
226 [[nodiscard]] AMREX_GPU_HOST_DEVICE
227 bool contains (int i, int j, int k) const noexcept {
228 Dim3 p3d{i, j, k};
229 return contains(p3d);
230 }
231
235 [[nodiscard]] AMREX_GPU_HOST_DEVICE
236 bool contains (const BoxND& b) const noexcept
237 {
238 BL_ASSERT(sameType(b));
239 return b.smallend.allGE(smallend) && b.bigend.allLE(bigend);
240 }
241
243 [[nodiscard]] AMREX_GPU_HOST_DEVICE
244 bool strictly_contains (const IntVectND<dim>& p) const noexcept {
245 return p.allGT(smallend) && p.allLT(bigend);
246 }
247
252 [[nodiscard]] AMREX_GPU_HOST_DEVICE
253 bool strictly_contains (const BoxND& b) const noexcept
254 {
255 BL_ASSERT(sameType(b));
256 return b.smallend.allGT(smallend) && b.bigend.allLT(bigend);
257 }
258
260 template <int N=dim, std::enable_if_t<( 1<=N && N<=3 ), int> = 0>
261 [[nodiscard]] AMREX_GPU_HOST_DEVICE
262 bool strictly_contains (const Dim3& p) const noexcept {
263 IntVectND<dim> piv{p};
264 return strictly_contains(piv);
265 }
266
268 template <int N=dim, std::enable_if_t<( 1<=N && N<=3 ), int> = 0>
269 [[nodiscard]] AMREX_GPU_HOST_DEVICE
270 bool strictly_contains (int i, int j, int k) const noexcept {
271 Dim3 p3d{i, j, k};
272 return strictly_contains(p3d);
273 }
274
279 [[nodiscard]] AMREX_GPU_HOST_DEVICE
280 bool intersects (const BoxND& b) const noexcept { BoxND isect(*this); isect &= b; return isect.ok(); }
281
286 [[nodiscard]] AMREX_GPU_HOST_DEVICE
287 bool sameSize (const BoxND& b) const noexcept {
288 BL_ASSERT(sameType(b));
289 return length() == b.length();
290 }
291
293 [[nodiscard]] AMREX_GPU_HOST_DEVICE
294 bool sameType (const BoxND &b) const noexcept { return btype == b.btype; }
295
297 [[nodiscard]] AMREX_GPU_HOST_DEVICE
298 bool operator== (const BoxND& b) const noexcept { return smallend == b.smallend && bigend == b.bigend && b.btype == btype; }
299
301 [[nodiscard]] AMREX_GPU_HOST_DEVICE
302 bool operator!= (const BoxND& b) const noexcept { return !operator==(b); }
303
304 [[nodiscard]] AMREX_GPU_HOST_DEVICE
305 bool operator< (const BoxND& rhs) const noexcept
306 {
307 return btype < rhs.btype ||
308 ((btype == rhs.btype) &&
309 ( (smallend < rhs.smallend) ||
310 ((smallend == rhs.smallend) && (bigend < rhs.bigend)) ));
311 }
312 [[nodiscard]] AMREX_GPU_HOST_DEVICE
313 bool operator <= (const BoxND& rhs) const noexcept {
314 return !(rhs < *this);
315 }
316 [[nodiscard]] AMREX_GPU_HOST_DEVICE
317 bool operator> (const BoxND& rhs) const noexcept {
318 return rhs < *this;
319 }
320 [[nodiscard]] AMREX_GPU_HOST_DEVICE
321 bool operator>= (const BoxND& rhs) const noexcept {
322 return !(*this < rhs);
323 }
324
326 [[nodiscard]] AMREX_GPU_HOST_DEVICE
327 bool cellCentered () const noexcept { return !btype.any(); }
328
331 void checkOverflow () const noexcept {
332 if (ok()) {
333 for (int i = 0; i < dim; ++i) {
334 auto lo = static_cast<Long>(smallend[i]);
335 auto hi = static_cast<Long>(bigend[i]);
336 Long len = hi - lo + 1;
337 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(len>=0 && len<std::numeric_limits<int>::max(),
338 "Overflow when computing length of box");
339 }
340 auto num_pts = static_cast<Long>(length(0));
341 for (int i = 1; i < dim; ++i) {
342 auto len = static_cast<Long>(length(i));
343 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(num_pts == 0 || len == 0 ||
344 num_pts <= std::numeric_limits<Long>::max() / len,
345 "Overflow when computing numPts of box");
346 num_pts *= len;
347 }
348 }
349 }
351
355 [[nodiscard]] AMREX_GPU_HOST_DEVICE
356 Long numPts () const noexcept {
357#if defined(AMREX_DEBUG) || defined(AMREX_USE_ASSERTION)
358 AMREX_IF_ON_HOST((checkOverflow();))
359#endif
360 if (ok()) {
361 auto num_pts = static_cast<Long>(length(0));
362 for (int i = 1; i < dim; ++i) {
363 num_pts *= static_cast<Long>(length(i));
364 }
365 return num_pts;
366 } else {
367 return Long(0);
368 }
369 }
370
375 [[nodiscard]] AMREX_GPU_HOST_DEVICE
376 double d_numPts () const noexcept {
377 if (ok()) {
378 auto num_pts = static_cast<double>(length(0));
379 for (int i = 1; i < dim; ++i) {
380 num_pts *= static_cast<double>(length(i));
381 }
382 return num_pts;
383 } else {
384 return 0.0;
385 }
386 }
387
393 [[nodiscard]] AMREX_GPU_HOST_DEVICE
394 Long volume () const noexcept {
395 if (ok()) {
396 auto num_pts = static_cast<Long>(length(0)-btype[0]);
397 for (int i = 1; i < dim; ++i) {
398 num_pts *= static_cast<Long>(length(i)-btype[i]);
399 }
400 return num_pts;
401 } else {
402 return Long(0);
403 }
404 }
405
410 [[nodiscard]] AMREX_GPU_HOST_DEVICE
411 int longside (int& dir) const noexcept {
412 int maxlen = length(0);
413 dir = 0;
414 for (int i = 1; i < dim; i++)
415 {
416 if (length(i) > maxlen)
417 {
418 maxlen = length(i);
419 dir = i;
420 }
421 }
422 return maxlen;
423 }
424
426 [[nodiscard]] AMREX_GPU_HOST_DEVICE
427 int longside () const noexcept {
428 int ignore = 0;
429 return longside(ignore);
430 }
431
436 [[nodiscard]] AMREX_GPU_HOST_DEVICE
437 int shortside (int& dir) const noexcept {
438 int minlen = length(0);
439 dir = 0;
440 for (int i = 1; i < dim; i++)
441 {
442 if (length(i) < minlen)
443 {
444 minlen = length(i);
445 dir = i;
446 }
447 }
448 return minlen;
449 }
450
452 [[nodiscard]] AMREX_GPU_HOST_DEVICE
453 int shortside () const noexcept {
454 int ignore = 0;
455 return shortside(ignore);
456 }
457
463 [[nodiscard]] AMREX_GPU_HOST_DEVICE
464 Long index (const IntVectND<dim>& v) const noexcept;
465
467 [[nodiscard]] AMREX_GPU_HOST_DEVICE
468 IntVectND<dim> atOffset (Long offset) const noexcept;
469
473 template <int N=dim, std::enable_if_t<( 1<=N && N<=3 ), int> = 0>
474 [[nodiscard]] AMREX_GPU_HOST_DEVICE
475 GpuArray<int,3> atOffset3d (Long offset) const noexcept;
476
478 AMREX_GPU_HOST_DEVICE
479 BoxND& setSmall (const IntVectND<dim>& sm) noexcept { smallend = sm; return *this; }
480
482 AMREX_GPU_HOST_DEVICE
483 BoxND& setSmall (int dir, int sm_index) noexcept { smallend.setVal(dir,sm_index); return *this; }
484
486 AMREX_GPU_HOST_DEVICE
487 BoxND& setBig (const IntVectND<dim>& bg) noexcept { bigend = bg; return *this; }
488
490 AMREX_GPU_HOST_DEVICE
491 BoxND& setBig (int dir, int bg_index) noexcept { bigend.setVal(dir,bg_index); return *this; }
492
498 AMREX_GPU_HOST_DEVICE
499 BoxND& setRange (int dir,
500 int sm_index,
501 int n_cells = 1) noexcept;
502
504 AMREX_GPU_HOST_DEVICE
505 BoxND& setType (const IndexTypeND<dim>& t) noexcept { btype = t; return *this; }
506
508 AMREX_GPU_HOST_DEVICE
509 BoxND& shift (int dir, int nzones) noexcept { smallend.shift(dir,nzones); bigend.shift(dir,nzones); return *this; }
510
512 AMREX_GPU_HOST_DEVICE
513 BoxND& shift (const IntVectND<dim>& iv) noexcept { smallend.shift(iv); bigend.shift(iv); return *this; }
514
524 AMREX_GPU_HOST_DEVICE
525 BoxND& shiftHalf (int dir, int num_halfs) noexcept;
526
528 AMREX_GPU_HOST_DEVICE
529 BoxND& shiftHalf (const IntVectND<dim>& iv) noexcept;
530
538 AMREX_GPU_HOST_DEVICE
539 BoxND& convert (IndexTypeND<dim> typ) noexcept;
540
548 AMREX_GPU_HOST_DEVICE
549 BoxND& convert (const IntVectND<dim>& typ) noexcept;
550
552 AMREX_GPU_HOST_DEVICE
553 BoxND& surroundingNodes () noexcept;
554
556 AMREX_GPU_HOST_DEVICE
557 BoxND& surroundingNodes (int dir) noexcept;
558
559 AMREX_GPU_HOST_DEVICE
560 BoxND& surroundingNodes (Direction d) noexcept { return surroundingNodes(static_cast<int>(d)); }
561
563 AMREX_GPU_HOST_DEVICE
564 BoxND& enclosedCells () noexcept;
565
567 AMREX_GPU_HOST_DEVICE
568 BoxND& enclosedCells (int dir) noexcept;
569
571 AMREX_GPU_HOST_DEVICE
572 BoxND& enclosedCells (Direction d) noexcept { return enclosedCells(static_cast<int>(d)); }
573
578 AMREX_GPU_HOST_DEVICE
579 BoxND operator& (const BoxND& rhs) const noexcept { BoxND lhs(*this); lhs &= rhs; return lhs; }
580
582 AMREX_GPU_HOST_DEVICE
583 BoxND& operator&= (const BoxND& rhs) noexcept
584 {
585 BL_ASSERT(sameType(rhs));
586 smallend.max(rhs.smallend);
587 bigend.min(rhs.bigend);
588 return *this;
589 }
590
596 AMREX_GPU_HOST_DEVICE
597 BoxND& minBox (const BoxND& b) noexcept {
598 // BoxArray may call this with not ok boxes. BL_ASSERT(b.ok() && ok());
599 BL_ASSERT(sameType(b));
600 smallend.min(b.smallend);
601 bigend.max(b.bigend);
602 return *this;
603 }
604
606 AMREX_GPU_HOST_DEVICE
607 BoxND& operator+= (const IntVectND<dim>& v) noexcept { smallend += v; bigend += v; return *this; }
608
610 AMREX_GPU_HOST_DEVICE
611 BoxND operator+ (const IntVectND<dim>& v) const noexcept { BoxND r(*this); r += v; return r; }
612
614 AMREX_GPU_HOST_DEVICE
615 BoxND& operator-= (const IntVectND<dim>& v) noexcept { smallend -= v; bigend -= v; return *this; }
616
618 AMREX_GPU_HOST_DEVICE
619 BoxND operator- (const IntVectND<dim>& v) const noexcept { BoxND r(*this); r -= v; return r; }
620
633 AMREX_GPU_HOST_DEVICE
634 BoxND chop (int dir, int chop_pnt) noexcept;
635
636 /*
637 * \brief Grow BoxND in all directions by given amount.
638 * NOTE: n_cell negative shrinks the BoxND by that number of cells.
639 */
640 AMREX_GPU_HOST_DEVICE
641 BoxND& grow (int i) noexcept { smallend.diagShift(-i); bigend.diagShift(i); return *this; }
642
644 AMREX_GPU_HOST_DEVICE
645 BoxND& grow (const IntVectND<dim>& v) noexcept { smallend -= v; bigend += v; return *this;}
646
651 AMREX_GPU_HOST_DEVICE
652 BoxND& grow (int idir, int n_cell) noexcept { smallend.shift(idir, -n_cell); bigend.shift(idir, n_cell); return *this; }
653
654 AMREX_GPU_HOST_DEVICE
655 BoxND& grow (Direction d, int n_cell) noexcept { return grow(static_cast<int>(d), n_cell); }
656
661 AMREX_GPU_HOST_DEVICE
662 BoxND& growLo (int idir, int n_cell = 1) noexcept { smallend.shift(idir, -n_cell); return *this; }
663
664 AMREX_GPU_HOST_DEVICE
665 BoxND& growLo (Direction d, int n_cell = 1) noexcept { return growLo(static_cast<int>(d), n_cell); }
666
672 AMREX_GPU_HOST_DEVICE
673 BoxND& growHi (int idir, int n_cell = 1) noexcept { bigend.shift(idir,n_cell); return *this; }
674
675 AMREX_GPU_HOST_DEVICE
676 BoxND& growHi (Direction d, int n_cell = 1) noexcept { return growHi(static_cast<int>(d), n_cell); }
677
679 AMREX_GPU_HOST_DEVICE
680 BoxND& grow (Orientation face, int n_cell = 1) noexcept {
681 int idir = face.coordDir();
682 if (face.isLow()) {
683 smallend.shift(idir, -n_cell);
684 } else {
685 bigend.shift(idir,n_cell);
686 }
687 return *this;
688 }
689
697 AMREX_GPU_HOST_DEVICE
698 BoxND& refine (int ref_ratio) noexcept {
699 return this->refine(IntVectND<dim>(ref_ratio));
700 }
701
702 /*
703 * \brief Refine BoxND by given (positive) refinement ratio.
704 * NOTE: if type(dir) = CELL centered: lo <- lo*ratio and
705 * hi <- (hi+1)*ratio - 1.
706 * NOTE: if type(dir) = NODE centered: lo <- lo*ratio and
707 * hi <- hi*ratio.
708 */
709 AMREX_GPU_HOST_DEVICE
710 BoxND& refine (const IntVectND<dim>& ref_ratio) noexcept;
711
721 AMREX_GPU_HOST_DEVICE
722 BoxND& coarsen (int ref_ratio) noexcept {
723 return this->coarsen(IntVectND<dim>(ref_ratio));
724 }
725
736 BoxND& coarsen (const IntVectND<dim>& ref_ratio) noexcept;
737
743 void next (IntVectND<dim> &) const noexcept;
744
754
755 [[nodiscard]] AMREX_GPU_HOST_DEVICE
756 bool isSquare() const noexcept;
757
769 [[nodiscard]] AMREX_GPU_HOST_DEVICE
770 bool coarsenable (const IntVectND<dim>& refrat, const IntVectND<dim>& min_width) const noexcept
771 {
772 if (!size().allGE(refrat*min_width)) {
773 return false;
774 } else {
775 BoxND testBox = *this;
776 testBox.coarsen(refrat);
777 testBox.refine (refrat);
778 return (*this == testBox);
779 }
780 }
781
793 [[nodiscard]] AMREX_GPU_HOST_DEVICE
794 bool coarsenable (int refrat, int min_width=1) const noexcept {
795 return coarsenable(IntVectND<dim>(refrat), IntVectND<dim>(min_width));
796 }
797
809 [[nodiscard]] AMREX_GPU_HOST_DEVICE
810 bool coarsenable (const IntVectND<dim>& refrat, int min_width=1) const noexcept
811 {
812 return coarsenable(refrat, IntVectND<dim>(min_width));
813 }
814
816 void normalize () noexcept
817 {
818 for (int idim=0; idim < dim; ++idim) {
819 if (this->length(idim) == 0) {
820 this->growHi(idim,1);
821 }
822 }
823 }
824
826 BoxND& makeSlab (int direction, int slab_index) noexcept
827 {
828 smallend[direction] = slab_index;
829 bigend[direction] = slab_index;
830 return *this;
831 }
832
844 [[nodiscard]] BoxIteratorND<dim> iterator () const noexcept {
845 return BoxIteratorND<dim>{*this};
846 }
847
849 static constexpr std::size_t ndims () noexcept {
850 return static_cast<std::size_t>(dim);
851 }
852
854 static constexpr int indims () noexcept {
855 return dim;
856 }
857
862 template<int new_dim>
864 BoxND<new_dim> shrink () const noexcept {
865 static_assert(new_dim <= dim);
866 auto lo = smallend.template shrink<new_dim>();
867 auto hi = bigend.template shrink<new_dim>();
868 auto typ = btype.template shrink<new_dim>();
869 return BoxND<new_dim>(lo, hi, typ);
870 }
871
877 template<int new_dim>
879 BoxND<new_dim> expand () const noexcept {
880 static_assert(new_dim >= dim);
881 auto lo = smallend.template expand<new_dim>(0);
882 auto hi = bigend.template expand<new_dim>(0);
883 auto typ = btype.template expand<new_dim>(IndexType::CellIndex::CELL);
884 return BoxND<new_dim>(lo, hi, typ);
885 }
886
891 template<int new_dim>
893 BoxND<new_dim> resize () const noexcept {
894 if constexpr (new_dim > dim) {
895 return expand<new_dim>();
896 } else {
897 return shrink<new_dim>();
898 }
899 }
900
901private:
905};
906
907template<int dim>
911BoxND<dim>::refine (const IntVectND<dim>& ref_ratio) noexcept
912{
913 if (ref_ratio != 1) {
914 IntVectND<dim> shft(1);
915 shft -= btype.ixType();
916 smallend *= ref_ratio;
917 bigend += shft;
918 bigend *= ref_ratio;
919 bigend -= shft;
920 }
921 return *this;
922}
923
924template<int dim>
928BoxND<dim>::coarsen (const IntVectND<dim>& ref_ratio) noexcept
929{
930 if (ref_ratio != 1)
931 {
932 smallend.coarsen(ref_ratio);
933
934 if (btype.any())
935 {
936 IntVectND<dim> off(0);
937 for (int dir = 0; dir < dim; dir++)
938 {
939 if (btype[dir]) {
940 if (bigend[dir]%ref_ratio[dir]) {
941 off.setVal(dir,1);
942 }
943 }
944 }
945 bigend.coarsen(ref_ratio);
946 bigend += off;
947 }
948 else
949 {
950 bigend.coarsen(ref_ratio);
951 }
952 }
953
954 return *this;
955}
956
957template<int dim>
962{
963 BL_ASSERT(typ.allGE(0) && typ.allLE(1));
964 IntVectND<dim> shft(typ - btype.ixType());
965 bigend += shft;
966 btype = IndexTypeND<dim>(typ);
967 return *this;
968}
969
970template<int dim>
975{
976 for (int dir = 0; dir < dim; dir++)
977 {
978 const auto typ = t[dir];
979 const auto bitval = btype[dir];
980 const int off = typ - bitval;
981 bigend.shift(dir,off);
982 btype.setType(dir, static_cast<IndexType::CellIndex>(typ));
983 }
984 return *this;
985}
986
987template<int dim>
992{
993 if (!(btype[dir]))
994 {
995 bigend.shift(dir,1);
996 //
997 // Set dir'th bit to 1 = IndexType::NODE.
998 //
999 btype.set(dir);
1000 }
1001 return *this;
1002}
1003
1004template<int dim>
1009{
1010 for (int i = 0; i < dim; ++i) {
1011 if ((btype[i] == 0)) {
1012 bigend.shift(i,1);
1013 }
1014 }
1015 btype.setall();
1016 return *this;
1017}
1018
1019template<int dim>
1024{
1025 if (btype[dir])
1026 {
1027 bigend.shift(dir,-1);
1028 //
1029 // Set dir'th bit to 0 = IndexType::CELL.
1030 //
1031 btype.unset(dir);
1032 }
1033 return *this;
1034}
1035
1036template<int dim>
1041{
1042 for (int i = 0 ; i < dim; ++i) {
1043 if (btype[i]) {
1044 bigend.shift(i,-1);
1045 }
1046 }
1047 btype.clear();
1048 return *this;
1049}
1050
1051template<int dim>
1054Long
1055BoxND<dim>::index (const IntVectND<dim>& v) const noexcept
1056{
1057 IntVectND<dim> vz = v - smallend;
1058 Long result = vz[0];
1059 Long mult = length(0);
1060 for (int i = 1 ; i < dim; ++i) {
1061 result += mult * vz[i];
1062 mult *= length(i);
1063 }
1064 return result;
1065}
1066
1067template<int dim>
1072{
1073 IntVectND<dim> result = smallend;
1074
1075 if constexpr (dim > 1) {
1076 GpuArray<Long, dim-1> mult{};
1077 mult[0] = length(0);
1078 for (int i = 1 ; i < dim-1; ++i) {
1079 mult[i] = mult[i-1] * length(i);
1080 }
1081 for (int i = dim-1 ; i > 0; --i) {
1082 Long idx = offset / mult[i-1];
1083 offset -= idx * mult[i-1];
1084 result[i] += static_cast<int>(idx);
1085 }
1086 }
1087
1088 result[0] += static_cast<int>(offset);
1089
1090 return result;
1091}
1092
1093template<int dim>
1094template <int N, std::enable_if_t<( 1<=N && N<=3 ), int>>
1095AMREX_GPU_HOST_DEVICE
1096AMREX_FORCE_INLINE
1097GpuArray<int,3>
1098BoxND<dim>::atOffset3d (Long offset) const noexcept
1099{
1100 Dim3 iv3d = atOffset(offset).dim3(0);
1101 return {{iv3d.x, iv3d.y, iv3d.z}};
1102}
1103
1104template<int dim>
1109 int sm_index,
1110 int n_cells) noexcept
1111{
1112 smallend.setVal(dir,sm_index);
1113 bigend.setVal(dir,sm_index+n_cells-1);
1114 return *this;
1115}
1116
1117template<int dim>
1120void
1121BoxND<dim>::next (IntVectND<dim>& p) const noexcept // NOLINT(readability-convert-member-functions-to-static)
1122{
1123 BL_ASSERT(contains(p));
1124
1125 ++p[0];
1126
1127 for (int i = 0 ; i < dim-1; ++i) {
1128 if (p[i] > bigend[i]) {
1129 p[i] = smallend[i];
1130 ++p[i+1];
1131 } else {
1132 break;
1133 }
1134 }
1135}
1136
1137template<int dim>
1140bool
1141BoxND<dim>::isSquare () const noexcept // NOLINT(readability-convert-member-functions-to-static)
1142{
1143 if constexpr (dim == 1) {
1144 return false; // can't build a square in 1-D
1145 } else {
1146 bool is_square = true;
1147 const IntVectND<dim>& sz = this->size();
1148 for (int i = 0 ; i < dim-1; ++i) {
1149 is_square = is_square && (sz[i] == sz[i+1]);
1150 }
1151 return is_square;
1152 }
1153}
1154
1155//
1156// Modified BoxND is low end, returned BoxND is high end.
1157// If CELL: chop_pnt included in high end.
1158// If NODE: chop_pnt included in both Boxes.
1159//
1160
1161template<int dim>
1163inline
1165BoxND<dim>::chop (int dir, int chop_pnt) noexcept
1166{
1167 //
1168 // Define new high end BoxND including chop_pnt.
1169 //
1170 IntVectND<dim> sm(smallend);
1171 IntVectND<dim> bg(bigend);
1172 sm.setVal(dir,chop_pnt);
1173 if (btype[dir])
1174 {
1175 //
1176 // NODE centered BoxND.
1177 //
1178 BL_ASSERT(chop_pnt > smallend[dir] && chop_pnt < bigend[dir]);
1179 //
1180 // Shrink original BoxND to just contain chop_pnt.
1181 //
1182 bigend.setVal(dir,chop_pnt);
1183 }
1184 else
1185 {
1186 //
1187 // CELL centered BoxND.
1188 //
1189 BL_ASSERT(chop_pnt > smallend[dir] && chop_pnt <= bigend[dir]);
1190 //
1191 // Shrink original BoxND to one below chop_pnt.
1192 //
1193 bigend.setVal(dir,chop_pnt-1);
1194 }
1195 return BoxND<dim>(sm,bg,btype);
1196}
1197
1198template<int dim>
1200inline
1202BoxND<dim>::shiftHalf (int dir, int num_halfs) noexcept
1203{
1204 const int nbit = (num_halfs<0 ? -num_halfs : num_halfs)%2;
1205 int nshift = num_halfs/2;
1206 //
1207 // Toggle btyp bit if num_halfs is odd.
1208 //
1209 const unsigned int bit_dir = btype[dir];
1210 if (nbit) {
1211 btype.flip(dir);
1212 }
1213 if (num_halfs < 0) {
1214 nshift -= (bit_dir ? nbit : 0);
1215 } else {
1216 nshift += (bit_dir ? 0 : nbit);
1217 }
1218 smallend.shift(dir,nshift);
1219 bigend.shift(dir,nshift);
1220 return *this;
1221}
1222
1223template<int dim>
1225inline
1228{
1229 for (int i = 0; i < dim; i++) {
1230 shiftHalf(i,iv[i]);
1231 }
1232 return *this;
1233}
1234
1236class BoxCommHelper
1237{
1238public:
1239
1240 explicit BoxCommHelper (const Box& bx, int* p_ = nullptr);
1241
1242 [[nodiscard]] int* data () const& noexcept { return p; }
1243 int* data () && = delete;
1244
1245 [[nodiscard]] Box make_box () const noexcept {
1246 return Box(IntVect(p), IntVect(p+AMREX_SPACEDIM), IntVect(p+2*AMREX_SPACEDIM));
1247 }
1248
1249 [[nodiscard]] static int size () noexcept { return 3*AMREX_SPACEDIM; }
1250
1251private:
1252 int* p;
1253 std::vector<int> v;
1254};
1256
1257class BoxConverter { // NOLINT
1258public:
1259 virtual Box doit (const Box& fine) const = 0; // NOLINT
1260 virtual BoxConverter* clone () const = 0; // NOLINT
1261 virtual ~BoxConverter () = default;
1262};
1263
1264void AllGatherBoxes (Vector<Box>& bxs, int n_extra_reserve=0);
1265
1276template<int dim>
1277[[nodiscard]]
1280BoxND<dim> grow (const BoxND<dim>& b, int i) noexcept
1281{
1282 BoxND<dim> result = b;
1283 result.grow(i);
1284 return result;
1285}
1286
1297template<int dim>
1298[[nodiscard]]
1301BoxND<dim> grow (const BoxND<dim>& b, const IntVectND<dim>& v) noexcept
1302{
1303 BoxND<dim> result = b;
1304 result.grow(v);
1305 return result;
1306}
1307
1319template<int dim>
1320[[nodiscard]]
1323BoxND<dim> grow (const BoxND<dim>& b, int idir, int n_cell) noexcept
1324{
1325 BoxND<dim> result = b;
1326 result.grow(idir, n_cell);
1327 return result;
1328}
1329
1341template<int dim>
1342[[nodiscard]]
1345BoxND<dim> grow (const BoxND<dim>& b, Direction d, int n_cell) noexcept
1346{
1347 return grow(b, static_cast<int>(d), n_cell);
1348}
1349
1350template<int dim>
1351[[nodiscard]]
1354BoxND<dim> growLo (const BoxND<dim>& b, int idir, int n_cell) noexcept
1355{
1356 BoxND<dim> result = b;
1357 result.growLo(idir, n_cell);
1358 return result;
1359}
1360
1361template<int dim>
1362[[nodiscard]]
1365BoxND<dim> growLo (const BoxND<dim>& b, Direction d, int n_cell) noexcept
1366{
1367 return growLo(b, static_cast<int>(d), n_cell);
1368}
1369
1370template<int dim>
1371[[nodiscard]]
1374BoxND<dim> growHi (const BoxND<dim>& b, int idir, int n_cell) noexcept
1375{
1376 BoxND<dim> result = b;
1377 result.growHi(idir, n_cell);
1378 return result;
1379}
1380
1381template<int dim>
1382[[nodiscard]]
1385BoxND<dim> growHi (const BoxND<dim>& b, Direction d, int n_cell) noexcept
1386{
1387 return growHi(b, static_cast<int>(d), n_cell);
1388}
1389
1405template<int dim>
1406[[nodiscard]]
1409BoxND<dim> coarsen (const BoxND<dim>& b, int ref_ratio) noexcept
1410{
1411 BoxND<dim> result = b;
1412 result.coarsen(IntVectND<dim>(ref_ratio));
1413 return result;
1414}
1415
1431template<int dim>
1432[[nodiscard]]
1435BoxND<dim> coarsen (const BoxND<dim>& b, const IntVectND<dim>& ref_ratio) noexcept
1436{
1437 BoxND<dim> result = b;
1438 result.coarsen(ref_ratio);
1439 return result;
1440}
1441
1455template<int dim>
1456[[nodiscard]]
1459BoxND<dim> refine (const BoxND<dim>& b, int ref_ratio) noexcept
1460{
1461 BoxND<dim> result = b;
1462 result.refine(IntVectND<dim>(ref_ratio));
1463 return result;
1464}
1465
1479template<int dim>
1480[[nodiscard]]
1483BoxND<dim> refine (const BoxND<dim>& b, const IntVectND<dim>& ref_ratio) noexcept
1484{
1485 BoxND<dim> result = b;
1486 result.refine(ref_ratio);
1487 return result;
1488}
1489
1491template<int dim>
1492[[nodiscard]]
1495BoxND<dim> shift (const BoxND<dim>& b, int dir, int nzones) noexcept
1496{
1497 BoxND<dim> result = b;
1498 result.shift(dir, nzones);
1499 return result;
1500}
1501
1502template<int dim>
1503[[nodiscard]]
1506BoxND<dim> shift (const BoxND<dim>& b, const IntVectND<dim>& nzones) noexcept
1507{
1508 BoxND<dim> result = b;
1509 result.shift(nzones);
1510 return result;
1511}
1512
1518template<int dim>
1519[[nodiscard]]
1522BoxND<dim> surroundingNodes (const BoxND<dim>& b, int dir) noexcept
1523{
1524 BoxND<dim> bx(b);
1525 bx.surroundingNodes(dir);
1526 return bx;
1527}
1528
1529template<int dim>
1530[[nodiscard]]
1534{
1535 return surroundingNodes(b, static_cast<int>(d));
1536}
1537
1542template<int dim>
1543[[nodiscard]]
1547{
1548 BoxND<dim> bx(b);
1549 bx.surroundingNodes();
1550 return bx;
1551}
1552
1554template<int dim>
1555[[nodiscard]]
1558BoxND<dim> convert (const BoxND<dim>& b, const IntVectND<dim>& typ) noexcept
1559{
1560 BoxND<dim> bx(b);
1561 bx.convert(typ);
1562 return bx;
1563}
1564
1565template<int dim>
1566[[nodiscard]]
1569BoxND<dim> convert (const BoxND<dim>& b, const IndexTypeND<dim>& typ) noexcept
1570{
1571 BoxND<dim> bx(b);
1572 bx.convert(typ);
1573 return bx;
1574}
1575
1582template<int dim>
1583[[nodiscard]]
1586BoxND<dim> enclosedCells (const BoxND<dim>& b, int dir) noexcept
1587{
1588 BoxND<dim> bx(b);
1589 bx.enclosedCells(dir);
1590 return bx;
1591}
1592
1593template<int dim>
1594[[nodiscard]]
1598{
1599 return enclosedCells(b, static_cast<int>(d));
1600}
1601
1606template<int dim>
1607[[nodiscard]]
1611{
1612 BoxND<dim> bx(b);
1613 bx.enclosedCells();
1614 return bx;
1615}
1616
1621template<int dim>
1622[[nodiscard]]
1625BoxND<dim> bdryLo (const BoxND<dim>& b, int dir, int len=1) noexcept
1626{
1627 IntVectND<dim> low(b.smallEnd());
1628 IntVectND<dim> hi(b.bigEnd());
1629 int sm = low[dir];
1630 low.setVal(dir,sm-len+1);
1631 hi.setVal(dir,sm);
1632 //
1633 // set dir'th bit to 1 = IndexType::NODE.
1634 //
1635 IndexTypeND<dim> typ(b.ixType());
1636 typ.set(dir);
1637 return BoxND<dim>(low,hi,typ);
1638}
1639
1644template<int dim>
1645[[nodiscard]]
1648BoxND<dim> bdryHi (const BoxND<dim>& b, int dir, int len=1) noexcept
1649{
1650 IntVectND<dim> low(b.smallEnd());
1651 IntVectND<dim> hi(b.bigEnd());
1652 auto const bitval = b.type()[dir];
1653 int bg = hi[dir] + 1 - bitval%2;
1654 low.setVal(dir,bg);
1655 hi.setVal(dir,bg+len-1);
1656 //
1657 // Set dir'th bit to 1 = IndexType::NODE.
1658 //
1659 IndexTypeND<dim> typ(b.ixType());
1660 typ.set(dir);
1661 return BoxND<dim>(low,hi,typ);
1662}
1663
1668template<int dim>
1669[[nodiscard]]
1672BoxND<dim> bdryNode (const BoxND<dim>& b, Orientation face, int len=1) noexcept
1673{
1674 int dir = face.coordDir();
1675 IntVectND<dim> low(b.smallEnd());
1676 IntVectND<dim> hi(b.bigEnd());
1677 if (face.isLow())
1678 {
1679 int sm = low[dir];
1680 low.setVal(dir,sm-len+1);
1681 hi.setVal(dir,sm);
1682 }
1683 else
1684 {
1685 int bitval = b.type()[dir];
1686 int bg = hi[dir] + 1 - bitval%2;
1687 low.setVal(dir,bg);
1688 hi.setVal(dir,bg+len-1);
1689 }
1690 //
1691 // Set dir'th bit to 1 = IndexType::NODE.
1692 //
1693 IndexTypeND<dim> typ(b.ixType());
1694 typ.set(dir);
1695 return BoxND<dim>(low,hi,typ);
1696}
1697
1710template<int dim>
1711[[nodiscard]]
1714BoxND<dim> adjCellLo (const BoxND<dim>& b, int dir, int len=1) noexcept
1715{
1716 BL_ASSERT(len > 0);
1717 IntVectND<dim> low(b.smallEnd());
1718 IntVectND<dim> hi(b.bigEnd());
1719 int sm = low[dir];
1720 low.setVal(dir,sm - len);
1721 hi.setVal(dir,sm - 1);
1722 //
1723 // Set dir'th bit to 0 = IndexType::CELL.
1724 //
1725 IndexTypeND<dim> typ(b.ixType());
1726 typ.unset(dir);
1727 return BoxND<dim>(low,hi,typ);
1728}
1729
1731template<int dim>
1732[[nodiscard]]
1735BoxND<dim> adjCellHi (const BoxND<dim>& b, int dir, int len=1) noexcept
1736{
1737 BL_ASSERT(len > 0);
1738 IntVectND<dim> low(b.smallEnd());
1739 IntVectND<dim> hi(b.bigEnd());
1740 int bitval = b.type()[dir];
1741 int bg = hi[dir] + 1 - bitval%2;
1742 low.setVal(dir,bg);
1743 hi.setVal(dir,bg + len - 1);
1744 //
1745 // Set dir'th bit to 0 = IndexType::CELL.
1746 //
1747 IndexTypeND<dim> typ(b.ixType());
1748 typ.unset(dir);
1749 return BoxND<dim>(low,hi,typ);
1750}
1751
1753template<int dim>
1754[[nodiscard]]
1757BoxND<dim> adjCell (const BoxND<dim>& b, Orientation face, int len=1) noexcept
1758{
1759 BL_ASSERT(len > 0);
1760 IntVectND<dim> low(b.smallEnd());
1761 IntVectND<dim> hi(b.bigEnd());
1762 int dir = face.coordDir();
1763 if (face.isLow())
1764 {
1765 int sm = low[dir];
1766 low.setVal(dir,sm - len);
1767 hi.setVal(dir,sm - 1);
1768 }
1769 else
1770 {
1771 int bitval = b.type()[dir];
1772 int bg = hi[dir] + 1 - bitval%2;
1773 low.setVal(dir,bg);
1774 hi.setVal(dir,bg + len - 1);
1775 }
1776 //
1777 // Set dir'th bit to 0 = IndexType::CELL.
1778 //
1779 IndexTypeND<dim> typ(b.ixType());
1780 typ.unset(dir);
1781 return BoxND<dim>(low,hi,typ);
1782}
1783
1789template<int dim>
1790[[nodiscard]]
1793BoxND<dim> minBox (const BoxND<dim>& b1, const BoxND<dim>& b2) noexcept
1794{
1795 BoxND<dim> result = b1;
1796 result.minBox(b2);
1797 return result;
1798}
1799
1801namespace detail {
1802 std::ostream& box_write (std::ostream& os, const int * smallend, const int * bigend,
1803 const int * type, int dim);
1804 std::istream& box_read (std::istream& is, int * smallend, int * bigend, int * type, int dim);
1805
1806 template<std::size_t...Ns, class T, class U>
1808 auto BoxSplit_imp (std::index_sequence<Ns...>,
1809 const T& lo, const T& hi, const U& typ) noexcept {
1810 return makeTuple(BoxND(get<Ns>(lo), get<Ns>(hi), get<Ns>(typ))...);
1811 }
1812}
1814
1816template<int dim>
1817std::ostream& operator<< (std::ostream& os, const BoxND<dim>& bx)
1818{
1819 IntVectND<dim> type = bx.type();
1820 return detail::box_write(os, bx.smallEnd().begin(), bx.bigEnd().begin(), type.begin(), dim);
1821}
1822
1824template<int dim>
1825std::istream& operator>> (std::istream& is, BoxND<dim>& bx) {
1826 IntVectND<dim> small;
1827 IntVectND<dim> big;
1828 IntVectND<dim> type;
1829 detail::box_read(is, small.begin(), big.begin(), type.begin(), dim);
1830 bx = BoxND<dim>{small, big, type};
1831 return is;
1832}
1833
1838template<int d, int...dims>
1841constexpr BoxND<detail::get_sum<d, dims...>()>
1842BoxCat (const BoxND<d>& bx, const BoxND<dims>&...boxes) noexcept {
1843 auto lo = IntVectCat(bx.smallEnd(), boxes.smallEnd()...);
1844 auto hi = IntVectCat(bx.bigEnd(), boxes.bigEnd()...);
1845 auto typ = IndexTypeCat(bx.ixType(), boxes.ixType()...);
1846 return BoxND<detail::get_sum<d, dims...>()>{lo, hi, typ};
1847}
1848
1853template<int d, int...dims>
1856constexpr GpuTuple<BoxND<d>, BoxND<dims>...>
1857BoxSplit (const BoxND<detail::get_sum<d, dims...>()>& bx) noexcept {
1858 auto lo = IntVectSplit<d, dims...>(bx.smallEnd());
1859 auto hi = IntVectSplit<d, dims...>(bx.bigEnd());
1860 auto typ = IndexTypeSplit<d, dims...>(bx.ixType());
1861 return detail::BoxSplit_imp(std::make_index_sequence<1 + sizeof...(dims)>(), lo, hi, typ);
1862}
1863
1868template<int new_dim, int old_dim>
1871constexpr BoxND<new_dim>
1872BoxShrink (const BoxND<old_dim>& bx) noexcept {
1873 return bx.template shrink<new_dim>();
1874}
1875
1881template<int new_dim, int old_dim>
1884constexpr BoxND<new_dim>
1885BoxExpand (const BoxND<old_dim>& bx) noexcept {
1886 return bx.template expand<new_dim>();
1887}
1888
1893template<int new_dim, int old_dim>
1896constexpr BoxND<new_dim>
1897BoxResize (const BoxND<old_dim>& bx) noexcept {
1898 return bx.template resize<new_dim>();
1899}
1900
1901template<int dim>
1902[[nodiscard]]
1906{
1907 return box.smallEnd();
1908}
1909
1910template<int dim>
1911[[nodiscard]]
1915{
1916 return box.bigEnd();
1917}
1918
1919template<int dim>
1920[[nodiscard]]
1924{
1925 return box.smallEnd();
1926}
1927
1928template<int dim>
1929[[nodiscard]]
1932IntVectND<dim> end_iv (BoxND<dim> const& box) noexcept
1933{
1934 return box.bigEnd() + 1;
1935}
1936
1937template<int dim>
1938[[nodiscard]]
1942{
1943 return box.bigEnd() - box.smallEnd() + 1;
1944}
1945
1946// Max of lower bound
1947template<int dim>
1948[[nodiscard]]
1951IntVectND<dim> max_lbound_iv (BoxND<dim> const& b1, BoxND<dim> const& b2) noexcept
1952{
1953 return max(b1.smallEnd(), b2.smallEnd());
1954}
1955
1956template<int dim>
1957[[nodiscard]]
1961{
1962 return max(b1.smallEnd(), lo);
1963}
1964
1965// Min of upper bound
1966template<int dim>
1967[[nodiscard]]
1970IntVectND<dim> min_ubound_iv (BoxND<dim> const& b1, BoxND<dim> const& b2) noexcept
1971{
1972 return min(b1.bigEnd(), b2.bigEnd());
1973}
1974
1975template<int dim>
1976[[nodiscard]]
1980{
1981 return min(b1.bigEnd(), hi);
1982}
1983
1984template<int dim, std::enable_if_t<( 1<=dim && dim<=3 ), int> = 0>
1985[[nodiscard]]
1986AMREX_GPU_HOST_DEVICE
1987AMREX_FORCE_INLINE
1988Dim3 lbound (BoxND<dim> const& box) noexcept
1989{
1990 return box.smallEnd().dim3();
1991}
1992
1993template<int dim, std::enable_if_t<( 1<=dim && dim<=3 ), int> = 0>
1994[[nodiscard]]
1995AMREX_GPU_HOST_DEVICE
1996AMREX_FORCE_INLINE
1997Dim3 ubound (BoxND<dim> const& box) noexcept
1998{
1999 return box.bigEnd().dim3();
2000}
2001
2002template<int dim, std::enable_if_t<( 1<=dim && dim<=3 ), int> = 0>
2003[[nodiscard]]
2004AMREX_GPU_HOST_DEVICE
2005AMREX_FORCE_INLINE
2006Dim3 begin (BoxND<dim> const& box) noexcept
2007{
2008 return box.smallEnd().dim3();
2009}
2010
2011template<int dim, std::enable_if_t<( 1<=dim && dim<=3 ), int> = 0>
2012[[nodiscard]]
2013AMREX_GPU_HOST_DEVICE
2014AMREX_FORCE_INLINE
2015Dim3 end (BoxND<dim> const& box) noexcept
2016{
2017 return (box.bigEnd() + 1).dim3(1);
2018}
2019
2020template<int dim, std::enable_if_t<( 1<=dim && dim<=3 ), int> = 0>
2021[[nodiscard]]
2022AMREX_GPU_HOST_DEVICE
2023AMREX_FORCE_INLINE
2024Dim3 length (BoxND<dim> const& box) noexcept
2025{
2026 return (box.bigEnd() - box.smallEnd() + 1).dim3(1);
2027}
2028
2029// Max of lower bound
2030template<int dim, std::enable_if_t<( 1<=dim && dim<=3 ), int> = 0>
2031[[nodiscard]]
2032AMREX_GPU_HOST_DEVICE
2033AMREX_FORCE_INLINE
2034Dim3 max_lbound (BoxND<dim> const& b1, BoxND<dim> const& b2) noexcept
2035{
2036 return max(b1.smallEnd(), b2.smallEnd()).dim3();
2037}
2038
2039template<int dim, std::enable_if_t<( 1<=dim && dim<=3 ), int> = 0>
2040[[nodiscard]]
2041AMREX_GPU_HOST_DEVICE
2042AMREX_FORCE_INLINE
2043Dim3 max_lbound (BoxND<dim> const& b1, Dim3 const& lo) noexcept
2044{
2045 return max(b1.smallEnd(), IntVectND<dim>(lo)).dim3();
2046}
2047
2048// Min of upper bound
2049template<int dim, std::enable_if_t<( 1<=dim && dim<=3 ), int> = 0>
2050[[nodiscard]]
2051AMREX_GPU_HOST_DEVICE
2052AMREX_FORCE_INLINE
2053Dim3 min_ubound (BoxND<dim> const& b1, BoxND<dim> const& b2) noexcept
2054{
2055 return min(b1.bigEnd(), b2.bigEnd()).dim3();
2056}
2057
2058template<int dim, std::enable_if_t<( 1<=dim && dim<=3 ), int> = 0>
2059[[nodiscard]]
2060AMREX_GPU_HOST_DEVICE
2061AMREX_FORCE_INLINE
2062Dim3 min_ubound (BoxND<dim> const& b1, Dim3 const& hi) noexcept
2063{
2064 return min(b1.bigEnd(), IntVectND<dim>(hi)).dim3();
2065}
2066
2067// Return a BoxND that covers all the argument Boxes in index
2068// space. The types are ignored. Thus, the arguments can have
2069// different index types, and the returned BoxND's index type has no
2070// meaning.
2071template<int dim>
2072[[nodiscard]]
2075{
2076 return b1;
2077}
2078
2079template<int dim>
2080[[nodiscard]]
2082BoxND<dim> getIndexBounds (BoxND<dim> const& b1, BoxND<dim> const& b2) noexcept
2083{
2084 BoxND<dim> b = b1;
2085 b.setType(b2.ixType());
2086 b.minBox(b2);
2087 return b;
2088}
2089
2090template <class T, class ... Ts>
2091[[nodiscard]]
2093auto getIndexBounds (T const& b1, T const& b2, Ts const& ... b3) noexcept
2094{
2095 return getIndexBounds(getIndexBounds(b1,b2),b3...);
2096}
2097
2098
2099template<int dim>
2100[[nodiscard]]
2103IntVectND<dim> getCell (BoxND<dim> const* boxes, int nboxes, Long icell) noexcept
2104{
2105 int ibox;
2106 Long ncells_subtotal = 0;
2107 Long offset = 0;
2108 for (ibox = 0; ibox < nboxes; ++ibox) {
2109 const Long n = boxes[ibox].numPts();
2110 ncells_subtotal += n;
2111 if (icell < ncells_subtotal) {
2112 offset = icell - (ncells_subtotal - n);
2113 break;
2114 }
2115 }
2116 return boxes[ibox].atOffset(offset);
2117}
2118
2119template<int dim>
2120[[nodiscard]]
2123BoxND<dim> makeSlab (BoxND<dim> const& b, int direction, int slab_index) noexcept
2124{
2125 BoxND<dim> r = b;
2126 r.makeSlab(direction,slab_index);
2127 return r;
2128}
2129
2130template<int dim=AMREX_SPACEDIM, std::enable_if_t<( 1<=dim && dim<=3 ), int> = 0>
2131[[nodiscard]]
2132AMREX_GPU_HOST_DEVICE
2133AMREX_FORCE_INLINE
2134BoxND<dim> makeSingleCellBox (int i, int j, int k, IndexTypeND<dim> typ = IndexTypeND<dim>::TheCellType())
2135{
2136 Dim3 p3d{i, j, k};
2137 IntVectND<dim> vect{p3d};
2138 return BoxND<dim>{vect, vect, typ};
2139}
2140
2141template<int dim>
2142[[nodiscard]]
2146{
2147 return BoxND<dim>{vect, vect, typ};
2148}
2149
2150template<int dim>
2152{
2153 std::uint64_t npts;
2156
2158 : npts(box.numPts()),
2159 lo (box.smallEnd())
2160 {
2161 std::uint64_t mult = 1;
2162 for (int i=0; i<dim-1; ++i) {
2163 mult *= box.length(i);
2164 fdm[i] = Math::FastDivmodU64(mult);
2165 }
2166 }
2167
2169 IntVectND<dim> intVect (std::uint64_t icell) const
2170 {
2171 IntVectND<dim> retval = lo;
2172
2173 for (int i=dim-1; i>0; --i) {
2174 std::uint64_t quotient, remainder;
2175 fdm[i-1](quotient, remainder, icell);
2176 retval[i] += quotient;
2177 icell = remainder;
2178 }
2179
2180 retval[0] += icell;
2181
2182 return retval;
2183 }
2184
2185 template <int N=dim, std::enable_if_t<( 1<=N && N<=3 ), int> = 0>
2186 [[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
2187 Dim3 operator() (std::uint64_t icell) const
2188 {
2189 return intVect(icell).dim3();
2190 }
2191
2192 [[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
2193 std::uint64_t numPts () const { return npts; }
2194};
2195
2196template<>
2198{
2199 std::uint64_t npts;
2200
2201 int lo;
2202
2204 : npts(box.numPts()),
2205 lo(box.smallEnd(0))
2206 {}
2207
2209 IntVectND<1> intVect (std::uint64_t icell) const
2210 {
2211 return IntVectND<1>{int(icell)+lo};
2212 }
2213
2215 Dim3 operator() (std::uint64_t icell) const
2216 {
2217 return {int(icell)+lo, 0, 0};
2218 }
2219
2221 std::uint64_t numPts () const { return npts; }
2222};
2223
2225
2226}
2227
2228#endif /*AMREX_BOX_H*/
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#define AMREX_FORCE_INLINE
Definition AMReX_Extension.H:119
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
int idir
Definition AMReX_HypreMLABecLap.cpp:1093
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1089
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Definition AMReX_Box.H:1257
virtual Box doit(const Box &fine) const =0
virtual ~BoxConverter()=default
virtual BoxConverter * clone() const =0
iterates through the IntVects of a Box
Definition AMReX_BoxIterator.H:51
A Rectangular Domain on an Integer Lattice.
Definition AMReX_Box.H:49
IntVectND< dim > smallend
Definition AMReX_Box.H:902
__host__ __device__ constexpr BoxND(const IntVectND< dim > &small, const IntVectND< dim > &big) noexcept
Construct cell-centered type BoxND.
Definition AMReX_Box.H:67
__host__ __device__ BoxND & grow(int i) noexcept
Definition AMReX_Box.H:641
__host__ __device__ const IntVectND< dim > & bigEnd() const &noexcept
Return the inclusive upper bound of the box.
Definition AMReX_Box.H:123
IndexTypeND< dim > btype
Definition AMReX_Box.H:904
__host__ __device__ BoxND(const IntVectND< dim > &small, const IntVectND< dim > &big, IndexTypeND< dim > t) noexcept
Construct dimension specific Boxes.
Definition AMReX_Box.H:94
__host__ static __device__ constexpr int indims() noexcept
Definition AMReX_Box.H:854
__host__ __device__ BoxND< new_dim > resize() const noexcept
Return a new BoxND of size new_dim by either shrinking or expanding this BoxND.
Definition AMReX_Box.H:893
__host__ __device__ BoxND & makeSlab(int direction, int slab_index) noexcept
Definition AMReX_Box.H:826
__host__ __device__ BoxND & minBox(const BoxND &b) noexcept
Modify BoxND to that of the minimum BoxND containing both the original BoxND and the argument....
Definition AMReX_Box.H:597
__host__ __device__ BoxND< new_dim > expand() const noexcept
Return a new BoxND of size new_dim and assigns all values of this BoxND to it and (small=0,...
Definition AMReX_Box.H:879
__host__ __device__ bool isSquare() const noexcept
Definition AMReX_Box.H:1141
__host__ __device__ BoxND(const IntVectND< dim > &small, const int *vec_len) noexcept
Construct BoxND with specified lengths.
Definition AMReX_Box.H:74
__host__ __device__ bool coarsenable(int refrat, int min_width=1) const noexcept
Return whether this Box is coarsenable.
Definition AMReX_Box.H:794
__host__ __device__ Long numPts() const noexcept
Return the number of points contained in the BoxND.
Definition AMReX_Box.H:356
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:154
__host__ __device__ BoxND & convert(IndexTypeND< dim > typ) noexcept
Convert the BoxND from the current type into the argument type. This may change the BoxND coordinates...
Definition AMReX_Box.H:974
__host__ __device__ BoxND & shift(int dir, int nzones) noexcept
Shift this BoxND nzones indexing positions in coordinate direction dir.
Definition AMReX_Box.H:509
BoxIteratorND< dim > iterator() const noexcept
Return a BoxIteratorND that can be used to loop over the IntVects contained by the BoxND.
Definition AMReX_Box.H:844
__host__ static __device__ constexpr std::size_t ndims() noexcept
Definition AMReX_Box.H:849
IntVectND< dim > bigend
Definition AMReX_Box.H:903
__host__ __device__ IntVectND< dim > atOffset(Long offset) const noexcept
Given the offset, compute IntVectND<dim>
Definition AMReX_Box.H:1071
__host__ __device__ BoxND & growLo(int idir, int n_cell=1) noexcept
Grow the BoxND on the low end by n_cell cells in direction idir. NOTE: n_cell negative shrinks the Bo...
Definition AMReX_Box.H:662
__host__ __device__ IndexTypeND< dim > ixType() const noexcept
Return the indexing type.
Definition AMReX_Box.H:135
__host__ __device__ BoxND & growHi(int idir, int n_cell=1) noexcept
Grow the BoxND on the high end by n_cell cells in direction idir. NOTE: n_cell negative shrinks the B...
Definition AMReX_Box.H:673
__host__ __device__ IntVectND< dim > type() const noexcept
Return the indexing type.
Definition AMReX_Box.H:139
__host__ __device__ BoxND & coarsen(int ref_ratio) noexcept
Coarsen BoxND by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo/rati...
Definition AMReX_Box.H:722
__host__ __device__ BoxND & refine(int ref_ratio) noexcept
Refine BoxND by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo*ratio...
Definition AMReX_Box.H:698
__host__ __device__ BoxND & setType(const IndexTypeND< dim > &t) noexcept
Set indexing type.
Definition AMReX_Box.H:505
__host__ __device__ constexpr BoxND() noexcept
Definition AMReX_Box.H:60
__host__ __device__ BoxND & enclosedCells() noexcept
Convert to CELL type in all directions.
Definition AMReX_Box.H:1040
__host__ __device__ void next(IntVectND< dim > &) const noexcept
Step through the rectangle. It is a runtime error to give a point not inside rectangle....
Definition AMReX_Box.H:1121
__host__ __device__ void normalize() noexcept
Definition AMReX_Box.H:816
__host__ __device__ BoxND< new_dim > shrink() const noexcept
Return a new BoxND of dimension new_dim and assigns the first new_dim dimension of this BoxND to it.
Definition AMReX_Box.H:864
__host__ __device__ BoxND & coarsen(const IntVectND< dim > &ref_ratio) noexcept
Coarsen BoxND by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo/rati...
Definition AMReX_Box.H:928
__host__ static __device__ BoxND TheUnitBox() noexcept
This static member function returns a constant reference to an object of type BoxND representing the ...
Definition AMReX_Box.H:751
__host__ __device__ bool coarsenable(const IntVectND< dim > &refrat, int min_width=1) const noexcept
Return whether this Box is coarsenable.
Definition AMReX_Box.H:810
friend class BoxCommHelper
Definition AMReX_Box.H:51
__host__ __device__ BoxND & surroundingNodes() noexcept
Convert to NODE type in all directions.
Definition AMReX_Box.H:1008
__host__ __device__ BoxND(const IntVectND< dim > &small, const IntVectND< dim > &big, const IntVectND< dim > &typ) noexcept
Construct BoxND with given type. small and big are expected to be consistent with given type.
Definition AMReX_Box.H:84
__host__ __device__ const IntVectND< dim > & smallEnd() const &noexcept
Return the inclusive lower bound of the box.
Definition AMReX_Box.H:111
Cell-Based or Node-Based Indices.
Definition AMReX_IndexType.H:36
__host__ __device__ void set(int dir) noexcept
Set IndexTypeND to be NODE based in direction dir.
Definition AMReX_IndexType.H:72
__host__ __device__ void unset(int dir) noexcept
Set IndexTypeND to be CELL based in direction dir.
Definition AMReX_IndexType.H:75
An Integer Vector in dim-Dimensional Space.
Definition AMReX_IntVect.H:57
__host__ __device__ IntVectND & setVal(int i, int val) noexcept
Set i'th coordinate of IntVectND to val.
Definition AMReX_IntVect.H:282
__host__ __device__ constexpr int * begin() noexcept
Returns a pointer to the first element of the IntVectND.
Definition AMReX_IntVect.H:266
__host__ __device__ IntVectND< dim > & coarsen(const IntVectND< dim > &p) noexcept
Modify IntVectND<dim> by component-wise integer projection.
Definition AMReX_IntVect.H:854
Encapsulation of the Orientation of the Faces of a Box.
Definition AMReX_Orientation.H:29
__host__ __device__ bool isLow() const noexcept
Returns true if Orientation is low.
Definition AMReX_Orientation.H:89
__host__ __device__ int coordDir() const noexcept
Returns the coordinate direction.
Definition AMReX_Orientation.H:83
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
amrex_long Long
Definition AMReX_INT.H:30
__host__ __device__ BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) coarsening ratio.
Definition AMReX_Box.H:1409
__host__ __device__ BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition AMReX_Box.H:1280
__host__ __device__ BoxND< dim > refine(const BoxND< dim > &b, int ref_ratio) noexcept
Definition AMReX_Box.H:1459
int MPI_Datatype
Definition AMReX_ccse-mpi.H:53
Definition AMReX_Amr.cpp:49
__host__ __device__ BoxND< dim > adjCellHi(const BoxND< dim > &b, int dir, int len=1) noexcept
Similar to adjCellLo but builds an adjacent BoxND on the high end.
Definition AMReX_Box.H:1735
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1558
__host__ __device__ BoxND< dim > makeSlab(BoxND< dim > const &b, int direction, int slab_index) noexcept
Definition AMReX_Box.H:2123
__host__ __device__ BoxND< dim > adjCellLo(const BoxND< dim > &b, int dir, int len=1) noexcept
Return the cell centered BoxND of length len adjacent to b on the low end along the coordinate direct...
Definition AMReX_Box.H:1714
__host__ __device__ BoxND< dim > surroundingNodes(const BoxND< dim > &b, int dir) noexcept
Return a BoxND with NODE based coordinates in direction dir that encloses BoxND b....
Definition AMReX_Box.H:1522
__host__ __device__ Dim3 length(Array4< T > const &a) noexcept
Definition AMReX_Array4.H:326
__host__ __device__ IntVectND< dim > end_iv(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:1932
__host__ __device__ IntVectND< dim > lbound_iv(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:1905
__host__ __device__ BoxND< dim > bdryLo(const BoxND< dim > &b, int dir, int len=1) noexcept
Return the edge-centered BoxND (in direction dir) defining the low side of BoxND b.
Definition AMReX_Box.H:1625
__host__ __device__ constexpr GpuTuple< IntVectND< d >, IntVectND< dims >... > IntVectSplit(const IntVectND< detail::get_sum< d, dims... >()> &v) noexcept
Returns a tuple of IntVectND obtained by splitting the input IntVectND according to the dimensions sp...
Definition AMReX_IntVect.H:1187
__host__ __device__ IntVectND< dim > min_ubound_iv(BoxND< dim > const &b1, BoxND< dim > const &b2) noexcept
Definition AMReX_Box.H:1970
__host__ __device__ BoxND< dim > minBox(const BoxND< dim > &b1, const BoxND< dim > &b2) noexcept
Modify BoxND to that of the minimum BoxND containing both the original BoxND and the argument....
Definition AMReX_Box.H:1793
__host__ __device__ BoxND< dim > makeSingleCellBox(int i, int j, int k, IndexTypeND< dim > typ=IndexTypeND< dim >::TheCellType())
Definition AMReX_Box.H:2134
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:27
__host__ __device__ constexpr BoxND< new_dim > BoxResize(const BoxND< old_dim > &bx) noexcept
Return a new BoxND of size new_dim by either shrinking or expanding this BoxND.
Definition AMReX_Box.H:1897
__host__ __device__ constexpr const T & min(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:21
__host__ __device__ BoxND< dim > adjCell(const BoxND< dim > &b, Orientation face, int len=1) noexcept
Similar to adjCellLo and adjCellHi; operates on given face.
Definition AMReX_Box.H:1757
__host__ __device__ IntVectND< dim > begin_iv(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:1923
__host__ __device__ BoxND< dim > shift(const BoxND< dim > &b, int dir, int nzones) noexcept
Return a BoxND with indices shifted by nzones in dir direction.
Definition AMReX_Box.H:1495
__host__ __device__ constexpr BoxND< new_dim > BoxShrink(const BoxND< old_dim > &bx) noexcept
Return a new BoxND of dimension new_dim and assigns the first new_dim dimension of this BoxND to it.
Definition AMReX_Box.H:1872
__host__ __device__ BoxND< dim > enclosedCells(const BoxND< dim > &b, int dir) noexcept
Return a BoxND with CELL based coordinates in direction dir that is enclosed by b....
Definition AMReX_Box.H:1586
__host__ __device__ constexpr IndexTypeND< detail::get_sum< d, dims... >()> IndexTypeCat(const IndexTypeND< d > &v, const IndexTypeND< dims > &...vects) noexcept
Returns a IndexTypeND obtained by concatenating the input IndexTypeNDs. The dimension of the return v...
Definition AMReX_IndexType.H:297
__host__ __device__ IntVectND< dim > length_iv(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:1941
Direction
Definition AMReX_Orientation.H:14
__host__ __device__ BoxND< dim > bdryNode(const BoxND< dim > &b, Orientation face, int len=1) noexcept
Similar to bdryLo and bdryHi except that it operates on the given face of box b.
Definition AMReX_Box.H:1672
__host__ __device__ constexpr const T & max(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:35
__host__ __device__ BoxND< dim > growLo(const BoxND< dim > &b, int idir, int n_cell) noexcept
Definition AMReX_Box.H:1354
__host__ __device__ IntVectND< dim > ubound_iv(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:1914
__host__ __device__ constexpr BoxND< new_dim > BoxExpand(const BoxND< old_dim > &bx) noexcept
Return a new BoxND of size new_dim and assigns all values of this BoxND to it and (small=0,...
Definition AMReX_Box.H:1885
__host__ __device__ IntVectND< dim > getCell(BoxND< dim > const *boxes, int nboxes, Long icell) noexcept
Definition AMReX_Box.H:2103
void AllGatherBoxes(Vector< Box > &bxs, int n_extra_reserve)
Definition AMReX_Box.cpp:124
__host__ __device__ constexpr IntVectND< detail::get_sum< d, dims... >()> IntVectCat(const IntVectND< d > &v, const IntVectND< dims > &...vects) noexcept
Returns a IntVectND obtained by concatenating the input IntVectNDs. The dimension of the return value...
Definition AMReX_IntVect.H:1171
__host__ __device__ BoxND< dim > bdryHi(const BoxND< dim > &b, int dir, int len=1) noexcept
Return the edge-centered BoxND (in direction dir) defining the high side of BoxND b.
Definition AMReX_Box.H:1648
__host__ __device__ BoxND< dim > growHi(const BoxND< dim > &b, int idir, int n_cell) noexcept
Definition AMReX_Box.H:1374
BoxND< dim > getIndexBounds(BoxND< dim > const &b1) noexcept
Definition AMReX_Box.H:2074
__host__ __device__ constexpr GpuTuple< BoxND< d >, BoxND< dims >... > BoxSplit(const BoxND< detail::get_sum< d, dims... >()> &bx) noexcept
Return a tuple of BoxNDs obtained by splitting the input BoxND according to the dimensions specified ...
Definition AMReX_Box.H:1857
__host__ __device__ constexpr GpuTuple< IndexTypeND< d >, IndexTypeND< dims >... > IndexTypeSplit(const IndexTypeND< detail::get_sum< d, dims... >()> &v) noexcept
Returns a tuple of IndexTypeND obtained by splitting the input IndexTypeND according to the dimension...
Definition AMReX_IndexType.H:317
__host__ __device__ IntVectND< dim > max_lbound_iv(BoxND< dim > const &b1, BoxND< dim > const &b2) noexcept
Definition AMReX_Box.H:1951
BoxIndexerND(BoxND< 1 > const &box)
Definition AMReX_Box.H:2203
std::uint64_t npts
Definition AMReX_Box.H:2199
__host__ __device__ std::uint64_t numPts() const
Definition AMReX_Box.H:2221
__host__ __device__ IntVectND< 1 > intVect(std::uint64_t icell) const
Definition AMReX_Box.H:2209
int lo
Definition AMReX_Box.H:2201
Definition AMReX_Box.H:2152
__host__ __device__ IntVectND< dim > intVect(std::uint64_t icell) const
Definition AMReX_Box.H:2169
BoxIndexerND(BoxND< dim > const &box)
Definition AMReX_Box.H:2157
IntVectND< dim > lo
Definition AMReX_Box.H:2155
std::uint64_t npts
Definition AMReX_Box.H:2153
CellIndex
The cell index type: one of CELL or NODE.
Definition AMReX_IndexType.H:21
@ CELL
Definition AMReX_IndexType.H:21
Definition AMReX_Dim3.H:12
Fixed-size array that can be used on GPU.
Definition AMReX_Array.H:40
Definition AMReX_Math.H:438