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>
101 requires ( 1<=dim && dim<=3 )
103 explicit BoxND (Array4<T> const& a) noexcept
104 : smallend(IntVectShrink<dim>(a.begin)),
105 bigend(IntVectShrink<dim>(a.end) - 1)
106 {}
107
108 // dtor, copy-ctor, copy-op=, move-ctor, and move-op= are compiler generated.
109
111 [[nodiscard]] AMREX_GPU_HOST_DEVICE
112 const IntVectND<dim>& smallEnd () const& noexcept { return smallend; }
113
115 [[nodiscard]] const IntVectND<dim>& smallEnd () && = delete;
117
119 [[nodiscard]] AMREX_GPU_HOST_DEVICE
120 int smallEnd (int dir) const& noexcept { return smallend[dir]; }
121
123 [[nodiscard]] AMREX_GPU_HOST_DEVICE
124 const IntVectND<dim>& bigEnd () const& noexcept { return bigend; }
125
127 [[nodiscard]] const IntVectND<dim>& bigEnd () && = delete;
129
131 [[nodiscard]] AMREX_GPU_HOST_DEVICE
132 int bigEnd (int dir) const noexcept { return bigend[dir]; }
133
135 [[nodiscard]] AMREX_GPU_HOST_DEVICE
136 IndexTypeND<dim> ixType () const noexcept { return btype; }
137
139 [[nodiscard]] AMREX_GPU_HOST_DEVICE
140 IntVectND<dim> type () const noexcept { return btype.ixType(); }
141
143 [[nodiscard]] AMREX_GPU_HOST_DEVICE
144 IndexType::CellIndex type (int dir) const noexcept { return btype.ixType(dir); }
145
147 [[nodiscard]] AMREX_GPU_HOST_DEVICE
148 IntVectND<dim> size () const noexcept
149 {
150 return bigend - smallend + 1;
151 }
152
154 [[nodiscard]] AMREX_GPU_HOST_DEVICE
155 IntVectND<dim> length () const noexcept
156 {
157 return bigend - smallend + 1;
158 }
159
161 [[nodiscard]] AMREX_GPU_HOST_DEVICE
162 int length (int dir) const noexcept { return bigend[dir] - smallend[dir] + 1; }
163
164 [[nodiscard]] AMREX_GPU_HOST_DEVICE
165 GpuArray<int,3> length3d () const noexcept
166 requires (1 <= dim && dim <= 3)
167 {
168 Dim3 len3d = length().dim3(1);
169 return {{len3d.x, len3d.y, len3d.z}};
170 }
171
172 [[nodiscard]] AMREX_GPU_HOST_DEVICE
173 GpuArray<int,3> loVect3d () const noexcept
174 requires (1 <= dim && dim <= 3)
175 {
176 Dim3 lo3d = smallend.dim3(0);
177 return {{lo3d.x, lo3d.y, lo3d.z}};
178 }
179
180 [[nodiscard]] AMREX_GPU_HOST_DEVICE
181 GpuArray<int,3> hiVect3d () const noexcept
182 requires (1 <= dim && dim <= 3)
183 {
184 Dim3 hi3d = bigend.dim3(0);
185 return {{hi3d.x, hi3d.y, hi3d.z}};
186 }
187
189 [[nodiscard]] AMREX_GPU_HOST_DEVICE
190 const int* loVect () const& noexcept { return smallend.getVect(); }
192 const int* loVect () && = delete;
194 [[nodiscard]] AMREX_GPU_HOST_DEVICE
195 const int* hiVect () const& noexcept { return bigend.getVect(); }
197 const int* hiVect () && = delete;
198
200 [[nodiscard]] AMREX_GPU_HOST_DEVICE
201 int operator[] (Orientation face) const noexcept {
202 const int dir = face.coordDir();
203 return face.isLow() ? smallend[dir] : bigend[dir];
204 }
205
207 [[nodiscard]] AMREX_GPU_HOST_DEVICE
208 bool isEmpty () const noexcept { return !ok(); }
209
211 [[nodiscard]] AMREX_GPU_HOST_DEVICE
212 bool ok () const noexcept { return bigend.allGE(smallend) && btype.ok(); }
213
215 [[nodiscard]] AMREX_GPU_HOST_DEVICE
216 bool contains (const IntVectND<dim>& p) const noexcept {
217 return p.allGE(smallend) && p.allLE(bigend);
218 }
219
221 [[nodiscard]] AMREX_GPU_HOST_DEVICE
222 bool contains (const Dim3& p) const noexcept
223 requires (1 <= dim && dim <= 3)
224 {
225 IntVectND<dim> piv{p};
226 return contains(piv);
227 }
228
230 [[nodiscard]] AMREX_GPU_HOST_DEVICE
231 bool contains (int i, int j, int k) const noexcept
232 requires (1 <= dim && dim <= 3)
233 {
234 Dim3 p3d{.x = i, .y = j, .z = k};
235 return contains(p3d);
236 }
237
241 [[nodiscard]] AMREX_GPU_HOST_DEVICE
242 bool contains (const BoxND& b) const noexcept
243 {
245 return b.smallend.allGE(smallend) && b.bigend.allLE(bigend);
246 }
247
249 [[nodiscard]] AMREX_GPU_HOST_DEVICE
250 bool strictly_contains (const IntVectND<dim>& p) const noexcept {
251 return p.allGT(smallend) && p.allLT(bigend);
252 }
253
258 [[nodiscard]] AMREX_GPU_HOST_DEVICE
259 bool strictly_contains (const BoxND& b) const noexcept
260 {
262 return b.smallend.allGT(smallend) && b.bigend.allLT(bigend);
263 }
264
266 [[nodiscard]] AMREX_GPU_HOST_DEVICE
267 bool strictly_contains (const Dim3& p) const noexcept
268 requires (1 <= dim && dim <= 3)
269 {
270 IntVectND<dim> piv{p};
271 return strictly_contains(piv);
272 }
273
275 [[nodiscard]] AMREX_GPU_HOST_DEVICE
276 bool strictly_contains (int i, int j, int k) const noexcept
277 requires (1 <= dim && dim <= 3)
278 {
279 Dim3 p3d{.x = i, .y = j, .z = k};
280 return strictly_contains(p3d);
281 }
282
287 [[nodiscard]] AMREX_GPU_HOST_DEVICE
288 bool intersects (const BoxND& b) const noexcept { BoxND isect(*this); isect &= b; return isect.ok(); }
289
294 [[nodiscard]] AMREX_GPU_HOST_DEVICE
295 bool sameSize (const BoxND& b) const noexcept {
297 return length() == b.length();
298 }
299
301 [[nodiscard]] AMREX_GPU_HOST_DEVICE
302 bool sameType (const BoxND &b) const noexcept { return btype == b.btype; }
303
305 [[nodiscard]] AMREX_GPU_HOST_DEVICE
306 bool operator== (const BoxND& b) const noexcept { return smallend == b.smallend && bigend == b.bigend && b.btype == btype; }
307
309 [[nodiscard]] AMREX_GPU_HOST_DEVICE
310 bool operator!= (const BoxND& b) const noexcept { return !operator==(b); }
311
312 [[nodiscard]] AMREX_GPU_HOST_DEVICE
313 bool operator< (const BoxND& 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 BoxND& rhs) const noexcept {
322 return !(rhs < *this);
323 }
324 [[nodiscard]] AMREX_GPU_HOST_DEVICE
325 bool operator> (const BoxND& rhs) const noexcept {
326 return rhs < *this;
327 }
328 [[nodiscard]] AMREX_GPU_HOST_DEVICE
329 bool operator>= (const BoxND& rhs) const noexcept {
330 return !(*this < rhs);
331 }
332
334 [[nodiscard]] AMREX_GPU_HOST_DEVICE
335 bool cellCentered () const noexcept { return !btype.any(); }
336
339 void checkOverflow () const {
340 if (ok()) {
341 for (int i = 0; i < dim; ++i) {
342 auto lo = static_cast<Long>(smallend[i]);
343 auto hi = static_cast<Long>(bigend[i]);
344 Long len = hi - lo + 1;
345 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(len>=0 && len<std::numeric_limits<int>::max(),
346 "Overflow when computing length of box");
347 }
348 auto num_pts = static_cast<Long>(length(0));
349 for (int i = 1; i < dim; ++i) {
350 auto len = static_cast<Long>(length(i));
351 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(num_pts == 0 || len == 0 ||
352 num_pts <= std::numeric_limits<Long>::max() / len,
353 "Overflow when computing numPts of box");
354 num_pts *= len;
355 }
356 }
357 }
359
363 [[nodiscard]] AMREX_GPU_HOST_DEVICE
364 Long numPts () const noexcept {
365#if defined(AMREX_DEBUG) || defined(AMREX_USE_ASSERTION)
366 AMREX_IF_ON_HOST((checkOverflow();))
367#endif
368 if (ok()) {
369 auto num_pts = static_cast<Long>(length(0));
370 for (int i = 1; i < dim; ++i) {
371 num_pts *= static_cast<Long>(length(i));
372 }
373 return num_pts;
374 } else {
375 return Long(0);
376 }
377 }
378
383 [[nodiscard]] AMREX_GPU_HOST_DEVICE
384 double d_numPts () const noexcept {
385 if (ok()) {
386 auto num_pts = static_cast<double>(length(0));
387 for (int i = 1; i < dim; ++i) {
388 num_pts *= static_cast<double>(length(i));
389 }
390 return num_pts;
391 } else {
392 return 0.0;
393 }
394 }
395
401 [[nodiscard]] AMREX_GPU_HOST_DEVICE
402 Long volume () const noexcept {
403 if (ok()) {
404 auto num_pts = static_cast<Long>(length(0)-btype[0]);
405 for (int i = 1; i < dim; ++i) {
406 num_pts *= static_cast<Long>(length(i)-btype[i]);
407 }
408 return num_pts;
409 } else {
410 return Long(0);
411 }
412 }
413
418 [[nodiscard]] AMREX_GPU_HOST_DEVICE
419 int longside (int& dir) const noexcept {
420 int maxlen = length(0);
421 dir = 0;
422 for (int i = 1; i < dim; i++)
423 {
424 if (length(i) > maxlen)
425 {
426 maxlen = length(i);
427 dir = i;
428 }
429 }
430 return maxlen;
431 }
432
434 [[nodiscard]] AMREX_GPU_HOST_DEVICE
435 int longside () const noexcept {
436 int ignore = 0;
437 return longside(ignore);
438 }
439
444 [[nodiscard]] AMREX_GPU_HOST_DEVICE
445 int shortside (int& dir) const noexcept {
446 int minlen = length(0);
447 dir = 0;
448 for (int i = 1; i < dim; i++)
449 {
450 if (length(i) < minlen)
451 {
452 minlen = length(i);
453 dir = i;
454 }
455 }
456 return minlen;
457 }
458
460 [[nodiscard]] AMREX_GPU_HOST_DEVICE
461 int shortside () const noexcept {
462 int ignore = 0;
463 return shortside(ignore);
464 }
465
471 [[nodiscard]] AMREX_GPU_HOST_DEVICE
472 Long index (const IntVectND<dim>& v) const noexcept;
473
475 [[nodiscard]] AMREX_GPU_HOST_DEVICE
477
481 template <int N=dim>
482 requires (1 <= N && N <= 3)
483 [[nodiscard]] AMREX_GPU_HOST_DEVICE
485
487 BoxND& setSmall (const IntVectND<dim>& sm) noexcept { smallend = sm; return *this; }
488
491 BoxND& setSmall (int dir, int sm_index) noexcept { smallend.setVal(dir,sm_index); return *this; }
492
495 BoxND& setBig (const IntVectND<dim>& bg) noexcept { bigend = bg; return *this; }
496
499 BoxND& setBig (int dir, int bg_index) noexcept { bigend.setVal(dir,bg_index); return *this; }
500
507 BoxND& setRange (int dir,
508 int sm_index,
509 int n_cells = 1) noexcept;
510
513 BoxND& setType (const IndexTypeND<dim>& t) noexcept { btype = t; return *this; }
514
517 BoxND& shift (int dir, int nzones) noexcept { smallend.shift(dir,nzones); bigend.shift(dir,nzones); return *this; }
518
521 BoxND& shift (const IntVectND<dim>& iv) noexcept { smallend.shift(iv); bigend.shift(iv); return *this; }
522
533 BoxND& shiftHalf (int dir, int num_halfs) noexcept;
534
537 BoxND& shiftHalf (const IntVectND<dim>& iv) noexcept;
538
548
557 BoxND& convert (const IntVectND<dim>& typ) noexcept;
558
562
565 BoxND& surroundingNodes (int dir) noexcept;
566
568 BoxND& surroundingNodes (Direction d) noexcept { return surroundingNodes(static_cast<int>(d)); }
569
572 BoxND& enclosedCells () noexcept;
573
576 BoxND& enclosedCells (int dir) noexcept;
577
580 BoxND& enclosedCells (Direction d) noexcept { return enclosedCells(static_cast<int>(d)); }
581
587 BoxND operator& (const BoxND& rhs) const noexcept { BoxND lhs(*this); lhs &= rhs; return lhs; }
588
591 BoxND& operator&= (const BoxND& rhs) noexcept
592 {
593 BL_ASSERT(sameType(rhs));
594 smallend.max(rhs.smallend);
595 bigend.min(rhs.bigend);
596 return *this;
597 }
598
605 BoxND& minBox (const BoxND& b) noexcept {
606 // BoxArray may call this with not ok boxes. BL_ASSERT(b.ok() && ok());
608 smallend.min(b.smallend);
609 bigend.max(b.bigend);
610 return *this;
611 }
612
615 BoxND& operator+= (const IntVectND<dim>& v) noexcept { smallend += v; bigend += v; return *this; }
616
619 BoxND operator+ (const IntVectND<dim>& v) const noexcept { BoxND r(*this); r += v; return r; }
620
623 BoxND& operator-= (const IntVectND<dim>& v) noexcept { smallend -= v; bigend -= v; return *this; }
624
627 BoxND operator- (const IntVectND<dim>& v) const noexcept { BoxND r(*this); r -= v; return r; }
628
642 BoxND chop (int dir, int chop_pnt) noexcept;
643
644 /*
645 * \brief Grow BoxND in all directions by given amount.
646 * NOTE: n_cell negative shrinks the BoxND by that number of cells.
647 */
649 BoxND& grow (int i) noexcept { smallend.diagShift(-i); bigend.diagShift(i); return *this; }
650
653 BoxND& grow (const IntVectND<dim>& v) noexcept { smallend -= v; bigend += v; return *this;}
654
660 BoxND& grow (int idir, int n_cell) noexcept { smallend.shift(idir, -n_cell); bigend.shift(idir, n_cell); return *this; }
661
663 BoxND& grow (Direction d, int n_cell) noexcept { return grow(static_cast<int>(d), n_cell); }
664
670 BoxND& growLo (int idir, int n_cell = 1) noexcept { smallend.shift(idir, -n_cell); return *this; }
671
673 BoxND& growLo (Direction d, int n_cell = 1) noexcept { return growLo(static_cast<int>(d), n_cell); }
674
681 BoxND& growHi (int idir, int n_cell = 1) noexcept { bigend.shift(idir,n_cell); return *this; }
682
684 BoxND& growHi (Direction d, int n_cell = 1) noexcept { return growHi(static_cast<int>(d), n_cell); }
685
688 BoxND& grow (Orientation face, int n_cell = 1) noexcept {
689 int idir = face.coordDir();
690 if (face.isLow()) {
691 smallend.shift(idir, -n_cell);
692 } else {
693 bigend.shift(idir,n_cell);
694 }
695 return *this;
696 }
697
706 BoxND& refine (int ref_ratio) noexcept {
707 return this->refine(IntVectND<dim>(ref_ratio));
708 }
709
710 /*
711 * \brief Refine BoxND by given (positive) refinement ratio.
712 * NOTE: if type(dir) = CELL centered: lo <- lo*ratio and
713 * hi <- (hi+1)*ratio - 1.
714 * NOTE: if type(dir) = NODE centered: lo <- lo*ratio and
715 * hi <- hi*ratio.
716 */
718 BoxND& refine (const IntVectND<dim>& ref_ratio) noexcept;
719
730 BoxND& coarsen (int ref_ratio) noexcept {
731 return this->coarsen(IntVectND<dim>(ref_ratio));
732 }
733
744 BoxND& coarsen (const IntVectND<dim>& ref_ratio) noexcept;
745
751 void next (IntVectND<dim> &) const noexcept;
752
762
763 [[nodiscard]] AMREX_GPU_HOST_DEVICE
764 bool isSquare() const noexcept;
765
777 [[nodiscard]] AMREX_GPU_HOST_DEVICE
778 bool coarsenable (const IntVectND<dim>& refrat, const IntVectND<dim>& min_width) const noexcept
779 {
780 if (!size().allGE(refrat*min_width)) {
781 return false;
782 } else {
783 BoxND testBox = *this;
784 testBox.coarsen(refrat);
785 testBox.refine (refrat);
786 return (*this == testBox);
787 }
788 }
789
801 [[nodiscard]] AMREX_GPU_HOST_DEVICE
802 bool coarsenable (int refrat, int min_width=1) const noexcept {
803 return coarsenable(IntVectND<dim>(refrat), IntVectND<dim>(min_width));
804 }
805
817 [[nodiscard]] AMREX_GPU_HOST_DEVICE
818 bool coarsenable (const IntVectND<dim>& refrat, int min_width=1) const noexcept
819 {
820 return coarsenable(refrat, IntVectND<dim>(min_width));
821 }
822
824 void normalize () noexcept
825 {
826 for (int idim=0; idim < dim; ++idim) {
827 if (this->length(idim) == 0) {
828 this->growHi(idim,1);
829 }
830 }
831 }
832
834 BoxND& makeSlab (int direction, int slab_index) noexcept
835 {
836 smallend[direction] = slab_index;
837 bigend[direction] = slab_index;
838 return *this;
839 }
840
852 [[nodiscard]] BoxIteratorND<dim> iterator () const noexcept {
853 return BoxIteratorND<dim>{*this};
854 }
855
857 static constexpr std::size_t ndims () noexcept {
858 return static_cast<std::size_t>(dim);
859 }
860
862 static constexpr int indims () noexcept {
863 return dim;
864 }
865
870 template<int new_dim>
872 BoxND<new_dim> shrink () const noexcept {
873 static_assert(new_dim <= dim);
874 auto lo = smallend.template shrink<new_dim>();
875 auto hi = bigend.template shrink<new_dim>();
876 auto typ = btype.template shrink<new_dim>();
877 return BoxND<new_dim>(lo, hi, typ);
878 }
879
885 template<int new_dim>
887 BoxND<new_dim> expand () const noexcept {
888 static_assert(new_dim >= dim);
889 auto lo = smallend.template expand<new_dim>(0);
890 auto hi = bigend.template expand<new_dim>(0);
891 auto typ = btype.template expand<new_dim>(IndexType::CellIndex::CELL);
892 return BoxND<new_dim>(lo, hi, typ);
893 }
894
899 template<int new_dim>
901 BoxND<new_dim> resize () const noexcept {
902 if constexpr (new_dim > dim) {
903 return expand<new_dim>();
904 } else {
905 return shrink<new_dim>();
906 }
907 }
908
909private:
910 IntVectND<dim> smallend;
911 IntVectND<dim> bigend;
912 IndexTypeND<dim> btype;
913};
914
915template<int dim>
918BoxND<dim>&
919BoxND<dim>::refine (const IntVectND<dim>& ref_ratio) noexcept
920{
921 if (ref_ratio != 1) {
922 IntVectND<dim> shft(1);
923 shft -= btype.ixType();
924 smallend *= ref_ratio;
925 bigend += shft;
926 bigend *= ref_ratio;
927 bigend -= shft;
928 }
929 return *this;
930}
931
932template<int dim>
936BoxND<dim>::coarsen (const IntVectND<dim>& ref_ratio) noexcept
937{
938 if (ref_ratio != 1)
939 {
940 smallend.coarsen(ref_ratio);
941
942 if (btype.any())
943 {
944 IntVectND<dim> off(0);
945 for (int dir = 0; dir < dim; dir++)
946 {
947 if (btype[dir]) {
948 if (bigend[dir]%ref_ratio[dir]) {
949 off.setVal(dir,1);
950 }
951 }
952 }
953 bigend.coarsen(ref_ratio);
954 bigend += off;
955 }
956 else
957 {
958 bigend.coarsen(ref_ratio);
959 }
960 }
961
962 return *this;
963}
964
965template<int dim>
970{
971 BL_ASSERT(typ.allGE(0) && typ.allLE(1));
972 IntVectND<dim> shft(typ - btype.ixType());
973 bigend += shft;
974 btype = IndexTypeND<dim>(typ);
975 return *this;
976}
977
978template<int dim>
983{
984 for (int dir = 0; dir < dim; dir++)
985 {
986 const auto typ = t[dir];
987 const auto bitval = btype[dir];
988 const int off = typ - bitval;
989 bigend.shift(dir,off);
990 btype.setType(dir, static_cast<IndexType::CellIndex>(typ));
991 }
992 return *this;
993}
994
995template<int dim>
1000{
1001 if (!(btype[dir]))
1002 {
1003 bigend.shift(dir,1);
1004 //
1005 // Set dir'th bit to 1 = IndexType::NODE.
1006 //
1007 btype.set(dir);
1008 }
1009 return *this;
1010}
1011
1012template<int dim>
1017{
1018 for (int i = 0; i < dim; ++i) {
1019 if ((btype[i] == 0)) {
1020 bigend.shift(i,1);
1021 }
1022 }
1023 btype.setall();
1024 return *this;
1025}
1026
1027template<int dim>
1032{
1033 if (btype[dir])
1034 {
1035 bigend.shift(dir,-1);
1036 //
1037 // Set dir'th bit to 0 = IndexType::CELL.
1038 //
1039 btype.unset(dir);
1040 }
1041 return *this;
1042}
1043
1044template<int dim>
1049{
1050 for (int i = 0 ; i < dim; ++i) {
1051 if (btype[i]) {
1052 bigend.shift(i,-1);
1053 }
1054 }
1055 btype.clear();
1056 return *this;
1057}
1058
1059template<int dim>
1062Long
1063BoxND<dim>::index (const IntVectND<dim>& v) const noexcept
1064{
1065 IntVectND<dim> vz = v - smallend;
1066 Long result = vz[0];
1067 Long mult = length(0);
1068 for (int i = 1 ; i < dim; ++i) {
1069 result += mult * vz[i];
1070 mult *= length(i);
1071 }
1072 return result;
1073}
1074
1075template<int dim>
1080{
1081 IntVectND<dim> result = smallend;
1082
1083 if constexpr (dim > 1) {
1084 GpuArray<Long, dim-1> mult{};
1085 mult[0] = length(0);
1086 for (int i = 1 ; i < dim-1; ++i) {
1087 mult[i] = mult[i-1] * length(i);
1088 }
1089 for (int i = dim-1 ; i > 0; --i) {
1090 Long idx = offset / mult[i-1];
1091 offset -= idx * mult[i-1];
1092 result[i] += static_cast<int>(idx);
1093 }
1094 }
1095
1096 result[0] += static_cast<int>(offset);
1097
1098 return result;
1099}
1100
1101template<int dim>
1102template <int N>
1103requires (1 <= N && N <= 3)
1108{
1109 Dim3 iv3d = atOffset(offset).dim3(0);
1110 return {{iv3d.x, iv3d.y, iv3d.z}};
1111}
1112
1113template<int dim>
1118 int sm_index,
1119 int n_cells) noexcept
1120{
1121 smallend.setVal(dir,sm_index);
1122 bigend.setVal(dir,sm_index+n_cells-1);
1123 return *this;
1124}
1125
1126template<int dim>
1129void
1130BoxND<dim>::next (IntVectND<dim>& p) const noexcept // NOLINT(readability-convert-member-functions-to-static)
1131{
1132 BL_ASSERT(contains(p));
1133
1134 ++p[0];
1135
1136 for (int i = 0 ; i < dim-1; ++i) {
1137 if (p[i] > bigend[i]) {
1138 p[i] = smallend[i];
1139 ++p[i+1];
1140 } else {
1141 break;
1142 }
1143 }
1144}
1145
1146template<int dim>
1149bool
1150BoxND<dim>::isSquare () const noexcept // NOLINT(readability-convert-member-functions-to-static)
1151{
1152 if constexpr (dim == 1) {
1153 return false; // can't build a square in 1-D
1154 } else {
1155 bool is_square = true;
1156 const IntVectND<dim>& sz = this->size();
1157 for (int i = 0 ; i < dim-1; ++i) {
1158 is_square = is_square && (sz[i] == sz[i+1]);
1159 }
1160 return is_square;
1161 }
1162}
1163
1164//
1165// Modified BoxND is low end, returned BoxND is high end.
1166// If CELL: chop_pnt included in high end.
1167// If NODE: chop_pnt included in both Boxes.
1168//
1169
1170template<int dim>
1172inline
1174BoxND<dim>::chop (int dir, int chop_pnt) noexcept
1175{
1176 //
1177 // Define new high end BoxND including chop_pnt.
1178 //
1179 IntVectND<dim> sm(smallend);
1180 IntVectND<dim> bg(bigend);
1181 sm.setVal(dir,chop_pnt);
1182 if (btype[dir])
1183 {
1184 //
1185 // NODE centered BoxND.
1186 //
1187 BL_ASSERT(chop_pnt > smallend[dir] && chop_pnt < bigend[dir]);
1188 //
1189 // Shrink original BoxND to just contain chop_pnt.
1190 //
1191 bigend.setVal(dir,chop_pnt);
1192 }
1193 else
1194 {
1195 //
1196 // CELL centered BoxND.
1197 //
1198 BL_ASSERT(chop_pnt > smallend[dir] && chop_pnt <= bigend[dir]);
1199 //
1200 // Shrink original BoxND to one below chop_pnt.
1201 //
1202 bigend.setVal(dir,chop_pnt-1);
1203 }
1204 return BoxND<dim>(sm,bg,btype);
1205}
1206
1207template<int dim>
1209inline
1211BoxND<dim>::shiftHalf (int dir, int num_halfs) noexcept
1212{
1213 const int nbit = (num_halfs<0 ? -num_halfs : num_halfs)%2;
1214 int nshift = num_halfs/2;
1215 //
1216 // Toggle btyp bit if num_halfs is odd.
1217 //
1218 const unsigned int bit_dir = btype[dir];
1219 if (nbit) {
1220 btype.flip(dir);
1221 }
1222 if (num_halfs < 0) {
1223 nshift -= (bit_dir ? nbit : 0);
1224 } else {
1225 nshift += (bit_dir ? 0 : nbit);
1226 }
1227 smallend.shift(dir,nshift);
1228 bigend.shift(dir,nshift);
1229 return *this;
1230}
1231
1232template<int dim>
1234inline
1237{
1238 for (int i = 0; i < dim; i++) {
1239 shiftHalf(i,iv[i]);
1240 }
1241 return *this;
1242}
1243
1245class BoxCommHelper
1246{
1247public:
1248
1249 explicit BoxCommHelper (const Box& bx, int* p_ = nullptr);
1250
1251 [[nodiscard]] int* data () const& noexcept { return p; }
1252 int* data () && = delete;
1253
1254 [[nodiscard]] Box make_box () const noexcept {
1255 return Box(IntVect(p), IntVect(p+AMREX_SPACEDIM), IntVect(p+2*AMREX_SPACEDIM));
1256 }
1257
1258 [[nodiscard]] static int size () noexcept { return 3*AMREX_SPACEDIM; }
1259
1260private:
1261 int* p;
1262 std::vector<int> v;
1263};
1265
1266class BoxConverter { // NOLINT
1267public:
1268 virtual Box doit (const Box& fine) const = 0; // NOLINT
1269 virtual BoxConverter* clone () const = 0; // NOLINT
1270 virtual ~BoxConverter () = default;
1271};
1272
1273void AllGatherBoxes (Vector<Box>& bxs, int n_extra_reserve=0);
1274
1285template<int dim>
1286[[nodiscard]]
1289BoxND<dim> grow (const BoxND<dim>& b, int i) noexcept
1290{
1291 BoxND<dim> result = b;
1292 result.grow(i);
1293 return result;
1294}
1295
1306template<int dim>
1307[[nodiscard]]
1310BoxND<dim> grow (const BoxND<dim>& b, const IntVectND<dim>& v) noexcept
1311{
1312 BoxND<dim> result = b;
1313 result.grow(v);
1314 return result;
1315}
1316
1328template<int dim>
1329[[nodiscard]]
1332BoxND<dim> grow (const BoxND<dim>& b, int idir, int n_cell) noexcept
1333{
1334 BoxND<dim> result = b;
1335 result.grow(idir, n_cell);
1336 return result;
1337}
1338
1350template<int dim>
1351[[nodiscard]]
1354BoxND<dim> grow (const BoxND<dim>& b, Direction d, int n_cell) noexcept
1355{
1356 return grow(b, static_cast<int>(d), n_cell);
1357}
1358
1359template<int dim>
1360[[nodiscard]]
1363BoxND<dim> growLo (const BoxND<dim>& b, int idir, int n_cell) noexcept
1364{
1365 BoxND<dim> result = b;
1366 result.growLo(idir, n_cell);
1367 return result;
1368}
1369
1370template<int dim>
1371[[nodiscard]]
1374BoxND<dim> growLo (const BoxND<dim>& b, Direction d, int n_cell) noexcept
1375{
1376 return growLo(b, static_cast<int>(d), n_cell);
1377}
1378
1379template<int dim>
1380[[nodiscard]]
1383BoxND<dim> growHi (const BoxND<dim>& b, int idir, int n_cell) noexcept
1384{
1385 BoxND<dim> result = b;
1386 result.growHi(idir, n_cell);
1387 return result;
1388}
1389
1390template<int dim>
1391[[nodiscard]]
1394BoxND<dim> growHi (const BoxND<dim>& b, Direction d, int n_cell) noexcept
1395{
1396 return growHi(b, static_cast<int>(d), n_cell);
1397}
1398
1414template<int dim>
1415[[nodiscard]]
1418BoxND<dim> coarsen (const BoxND<dim>& b, int ref_ratio) noexcept
1419{
1420 BoxND<dim> result = b;
1421 result.coarsen(IntVectND<dim>(ref_ratio));
1422 return result;
1423}
1424
1440template<int dim>
1441[[nodiscard]]
1444BoxND<dim> coarsen (const BoxND<dim>& b, const IntVectND<dim>& ref_ratio) noexcept
1445{
1446 BoxND<dim> result = b;
1447 result.coarsen(ref_ratio);
1448 return result;
1449}
1450
1464template<int dim>
1465[[nodiscard]]
1468BoxND<dim> refine (const BoxND<dim>& b, int ref_ratio) noexcept
1469{
1470 BoxND<dim> result = b;
1471 result.refine(IntVectND<dim>(ref_ratio));
1472 return result;
1473}
1474
1488template<int dim>
1489[[nodiscard]]
1492BoxND<dim> refine (const BoxND<dim>& b, const IntVectND<dim>& ref_ratio) noexcept
1493{
1494 BoxND<dim> result = b;
1495 result.refine(ref_ratio);
1496 return result;
1497}
1498
1500template<int dim>
1501[[nodiscard]]
1504BoxND<dim> shift (const BoxND<dim>& b, int dir, int nzones) noexcept
1505{
1506 BoxND<dim> result = b;
1507 result.shift(dir, nzones);
1508 return result;
1509}
1510
1511template<int dim>
1512[[nodiscard]]
1515BoxND<dim> shift (const BoxND<dim>& b, const IntVectND<dim>& nzones) noexcept
1516{
1517 BoxND<dim> result = b;
1518 result.shift(nzones);
1519 return result;
1520}
1521
1527template<int dim>
1528[[nodiscard]]
1531BoxND<dim> surroundingNodes (const BoxND<dim>& b, int dir) noexcept
1532{
1533 BoxND<dim> bx(b);
1534 bx.surroundingNodes(dir);
1535 return bx;
1536}
1537
1538template<int dim>
1539[[nodiscard]]
1543{
1544 return surroundingNodes(b, static_cast<int>(d));
1545}
1546
1551template<int dim>
1552[[nodiscard]]
1556{
1557 BoxND<dim> bx(b);
1558 bx.surroundingNodes();
1559 return bx;
1560}
1561
1563template<int dim>
1564[[nodiscard]]
1567BoxND<dim> convert (const BoxND<dim>& b, const IntVectND<dim>& typ) noexcept
1568{
1569 BoxND<dim> bx(b);
1570 bx.convert(typ);
1571 return bx;
1572}
1573
1574template<int dim>
1575[[nodiscard]]
1578BoxND<dim> convert (const BoxND<dim>& b, const IndexTypeND<dim>& typ) noexcept
1579{
1580 BoxND<dim> bx(b);
1581 bx.convert(typ);
1582 return bx;
1583}
1584
1591template<int dim>
1592[[nodiscard]]
1595BoxND<dim> enclosedCells (const BoxND<dim>& b, int dir) noexcept
1596{
1597 BoxND<dim> bx(b);
1598 bx.enclosedCells(dir);
1599 return bx;
1600}
1601
1602template<int dim>
1603[[nodiscard]]
1607{
1608 return enclosedCells(b, static_cast<int>(d));
1609}
1610
1615template<int dim>
1616[[nodiscard]]
1620{
1621 BoxND<dim> bx(b);
1622 bx.enclosedCells();
1623 return bx;
1624}
1625
1630template<int dim>
1631[[nodiscard]]
1634BoxND<dim> bdryLo (const BoxND<dim>& b, int dir, int len=1) noexcept
1635{
1636 IntVectND<dim> low(b.smallEnd());
1637 IntVectND<dim> hi(b.bigEnd());
1638 int sm = low[dir];
1639 low.setVal(dir,sm-len+1);
1640 hi.setVal(dir,sm);
1641 //
1642 // set dir'th bit to 1 = IndexType::NODE.
1643 //
1644 IndexTypeND<dim> typ(b.ixType());
1645 typ.set(dir);
1646 return BoxND<dim>(low,hi,typ);
1647}
1648
1653template<int dim>
1654[[nodiscard]]
1657BoxND<dim> bdryHi (const BoxND<dim>& b, int dir, int len=1) noexcept
1658{
1659 IntVectND<dim> low(b.smallEnd());
1660 IntVectND<dim> hi(b.bigEnd());
1661 auto const bitval = b.type()[dir];
1662 int bg = hi[dir] + 1 - bitval%2;
1663 low.setVal(dir,bg);
1664 hi.setVal(dir,bg+len-1);
1665 //
1666 // Set dir'th bit to 1 = IndexType::NODE.
1667 //
1668 IndexTypeND<dim> typ(b.ixType());
1669 typ.set(dir);
1670 return BoxND<dim>(low,hi,typ);
1671}
1672
1677template<int dim>
1678[[nodiscard]]
1681BoxND<dim> bdryNode (const BoxND<dim>& b, Orientation face, int len=1) noexcept
1682{
1683 int dir = face.coordDir();
1684 IntVectND<dim> low(b.smallEnd());
1685 IntVectND<dim> hi(b.bigEnd());
1686 if (face.isLow())
1687 {
1688 int sm = low[dir];
1689 low.setVal(dir,sm-len+1);
1690 hi.setVal(dir,sm);
1691 }
1692 else
1693 {
1694 int bitval = b.type()[dir];
1695 int bg = hi[dir] + 1 - bitval%2;
1696 low.setVal(dir,bg);
1697 hi.setVal(dir,bg+len-1);
1698 }
1699 //
1700 // Set dir'th bit to 1 = IndexType::NODE.
1701 //
1702 IndexTypeND<dim> typ(b.ixType());
1703 typ.set(dir);
1704 return BoxND<dim>(low,hi,typ);
1705}
1706
1719template<int dim>
1720[[nodiscard]]
1723BoxND<dim> adjCellLo (const BoxND<dim>& b, int dir, int len=1) noexcept
1724{
1725 BL_ASSERT(len > 0);
1726 IntVectND<dim> low(b.smallEnd());
1727 IntVectND<dim> hi(b.bigEnd());
1728 int sm = low[dir];
1729 low.setVal(dir,sm - len);
1730 hi.setVal(dir,sm - 1);
1731 //
1732 // Set dir'th bit to 0 = IndexType::CELL.
1733 //
1734 IndexTypeND<dim> typ(b.ixType());
1735 typ.unset(dir);
1736 return BoxND<dim>(low,hi,typ);
1737}
1738
1740template<int dim>
1741[[nodiscard]]
1744BoxND<dim> adjCellHi (const BoxND<dim>& b, int dir, int len=1) noexcept
1745{
1746 BL_ASSERT(len > 0);
1747 IntVectND<dim> low(b.smallEnd());
1748 IntVectND<dim> hi(b.bigEnd());
1749 int bitval = b.type()[dir];
1750 int bg = hi[dir] + 1 - bitval%2;
1751 low.setVal(dir,bg);
1752 hi.setVal(dir,bg + len - 1);
1753 //
1754 // Set dir'th bit to 0 = IndexType::CELL.
1755 //
1756 IndexTypeND<dim> typ(b.ixType());
1757 typ.unset(dir);
1758 return BoxND<dim>(low,hi,typ);
1759}
1760
1762template<int dim>
1763[[nodiscard]]
1766BoxND<dim> adjCell (const BoxND<dim>& b, Orientation face, int len=1) noexcept
1767{
1768 BL_ASSERT(len > 0);
1769 IntVectND<dim> low(b.smallEnd());
1770 IntVectND<dim> hi(b.bigEnd());
1771 int dir = face.coordDir();
1772 if (face.isLow())
1773 {
1774 int sm = low[dir];
1775 low.setVal(dir,sm - len);
1776 hi.setVal(dir,sm - 1);
1777 }
1778 else
1779 {
1780 int bitval = b.type()[dir];
1781 int bg = hi[dir] + 1 - bitval%2;
1782 low.setVal(dir,bg);
1783 hi.setVal(dir,bg + len - 1);
1784 }
1785 //
1786 // Set dir'th bit to 0 = IndexType::CELL.
1787 //
1788 IndexTypeND<dim> typ(b.ixType());
1789 typ.unset(dir);
1790 return BoxND<dim>(low,hi,typ);
1791}
1792
1798template<int dim>
1799[[nodiscard]]
1802BoxND<dim> minBox (const BoxND<dim>& b1, const BoxND<dim>& b2) noexcept
1803{
1804 BoxND<dim> result = b1;
1805 result.minBox(b2);
1806 return result;
1807}
1808
1810namespace detail {
1811 std::ostream& box_write (std::ostream& os, const int * smallend, const int * bigend,
1812 const int * type, int dim);
1813 std::istream& box_read (std::istream& is, int * smallend, int * bigend, int * type, int dim);
1814
1815 template<std::size_t...Ns, class T, class U>
1817 auto BoxSplit_imp (std::index_sequence<Ns...>,
1818 const T& lo, const T& hi, const U& typ) noexcept {
1819 return makeTuple(BoxND(get<Ns>(lo), get<Ns>(hi), get<Ns>(typ))...);
1820 }
1821}
1823
1825template<int dim>
1826std::ostream& operator<< (std::ostream& os, const BoxND<dim>& bx)
1827{
1828 IntVectND<dim> type = bx.type();
1829 return detail::box_write(os, bx.smallEnd().begin(), bx.bigEnd().begin(), type.begin(), dim);
1830}
1831
1833template<int dim>
1834std::istream& operator>> (std::istream& is, BoxND<dim>& bx) {
1835 IntVectND<dim> small;
1836 IntVectND<dim> big;
1837 IntVectND<dim> type;
1838 detail::box_read(is, small.begin(), big.begin(), type.begin(), dim);
1839 bx = BoxND<dim>{small, big, type};
1840 return is;
1841}
1842
1847template<int d, int...dims>
1850constexpr BoxND<detail::get_sum<d, dims...>()>
1851BoxCat (const BoxND<d>& bx, const BoxND<dims>&...boxes) noexcept {
1852 auto lo = IntVectCat(bx.smallEnd(), boxes.smallEnd()...);
1853 auto hi = IntVectCat(bx.bigEnd(), boxes.bigEnd()...);
1854 auto typ = IndexTypeCat(bx.ixType(), boxes.ixType()...);
1855 return BoxND<detail::get_sum<d, dims...>()>{lo, hi, typ};
1856}
1857
1862template<int d, int...dims>
1865constexpr GpuTuple<BoxND<d>, BoxND<dims>...>
1866BoxSplit (const BoxND<detail::get_sum<d, dims...>()>& bx) noexcept {
1867 auto lo = IntVectSplit<d, dims...>(bx.smallEnd());
1868 auto hi = IntVectSplit<d, dims...>(bx.bigEnd());
1869 auto typ = IndexTypeSplit<d, dims...>(bx.ixType());
1870 return detail::BoxSplit_imp(std::make_index_sequence<1 + sizeof...(dims)>(), lo, hi, typ);
1871}
1872
1877template<int new_dim, int old_dim>
1880constexpr BoxND<new_dim>
1881BoxShrink (const BoxND<old_dim>& bx) noexcept {
1882 return bx.template shrink<new_dim>();
1883}
1884
1890template<int new_dim, int old_dim>
1893constexpr BoxND<new_dim>
1894BoxExpand (const BoxND<old_dim>& bx) noexcept {
1895 return bx.template expand<new_dim>();
1896}
1897
1902template<int new_dim, int old_dim>
1905constexpr BoxND<new_dim>
1906BoxResize (const BoxND<old_dim>& bx) noexcept {
1907 return bx.template resize<new_dim>();
1908}
1909
1910template<int dim>
1911[[nodiscard]]
1915{
1916 return box.smallEnd();
1917}
1918
1919template<int dim>
1920[[nodiscard]]
1924{
1925 return box.bigEnd();
1926}
1927
1928template<int dim>
1929[[nodiscard]]
1933{
1934 return box.smallEnd();
1935}
1936
1937template<int dim>
1938[[nodiscard]]
1941IntVectND<dim> end_iv (BoxND<dim> const& box) noexcept
1942{
1943 return box.bigEnd() + 1;
1944}
1945
1946template<int dim>
1947[[nodiscard]]
1951{
1952 return box.bigEnd() - box.smallEnd() + 1;
1953}
1954
1955// Max of lower bound
1956template<int dim>
1957[[nodiscard]]
1960IntVectND<dim> max_lbound_iv (BoxND<dim> const& b1, BoxND<dim> const& b2) noexcept
1961{
1962 return max(b1.smallEnd(), b2.smallEnd());
1963}
1964
1965template<int dim>
1966[[nodiscard]]
1970{
1971 return max(b1.smallEnd(), lo);
1972}
1973
1974// Min of upper bound
1975template<int dim>
1976[[nodiscard]]
1979IntVectND<dim> min_ubound_iv (BoxND<dim> const& b1, BoxND<dim> const& b2) noexcept
1980{
1981 return min(b1.bigEnd(), b2.bigEnd());
1982}
1983
1984template<int dim>
1985[[nodiscard]]
1989{
1990 return min(b1.bigEnd(), hi);
1991}
1992
1993template<int dim>
1994requires ( 1<=dim && dim<=3 )
1995[[nodiscard]]
1998Dim3 lbound (BoxND<dim> const& box) noexcept
1999{
2000 return box.smallEnd().dim3();
2001}
2002
2003template<int dim>
2004requires ( 1<=dim && dim<=3 )
2005[[nodiscard]]
2008Dim3 ubound (BoxND<dim> const& box) noexcept
2009{
2010 return box.bigEnd().dim3();
2011}
2012
2013template<int dim>
2014requires ( 1<=dim && dim<=3 )
2015[[nodiscard]]
2018Dim3 begin (BoxND<dim> const& box) noexcept
2019{
2020 return box.smallEnd().dim3();
2021}
2022
2023template<int dim>
2024requires ( 1<=dim && dim<=3 )
2025[[nodiscard]]
2028Dim3 end (BoxND<dim> const& box) noexcept
2029{
2030 return (box.bigEnd() + 1).dim3(1);
2031}
2032
2033template<int dim>
2034requires ( 1<=dim && dim<=3 )
2035[[nodiscard]]
2038Dim3 length (BoxND<dim> const& box) noexcept
2039{
2040 return (box.bigEnd() - box.smallEnd() + 1).dim3(1);
2041}
2042
2043// Max of lower bound
2044template<int dim>
2045requires ( 1<=dim && dim<=3 )
2046[[nodiscard]]
2049Dim3 max_lbound (BoxND<dim> const& b1, BoxND<dim> const& b2) noexcept
2050{
2051 return max(b1.smallEnd(), b2.smallEnd()).dim3();
2052}
2053
2054template<int dim>
2055requires ( 1<=dim && dim<=3 )
2056[[nodiscard]]
2059Dim3 max_lbound (BoxND<dim> const& b1, Dim3 const& lo) noexcept
2060{
2061 return max(b1.smallEnd(), IntVectND<dim>(lo)).dim3();
2062}
2063
2064// Min of upper bound
2065template<int dim>
2066requires ( 1<=dim && dim<=3 )
2067[[nodiscard]]
2070Dim3 min_ubound (BoxND<dim> const& b1, BoxND<dim> const& b2) noexcept
2071{
2072 return min(b1.bigEnd(), b2.bigEnd()).dim3();
2073}
2074
2075template<int dim>
2076requires ( 1<=dim && dim<=3 )
2077[[nodiscard]]
2080Dim3 min_ubound (BoxND<dim> const& b1, Dim3 const& hi) noexcept
2081{
2082 return min(b1.bigEnd(), IntVectND<dim>(hi)).dim3();
2083}
2084
2085// Return a BoxND that covers all the argument Boxes in index
2086// space. The types are ignored. Thus, the arguments can have
2087// different index types, and the returned BoxND's index type has no
2088// meaning.
2089template<int dim>
2090[[nodiscard]]
2093{
2094 return b1;
2095}
2096
2097template<int dim>
2098[[nodiscard]]
2100BoxND<dim> getIndexBounds (BoxND<dim> const& b1, BoxND<dim> const& b2) noexcept
2101{
2102 BoxND<dim> b = b1;
2103 b.setType(b2.ixType());
2104 b.minBox(b2);
2105 return b;
2106}
2107
2108template <class T, class ... Ts>
2109[[nodiscard]]
2111auto getIndexBounds (T const& b1, T const& b2, Ts const& ... b3) noexcept
2112{
2113 return getIndexBounds(getIndexBounds(b1,b2),b3...);
2114}
2115
2116template<int dim>
2117[[nodiscard]]
2120IntVectND<dim> getCell (BoxND<dim> const* boxes, int nboxes, Long icell) noexcept
2121{
2122 int ibox;
2123 Long ncells_subtotal = 0;
2124 Long offset = 0;
2125 for (ibox = 0; ibox < nboxes; ++ibox) {
2126 const Long n = boxes[ibox].numPts();
2127 ncells_subtotal += n;
2128 if (icell < ncells_subtotal) {
2129 offset = icell - (ncells_subtotal - n);
2130 break;
2131 }
2132 }
2133 return boxes[ibox].atOffset(offset);
2134}
2135
2136template<int dim>
2137[[nodiscard]]
2140BoxND<dim> makeSlab (BoxND<dim> const& b, int direction, int slab_index) noexcept
2141{
2142 BoxND<dim> r = b;
2143 r.makeSlab(direction,slab_index);
2144 return r;
2145}
2146
2147template<int dim=AMREX_SPACEDIM>
2148requires ( 1<=dim && dim<=3 )
2149[[nodiscard]]
2153{
2154 Dim3 p3d{.x = i, .y = j, .z = k};
2155 IntVectND<dim> vect{p3d};
2156 return BoxND<dim>{vect, vect, typ};
2157}
2158
2159template<int dim>
2160[[nodiscard]]
2164{
2165 return BoxND<dim>{vect, vect, typ};
2166}
2167
2168template<int dim>
2170{
2171 std::uint64_t npts;
2174
2176 : npts(box.numPts()),
2177 lo (box.smallEnd())
2178 {
2179 std::uint64_t mult = 1;
2180 for (int i=0; i<dim-1; ++i) {
2181 mult *= box.length(i);
2182 fdm[i] = Math::FastDivmodU64(mult);
2183 }
2184 }
2185
2187 IntVectND<dim> intVect (std::uint64_t icell) const
2188 {
2189 IntVectND<dim> retval = lo;
2190
2191 for (int i=dim-1; i>0; --i) {
2192 std::uint64_t quotient, remainder;
2193 fdm[i-1](quotient, remainder, icell);
2194 retval[i] += quotient;
2195 icell = remainder;
2196 }
2197
2198 retval[0] += icell;
2199
2200 return retval;
2201 }
2202
2204 Dim3 operator()(std::uint64_t icell) const
2205 requires (1 <= dim && dim <= 3)
2206 {
2207 return intVect(icell).dim3();
2208 }
2209
2211 std::uint64_t numPts () const { return npts; }
2212};
2213
2214template<>
2216{
2217 std::uint64_t npts;
2218
2219 int lo;
2220
2222 : npts(box.numPts()),
2223 lo(box.smallEnd(0))
2224 {}
2225
2227 IntVectND<1> intVect (std::uint64_t icell) const
2228 {
2229 return IntVectND<1>{int(icell)+lo};
2230 }
2231
2233 Dim3 operator() (std::uint64_t icell) const
2234 {
2235 return Dim3{.x = int(icell)+lo, .y = 0, .z = 0};
2236 }
2237
2239 std::uint64_t numPts () const { return npts; }
2240};
2241
2243
2244}
2245
2246#endif /*AMREX_BOX_H*/
#define AMREX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:49
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#define AMREX_FORCE_INLINE
Definition AMReX_Extension.H:119
#define AMREX_IF_ON_HOST(CODE)
Definition AMReX_GpuQualifiers.H:58
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
int idir
Definition AMReX_HypreMLABecLap.cpp:1143
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1139
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Definition AMReX_Box.H:1266
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
__host__ __device__ const int * hiVect() &&=delete
__host__ __device__ BoxND & grow(const IntVectND< dim > &v) noexcept
Grow BoxND in each direction by specified amount.
Definition AMReX_Box.H:653
__host__ __device__ GpuArray< int, 3 > loVect3d() const noexcept
Definition AMReX_Box.H:173
__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__ bool operator<=(const BoxND &rhs) const noexcept
Definition AMReX_Box.H:321
__host__ __device__ BoxND & grow(int i) noexcept
Definition AMReX_Box.H:649
__host__ __device__ GpuArray< int, 3 > length3d() const noexcept
Definition AMReX_Box.H:165
__host__ __device__ const IntVectND< dim > & bigEnd() const &noexcept
Return the inclusive upper bound of the box.
Definition AMReX_Box.H:124
__host__ __device__ const int * hiVect() const &noexcept
Return a constant pointer the array of high end coordinates. Useful for calls to FORTRAN.
Definition AMReX_Box.H:195
__host__ __device__ BoxND & operator&=(const BoxND &rhs) noexcept
Intersect this BoxND with its argument. The Boxes MUST be of the same type.
Definition AMReX_Box.H:591
__host__ __device__ bool isEmpty() const noexcept
Checks if it is an empty BoxND.
Definition AMReX_Box.H:208
__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:862
__host__ __device__ int bigEnd(int dir) const noexcept
Return the inclusive upper bound in the given direction.
Definition AMReX_Box.H:132
__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:901
__host__ __device__ BoxND & operator-=(const IntVectND< dim > &v) noexcept
Shift BoxND (relative) by given IntVectND<dim>.
Definition AMReX_Box.H:623
__host__ __device__ BoxND & surroundingNodes(Direction d) noexcept
Definition AMReX_Box.H:568
__host__ __device__ BoxND & makeSlab(int direction, int slab_index) noexcept
Definition AMReX_Box.H:834
__host__ __device__ BoxND & setBig(const IntVectND< dim > &bg) noexcept
Redefine the big end of the BoxND.
Definition AMReX_Box.H:495
__host__ __device__ int shortside() const noexcept
Return length of shortest side. Ignores type.
Definition AMReX_Box.H:461
__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:605
__host__ __device__ Long volume() const noexcept
Return the volume, in indexing space, of region enclosed by this BoxND. This is identical to numPts()...
Definition AMReX_Box.H:402
__host__ __device__ BoxND & shift(const IntVectND< dim > &iv) noexcept
Equivalent to b.shift(0,iv[0]).shift(1,iv[1]) ....
Definition AMReX_Box.H:521
__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:887
__host__ __device__ bool isSquare() const noexcept
Definition AMReX_Box.H:1150
__host__ __device__ bool contains(const Dim3 &p) const noexcept
Return true if argument is contained within BoxND.
Definition AMReX_Box.H:222
__host__ __device__ bool cellCentered() const noexcept
Return true if BoxND is cell-centered in all indexing directions.
Definition AMReX_Box.H:335
__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:802
__host__ __device__ BoxND & operator+=(const IntVectND< dim > &v) noexcept
Shift BoxND (relative) by given IntVectND<dim>.
Definition AMReX_Box.H:615
__host__ __device__ const int * loVect() &&=delete
__host__ __device__ Long numPts() const noexcept
Return the number of points contained in the BoxND.
Definition AMReX_Box.H:364
__host__ __device__ bool operator>=(const BoxND &rhs) const noexcept
Definition AMReX_Box.H:329
__host__ __device__ bool sameType(const BoxND &b) const noexcept
Return true if Boxes have same type.
Definition AMReX_Box.H:302
__host__ __device__ bool operator!=(const BoxND &b) const noexcept
Return true if Boxes differ (including type).
Definition AMReX_Box.H:310
__host__ __device__ BoxND & setBig(int dir, int bg_index) noexcept
Redefine the big end of the BoxND.
Definition AMReX_Box.H:499
__host__ __device__ BoxND operator-(const IntVectND< dim > &v) const noexcept
Shift BoxND (relative) by given IntVectND<dim>.
Definition AMReX_Box.H:627
__host__ __device__ BoxND & refine(const IntVectND< dim > &ref_ratio) noexcept
Definition AMReX_Box.H:919
__host__ __device__ bool sameSize(const BoxND &b) const noexcept
Return true is Boxes same size, ie translates of each other,. It is an error if they have different t...
Definition AMReX_Box.H:295
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:155
__host__ __device__ GpuArray< int, 3 > hiVect3d() const noexcept
Definition AMReX_Box.H:181
__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:982
__host__ __device__ bool contains(const IntVectND< dim > &p) const noexcept
Return true if argument is contained within BoxND.
Definition AMReX_Box.H:216
__host__ __device__ int shortside(int &dir) const noexcept
Return length of shortest side. dir is modified to give direction with shortest side: 0....
Definition AMReX_Box.H:445
__host__ __device__ IntVectND< dim > size() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:148
__host__ __device__ int operator[](Orientation face) const noexcept
Return the coordinate normal to given face.
Definition AMReX_Box.H:201
__host__ __device__ bool strictly_contains(const Dim3 &p) const noexcept
Return true if argument is strictly contained within BoxND.
Definition AMReX_Box.H:267
__host__ __device__ BoxND & shift(int dir, int nzones) noexcept
Shift this BoxND nzones indexing positions in coordinate direction dir.
Definition AMReX_Box.H:517
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:852
__host__ static __device__ constexpr std::size_t ndims() noexcept
Definition AMReX_Box.H:857
__host__ __device__ bool coarsenable(const IntVectND< dim > &refrat, const IntVectND< dim > &min_width) const noexcept
Return whether this Box is coarsenable.
Definition AMReX_Box.H:778
__host__ __device__ BoxND & setSmall(int dir, int sm_index) noexcept
Redefine the small end of the BoxND.
Definition AMReX_Box.H:491
__host__ __device__ double d_numPts() const noexcept
Return the number of points contained in the BoxND. This is intended for use only in diagnostic messa...
Definition AMReX_Box.H:384
__host__ __device__ IntVectND< dim > atOffset(Long offset) const noexcept
Given the offset, compute IntVectND<dim>
Definition AMReX_Box.H:1079
__host__ __device__ BoxND operator+(const IntVectND< dim > &v) const noexcept
Shift BoxND (relative) by given IntVectND<dim>.
Definition AMReX_Box.H:619
__host__ __device__ BoxND & grow(int idir, int n_cell) noexcept
Grow the BoxND on the low and high end by n_cell cells in direction idir.
Definition AMReX_Box.H:660
__host__ __device__ const int * loVect() const &noexcept
Return a constant pointer the array of low end coordinates. Useful for calls to FORTRAN.
Definition AMReX_Box.H:190
__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:670
__host__ __device__ bool contains(const BoxND &b) const noexcept
Return true if argument is contained within BoxND. It is an error if the Boxes have different types.
Definition AMReX_Box.H:242
__host__ __device__ int longside(int &dir) const noexcept
Return length of longest side. dir is modified to give direction with longest side: 0....
Definition AMReX_Box.H:419
__host__ __device__ IndexTypeND< dim > ixType() const noexcept
Return the indexing type.
Definition AMReX_Box.H:136
__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:681
__host__ __device__ bool strictly_contains(const IntVectND< dim > &p) const noexcept
Return true if argument is strictly contained within BoxND.
Definition AMReX_Box.H:250
__host__ __device__ IntVectND< dim > type() const noexcept
Return the indexing type.
Definition AMReX_Box.H:140
__host__ __device__ BoxND & growHi(Direction d, int n_cell=1) noexcept
Definition AMReX_Box.H:684
__host__ __device__ BoxND & growLo(Direction d, int n_cell=1) noexcept
Definition AMReX_Box.H:673
__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:730
__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:706
__host__ __device__ BoxND & setRange(int dir, int sm_index, int n_cells=1) noexcept
Set the entire range in a given direction, starting at sm_index with length n_cells....
Definition AMReX_Box.H:1117
__host__ __device__ bool strictly_contains(int i, int j, int k) const noexcept
Return true if argument is strictly contained within BoxND.
Definition AMReX_Box.H:276
__host__ __device__ Long index(const IntVectND< dim > &v) const noexcept
Return offset of point from smallend; i.e. index(smallend) -> 0, bigend would return numPts()-1....
Definition AMReX_Box.H:1063
__host__ __device__ bool contains(int i, int j, int k) const noexcept
Return true if argument is contained within BoxND.
Definition AMReX_Box.H:231
__host__ __device__ BoxND & setType(const IndexTypeND< dim > &t) noexcept
Set indexing type.
Definition AMReX_Box.H:513
__host__ __device__ BoxND & grow(Direction d, int n_cell) noexcept
Definition AMReX_Box.H:663
__host__ __device__ int longside() const noexcept
Return length of longest side. Ignores type.
Definition AMReX_Box.H:435
__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:1048
__host__ __device__ GpuArray< int, 3 > atOffset3d(Long offset) const noexcept
Definition AMReX_Box.H:1107
__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:1130
__host__ __device__ void normalize() noexcept
Definition AMReX_Box.H:824
__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:872
__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:936
__host__ __device__ BoxND & shiftHalf(int dir, int num_halfs) noexcept
This member shifts the BoxND by "half" indices, thereby converting the BoxND from type CELL to NODE a...
Definition AMReX_Box.H:1211
__host__ __device__ bool operator==(const BoxND &b) const noexcept
Return true if Boxes are identical (including type).
Definition AMReX_Box.H:306
__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:759
__host__ __device__ int smallEnd(int dir) const &noexcept
Return the inclusive lower bound in the given direction.
Definition AMReX_Box.H:120
__host__ __device__ BoxND & convert(const IntVectND< dim > &typ) noexcept
Convert the BoxND from the current type into the argument type. This may change the BoxND coordinates...
Definition AMReX_Box.H:969
__host__ __device__ BoxND & shiftHalf(const IntVectND< dim > &iv) noexcept
Equivalent to b.shiftHalf(0,iv[0]).shiftHalf(1,iv[1]) ...
Definition AMReX_Box.H:1236
__host__ __device__ bool coarsenable(const IntVectND< dim > &refrat, int min_width=1) const noexcept
Return whether this Box is coarsenable.
Definition AMReX_Box.H:818
__host__ __device__ BoxND & enclosedCells(Direction d) noexcept
Convert to CELL type in given Direction d.
Definition AMReX_Box.H:580
__host__ __device__ bool operator<(const BoxND &rhs) const noexcept
Definition AMReX_Box.H:313
__host__ __device__ bool intersects(const BoxND &b) const noexcept
Return true if Boxes have non-null intersections. It is an error if the Boxes have different types.
Definition AMReX_Box.H:288
__host__ __device__ IndexType::CellIndex type(int dir) const noexcept
Return the indexing type in the specified direction.
Definition AMReX_Box.H:144
__host__ __device__ BoxND chop(int dir, int chop_pnt) noexcept
Chop the BoxND at chop_pnt in the dir direction return one BoxND and modify the object BoxND....
Definition AMReX_Box.H:1174
__host__ __device__ bool ok() const noexcept
Checks if it is a proper BoxND (including a valid type).
Definition AMReX_Box.H:212
__host__ __device__ int length(int dir) const noexcept
Return the length of the BoxND in given direction.
Definition AMReX_Box.H:162
friend class BoxCommHelper
Definition AMReX_Box.H:51
__host__ __device__ BoxND & surroundingNodes() noexcept
Convert to NODE type in all directions.
Definition AMReX_Box.H:1016
__host__ __device__ BoxND & grow(Orientation face, int n_cell=1) noexcept
Grow in the direction of the given face.
Definition AMReX_Box.H:688
__host__ __device__ BoxND operator&(const BoxND &rhs) const noexcept
Return BoxND that is intersection of this BoxND and argument. The Boxes MUST be of same type.
Definition AMReX_Box.H:587
__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:112
__host__ __device__ bool operator>(const BoxND &rhs) const noexcept
Definition AMReX_Box.H:325
__host__ __device__ BoxND & setSmall(const IntVectND< dim > &sm) noexcept
Definition AMReX_Box.H:487
__host__ __device__ bool strictly_contains(const BoxND &b) const noexcept
Return true if argument is strictly contained within BoxND. It is an error if the Boxes have differen...
Definition AMReX_Box.H:259
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:70
__host__ __device__ void unset(int dir) noexcept
Set IndexTypeND to be CELL based in direction dir.
Definition AMReX_IndexType.H:73
An Integer Vector in dim-Dimensional Space.
Definition AMReX_IntVect.H:149
__host__ __device__ IntVectND & setVal(int i, int val) noexcept
Set i'th coordinate of IntVectND to val.
Definition AMReX_IntVect.H:373
__host__ __device__ constexpr int * begin() noexcept
Returns a pointer to the first element of the IntVectND.
Definition AMReX_IntVect.H:357
__host__ __device__ IntVectND< dim > & coarsen(const IntVectND< dim > &p) noexcept
Modify IntVectND<dim> by component-wise integer projection.
Definition AMReX_IntVect.H:941
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:29
amrex_long Long
Definition AMReX_INT.H:30
__host__ __device__ Dim3 ubound(Array4< T > const &a) noexcept
Return the inclusive upper bounds of an Array4 in Dim3 form.
Definition AMReX_Array4.H:1359
__host__ __device__ Dim3 length(Array4< T > const &a) noexcept
Return the spatial extents of an Array4 in Dim3 form.
Definition AMReX_Array4.H:1373
__host__ __device__ Dim3 lbound(Array4< T > const &a) noexcept
Return the inclusive lower bounds of an Array4 in Dim3 form.
Definition AMReX_Array4.H:1345
__host__ __device__ BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) coarsening ratio.
Definition AMReX_Box.H:1418
__host__ __device__ BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition AMReX_Box.H:1289
__host__ __device__ BoxND< dim > refine(const BoxND< dim > &b, int ref_ratio) noexcept
Definition AMReX_Box.H:1468
int MPI_Datatype
Definition AMReX_ccse-mpi.H:53
Definition AMReX_Amr.cpp:50
std::ostream & operator<<(std::ostream &os, AmrMesh const &amr_mesh)
Stream helper; forwards to the friend declared inside AmrMesh.
Definition AMReX_AmrMesh.cpp:1306
__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:1744
__host__ __device__ BoxND< dim > makeSingleCellBox(int i, int j, int k, IndexTypeND< dim > typ=IndexTypeND< dim >::TheCellType())
Definition AMReX_Box.H:2152
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1567
__host__ __device__ BoxND< dim > makeSlab(BoxND< dim > const &b, int direction, int slab_index) noexcept
Definition AMReX_Box.H:2140
__host__ __device__ constexpr GpuTuple< detail::tuple_decay_t< Ts >... > makeTuple(Ts &&... args)
Definition AMReX_Tuple.H:269
__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:1723
__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:1531
__host__ __device__ IntVectND< dim > end_iv(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:1941
__host__ __device__ Dim3 begin(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2018
__host__ __device__ IntVectND< dim > lbound_iv(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:1914
__host__ __device__ Dim3 min_ubound(BoxND< dim > const &b1, BoxND< dim > const &b2) noexcept
Definition AMReX_Box.H:2070
__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:1634
__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:1288
__host__ __device__ IntVectND< dim > min_ubound_iv(BoxND< dim > const &b1, BoxND< dim > const &b2) noexcept
Definition AMReX_Box.H:1979
__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:1802
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
__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:1906
__host__ __device__ Dim3 max_lbound(BoxND< dim > const &b1, BoxND< dim > const &b2) noexcept
Definition AMReX_Box.H:2049
__host__ __device__ constexpr const T & min(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:25
__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:1766
__host__ __device__ IntVectND< dim > begin_iv(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:1932
__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:1504
__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:1881
__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:1595
__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:293
__host__ __device__ IntVectND< dim > length_iv(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:1950
Direction
Definition AMReX_Orientation.H:14
__host__ __device__ constexpr BoxND< detail::get_sum< d, dims... >()> BoxCat(const BoxND< d > &bx, const BoxND< dims > &...boxes) noexcept
Return a BoxND obtained by concatenating the input BoxNDs. The dimension of the return value equals t...
Definition AMReX_Box.H:1851
__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:1681
__host__ __device__ constexpr const T & max(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:45
__host__ __device__ BoxND< dim > growLo(const BoxND< dim > &b, int idir, int n_cell) noexcept
Definition AMReX_Box.H:1363
__host__ __device__ IntVectND< dim > ubound_iv(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:1923
__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:1894
__host__ __device__ IntVectND< dim > getCell(BoxND< dim > const *boxes, int nboxes, Long icell) noexcept
Definition AMReX_Box.H:2120
const int[]
Definition AMReX_BLProfiler.cpp:1664
void AllGatherBoxes(Vector< Box > &bxs, int n_extra_reserve)
Definition AMReX_Box.cpp:124
std::istream & operator>>(std::istream &is, BoxND< dim > &bx)
Read from istream.
Definition AMReX_Box.H:1834
__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:1272
__host__ __device__ Dim3 end(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2028
__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:1657
__host__ __device__ BoxND< dim > growHi(const BoxND< dim > &b, int idir, int n_cell) noexcept
Definition AMReX_Box.H:1383
BoxND< dim > getIndexBounds(BoxND< dim > const &b1) noexcept
Definition AMReX_Box.H:2092
__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:1866
__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:313
__host__ __device__ IntVectND< dim > max_lbound_iv(BoxND< dim > const &b1, BoxND< dim > const &b2) noexcept
Definition AMReX_Box.H:1960
A multidimensional array accessor.
Definition AMReX_Array4.H:285
IntVectND< N > end
Exclusive upper bounds.
Definition AMReX_Array4.H:299
IntVectND< N > begin
Inclusive lower bounds.
Definition AMReX_Array4.H:298
BoxIndexerND(BoxND< 1 > const &box)
Definition AMReX_Box.H:2221
std::uint64_t npts
Definition AMReX_Box.H:2217
__host__ __device__ std::uint64_t numPts() const
Definition AMReX_Box.H:2239
__host__ __device__ IntVectND< 1 > intVect(std::uint64_t icell) const
Definition AMReX_Box.H:2227
int lo
Definition AMReX_Box.H:2219
Definition AMReX_Box.H:2170
__host__ __device__ IntVectND< dim > intVect(std::uint64_t icell) const
Definition AMReX_Box.H:2187
BoxIndexerND(BoxND< dim > const &box)
Definition AMReX_Box.H:2175
IntVectND< dim > lo
Definition AMReX_Box.H:2173
__host__ __device__ std::uint64_t numPts() const
Definition AMReX_Box.H:2211
Math::FastDivmodU64 fdm[dim-1]
Definition AMReX_Box.H:2172
__host__ __device__ Dim3 operator()(std::uint64_t icell) const
Definition AMReX_Box.H:2204
std::uint64_t npts
Definition AMReX_Box.H:2171
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:13
int x
Definition AMReX_Dim3.H:13
int z
Definition AMReX_Dim3.H:13
int y
Definition AMReX_Dim3.H:13
Fixed-size array that can be used on GPU.
Definition AMReX_Array.H:43
Definition AMReX_Math.H:493