Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_Geometry.H
Go to the documentation of this file.
1#ifndef AMREX_GEOMETRY_H_
2#define AMREX_GEOMETRY_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_Array.H>
6#include <AMReX_CoordSys.H>
8#include <AMReX_RealBox.H>
9#include <AMReX_Periodicity.H>
10
11#ifdef AMREX_USE_OMP
12#include <omp.h>
13#endif
14
15#include <iosfwd>
16#include <map>
17#include <numbers>
18
19namespace amrex {
20
21class MultiFab;
22class DistributionMapping;
23class BoxArray;
24
26{
29 const Real* CellSize () const noexcept { return dx; }
30 //
32 Real CellSize (int dir) const noexcept { return dx[dir]; }
35 const Real* ProbLo () const noexcept { return prob_domain.lo(); }
36 //
38 Real ProbLo (int dir) const noexcept { return prob_domain.lo(dir); }
41 const Real* ProbHi () const noexcept { return prob_domain.hi(); }
42 //
44 Real ProbHi (int dir) const noexcept { return prob_domain.hi(dir); }
47 const Box& Domain () const noexcept { return domain; }
50 int isPeriodic (const int i) const noexcept { return is_periodic[i]; }
53 int Coord () const noexcept { return coord; }
54
55public:
58 Real dx[AMREX_SPACEDIM];
59 int is_periodic[AMREX_SPACEDIM];
60 int coord;
61};
62
73 :
74 public CoordSys
75{
76public:
81 Geometry () noexcept;
82
97 explicit Geometry (const Box& dom,
98 const RealBox* rb = nullptr,
99 int coord = -1,
100 int const* is_per = nullptr);
115 Geometry (const Box& dom, const RealBox& rb, int coord,
116 Array<int,AMREX_SPACEDIM> const& is_per);
117
118 ~Geometry () = default;
119 Geometry (const Geometry& rhs) = default;
120 Geometry (Geometry&& rhs) noexcept = default;
121 Geometry& operator= (const Geometry& rhs) = default;
122 Geometry& operator= (Geometry&& rhs) noexcept = default;
123
125 [[nodiscard]] GeometryData data() const noexcept {
126 return {
127 .prob_domain = prob_domain,
128 .domain = domain,
129 .dx = {AMREX_D_DECL(dx[0],dx[1],dx[2])},
130 .is_periodic = {AMREX_D_DECL(is_periodic[0], is_periodic[1], is_periodic[2])},
131 .coord = static_cast<int>(c_sys)
132 };
133 }
134
136 static void Setup (const RealBox* rb = nullptr, int coord = -1, int const* isper = nullptr);
137 static void ResetDefaultProbDomain (const RealBox& rb) noexcept;
138 static void ResetDefaultPeriodicity (const Array<int,AMREX_SPACEDIM>& is_per) noexcept;
139 static void ResetDefaultCoord (int coord) noexcept;
140
156 void define (const Box& dom, const RealBox* rb = nullptr, int coord = -1, int const* is_per = nullptr);
157
173 void define (const Box& dom, const RealBox& rb, int coord, Array<int,AMREX_SPACEDIM> const& is_per);
174
176 [[nodiscard]] const RealBox& ProbDomain () const noexcept { return prob_domain; }
178 void ProbDomain (const RealBox& rb)
179 {
180 prob_domain = rb;
182 }
184 [[nodiscard]] const Real* ProbLo () const noexcept { return prob_domain.lo(); }
186 [[nodiscard]] const Real* ProbHi () const noexcept { return prob_domain.hi(); }
188 [[nodiscard]] Real ProbLo (int dir) const noexcept { return prob_domain.lo(dir); }
190 [[nodiscard]] Real ProbHi (int dir) const noexcept { return prob_domain.hi(dir); }
191
192 [[nodiscard]] GpuArray<Real,AMREX_SPACEDIM> ProbLoArray () const noexcept {
193 return {{AMREX_D_DECL(prob_domain.lo(0),prob_domain.lo(1),prob_domain.lo(2))}};
194 }
195
196 [[nodiscard]] GpuArray<Real,AMREX_SPACEDIM> ProbHiArray () const noexcept {
197 return {{AMREX_D_DECL(prob_domain.hi(0),prob_domain.hi(1),prob_domain.hi(2))}};
198 }
199
201 return roundoff_lo;
202 }
203
205 return roundoff_hi;
206 }
207
209 [[nodiscard]] Real ProbSize () const noexcept
210 {
211 return AMREX_D_TERM(prob_domain.length(0),*prob_domain.length(1),*prob_domain.length(2));
212 }
214 [[nodiscard]] Real ProbLength (int dir) const noexcept { return prob_domain.length(dir); }
216 [[nodiscard]] const Box& Domain () const noexcept { return domain; }
218 void Domain (const Box& bx)
219 {
221 domain = bx;
223 }
225 void GetVolume (MultiFab& vol,
226 const BoxArray& grds,
227 const DistributionMapping& dm,
228 int grow) const;
230 void GetVolume (MultiFab& vol) const;
231
232 void GetVolume (FArrayBox& vol,
233 const BoxArray& grds,
234 int idx,
235 int grow) const;
236
239 static Real Volume (const IntVect& point, const GeometryData& geomdata)
240 {
241 const auto *dx = geomdata.CellSize();
242
243 Real vol;
244
245#if AMREX_SPACEDIM == 1
246
247 auto coord = geomdata.Coord();
248
249 if (coord == CoordSys::cartesian) {
250 // Cartesian
251
252 vol = dx[0];
253 }
254 else {
255 // Spherical
256
257 Real rl = geomdata.ProbLo()[0] + static_cast<Real>(point[0]) * dx[0];
258 Real rr = rl + dx[0];
259
260 constexpr Real pi = std::numbers::pi_v<Real>;
261 vol = (4.0_rt / 3.0_rt) * pi * dx[0] * (rl * rl + rl * rr + rr * rr);
262 }
263
264#elif AMREX_SPACEDIM == 2
265
266 auto coord = geomdata.Coord();
267
268 if (coord == CoordSys::cartesian) {
269 // Cartesian
270
271 vol = dx[0] * dx[1];
272 }
273 else if (coord == CoordSys::RZ) {
274 // Cylindrical
275
276 Real r_l = geomdata.ProbLo()[0] + static_cast<Real>(point[0]) * dx[0];
277 Real r_r = geomdata.ProbLo()[0] + static_cast<Real>(point[0]+1) * dx[0];
278
279 constexpr Real pi = std::numbers::pi_v<Real>;
280 vol = pi * (r_l + r_r) * dx[0] * dx[1];
281 }
282 else {
283 // Spherical
284
285 Real r_l = geomdata.ProbLo()[0] + static_cast<Real>(point[0]) * dx[0];
286 Real r_r = geomdata.ProbLo()[0] + static_cast<Real>(point[0]+1) * dx[0];
287
288 Real theta_l = geomdata.ProbLo()[1] + static_cast<Real>(point[1]) * dx[1];
289 Real theta_r = geomdata.ProbLo()[1] + static_cast<Real>(point[1]+1) * dx[1];
290
291 constexpr Real twoThirdsPi = static_cast<Real>((2.0 / 3.0) * std::numbers::pi_v<Real>);
292 vol = twoThirdsPi * (std::cos(theta_l) - std::cos(theta_r)) * dx[0] *
293 (r_r*r_r + r_r*r_l + r_l*r_l);
294 }
295
296#else
297
299
300 // Cartesian
301
302 vol = dx[0] * dx[1] * dx[2];
303
304#endif
305
306 return vol;
307 }
308
313 void GetDLogA (MultiFab& dloga,
314 const BoxArray& grds,
315 const DistributionMapping& dm,
316 int dir,
317 int grow) const;
322 void GetFaceArea (MultiFab& area,
323 const BoxArray& grds,
324 const DistributionMapping& dm,
325 int dir,
326 int grow) const;
328 void GetFaceArea (MultiFab& area,
329 int dir) const;
330
331 void GetFaceArea (FArrayBox& area,
332 const BoxArray& grds,
333 int idx,
334 int dir,
335 int grow) const;
337 [[nodiscard]] bool isPeriodic (int dir) const noexcept { return is_periodic[dir]; }
339 [[nodiscard]] bool isAnyPeriodic () const noexcept
340 {
341 return AMREX_D_TERM(isPeriodic(0),||isPeriodic(1),||isPeriodic(2));
342 }
344 [[nodiscard]] bool isAllPeriodic () const noexcept
345 {
346 return AMREX_D_TERM(isPeriodic(0),&&isPeriodic(1),&&isPeriodic(2));
347 }
348 [[nodiscard]] Array<int,AMREX_SPACEDIM> isPeriodic () const noexcept {
349 return {{AMREX_D_DECL(static_cast<int>(is_periodic[0]),
350 static_cast<int>(is_periodic[1]),
351 static_cast<int>(is_periodic[2]))}};
352 }
353 [[nodiscard]] GpuArray<int,AMREX_SPACEDIM> isPeriodicArray () const noexcept {
354 return {{AMREX_D_DECL(static_cast<int>(is_periodic[0]),
355 static_cast<int>(is_periodic[1]),
356 static_cast<int>(is_periodic[2]))}};
357 }
359 [[nodiscard]] int period (int dir) const noexcept { BL_ASSERT(is_periodic[dir]); return domain.length(dir); }
360
361 [[nodiscard]] Periodicity periodicity () const noexcept {
362 return Periodicity(IntVect(AMREX_D_DECL(domain.length(0) * is_periodic[0],
363 domain.length(1) * is_periodic[1],
364 domain.length(2) * is_periodic[2])));
365 }
366
367 [[nodiscard]] Periodicity periodicity (const Box& b) const noexcept {
368 AMREX_ASSERT(b.cellCentered());
369 return Periodicity(IntVect(AMREX_D_DECL(b.length(0) * is_periodic[0],
370 b.length(1) * is_periodic[1],
371 b.length(2) * is_periodic[2])));
372 }
373
382 void periodicShift (const Box& target,
383 const Box& src,
384 Vector<IntVect>& out) const noexcept;
385
387 [[nodiscard]] Box growNonPeriodicDomain (IntVect const& ngrow) const noexcept;
389 [[nodiscard]] Box growNonPeriodicDomain (int ngrow) const noexcept;
391 [[nodiscard]] Box growPeriodicDomain (IntVect const& ngrow) const noexcept;
393 [[nodiscard]] Box growPeriodicDomain (int ngrow) const noexcept;
394
400 Array<int,AMREX_SPACEDIM> r{{AMREX_D_DECL(is_periodic[0],
401 is_periodic[1],
402 is_periodic[2])}};
403 AMREX_D_TERM(is_periodic[0] = period[0];,
404 is_periodic[1] = period[1];,
405 is_periodic[2] = period[2];);
406 return r;
407 }
408
415 void coarsen (IntVect const& rr) {
416 domain.coarsen(rr);
417 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
418 dx[i] = (ProbHi(i) - ProbLo(i)) / static_cast<Real>(domain.length(i));
419 inv_dx[i] = Real(1.)/dx[i];
420 }
421 }
422
429 void refine (IntVect const& rr) {
430 domain.refine(rr);
431 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
432 dx[i] = (ProbHi(i) - ProbLo(i)) / static_cast<Real>(domain.length(i));
433 inv_dx[i] = Real(1.)/dx[i];
434 }
435 }
436
444
452
454 void computeRoundoffDomain ();
455
457 [[nodiscard]] GpuArray<ParticleReal, AMREX_SPACEDIM> const&
458 RoundOffLo () const { return roundoff_lo; }
459
461 [[nodiscard]] GpuArray<ParticleReal, AMREX_SPACEDIM> const&
462 RoundOffHi () const { return roundoff_hi; }
463
464private:
465 void read_params ();
466
467 // is_periodic and RealBox used to be static
468 bool is_periodic[AMREX_SPACEDIM] = {AMREX_D_DECL(false,false,false)};
469 RealBox prob_domain;
470
471 // Due to round-off errors, not all floating point numbers for which plo >= x < phi
472 // will map to a cell that is inside "domain". "roundoff_{lo,hi}" store
473 // positions that are very close to that in prob_domain, and for which all doubles and floats
474 // between them will map to a cell inside domain.
475 GpuArray<ParticleReal , AMREX_SPACEDIM> roundoff_lo, roundoff_hi;
476
477 //
478 Box domain;
479
480 friend std::istream& operator>> (std::istream&, Geometry&);
481};
482
483
485std::ostream& operator<< (std::ostream&, const Geometry&);
487std::istream& operator>> (std::istream&, Geometry&);
488
495inline
496Geometry
497coarsen (Geometry const& fine, IntVect const& rr) {
498 Geometry r{fine};
499 r.coarsen(rr);
500 return r;
501}
502
509inline
510Geometry
511coarsen (Geometry const& fine, int rr) { return coarsen(fine, IntVect(rr)); }
512
519inline
520Geometry
521refine (Geometry const& crse, IntVect const& rr) {
522 Geometry r{crse};
523 r.refine(rr);
524 return r;
525}
526
533inline
534Geometry
535refine (Geometry const& crse, int rr) { return refine(crse, IntVect(rr)); }
536
537inline
538const Geometry&
540 return *(AMReX::top()->getDefaultGeometry());
541}
542
543}
544
545#endif /*_GEOMETRY_H_*/
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_FORCE_INLINE
Definition AMReX_Extension.H:119
#define AMREX_INLINE
Definition AMReX_Extension.H:127
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Array4< Real const > crse
Definition AMReX_InterpFaceRegister.cpp:92
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
Geometry * getDefaultGeometry() noexcept
Definition AMReX.H:324
static AMReX * top() noexcept
Definition AMReX.H:311
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:564
__host__ __device__ bool cellCentered() const noexcept
Return true if BoxND is cell-centered in all indexing directions.
Definition AMReX_Box.H:327
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:154
__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
Coordinate System.
Definition AMReX_CoordSys.H:24
Real inv_dx[3]
Definition AMReX_CoordSys.H:244
CoordType c_sys
Definition AMReX_CoordSys.H:239
Real dx[3]
Definition AMReX_CoordSys.H:242
@ RZ
Definition AMReX_CoordSys.H:27
@ cartesian
Definition AMReX_CoordSys.H:27
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
A Fortran Array of REALs.
Definition AMReX_FArrayBox.H:231
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:75
Geometry(Geometry &&rhs) noexcept=default
Real ProbLo(int dir) const noexcept
Returns the lo end of the problem domain in specified direction.
Definition AMReX_Geometry.H:188
~Geometry()=default
GpuArray< ParticleReal, 3 > ProbLoArrayInParticleReal() const noexcept
Definition AMReX_Geometry.H:200
void ProbDomain(const RealBox &rb)
Sets the problem domain.
Definition AMReX_Geometry.H:178
bool outsideRoundoffDomain(ParticleReal x, ParticleReal y, ParticleReal z) const
Returns true if a point is outside the roundoff domain. All particles with positions inside the round...
Definition AMReX_Geometry.cpp:708
void GetVolume(MultiFab &vol, const BoxArray &grds, const DistributionMapping &dm, int grow) const
Define a multifab of areas and volumes with given grow factor.
Definition AMReX_Geometry.cpp:207
void refine(IntVect const &rr)
Refine the Geometry by rr.
Definition AMReX_Geometry.H:429
void coarsen(IntVect const &rr)
Coarsen the Geometry by rr.
Definition AMReX_Geometry.H:415
GpuArray< ParticleReal, 3 > ProbHiArrayInParticleReal() const noexcept
Definition AMReX_Geometry.H:204
const Real * ProbHi() const noexcept
Returns the hi end of the problem domain in each dimension.
Definition AMReX_Geometry.H:186
__host__ static __device__ Real Volume(const IntVect &point, const GeometryData &geomdata)
Return the volume of the specified cell.
Definition AMReX_Geometry.H:239
bool isAnyPeriodic() const noexcept
Is domain periodic in any direction?
Definition AMReX_Geometry.H:339
Array< int, 3 > setPeriodicity(Array< int, 3 > const &period) noexcept
Definition AMReX_Geometry.H:399
friend std::istream & operator>>(std::istream &, Geometry &)
Nice ASCII input.
Definition AMReX_Geometry.cpp:25
static void ResetDefaultProbDomain(const RealBox &rb) noexcept
Definition AMReX_Geometry.cpp:183
void GetFaceArea(MultiFab &area, const BoxArray &grds, const DistributionMapping &dm, int dir, int grow) const
Compute area of cell faces in given region and stuff stuff the results into the passed MultiFab.
Definition AMReX_Geometry.cpp:302
Real ProbHi(int dir) const noexcept
Returns the hi end of the problem domain in specified direction.
Definition AMReX_Geometry.H:190
Geometry(const Geometry &rhs)=default
int period(int dir) const noexcept
What's period in specified direction?
Definition AMReX_Geometry.H:359
Geometry & operator=(const Geometry &rhs)=default
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:216
void Domain(const Box &bx)
Sets our rectangular domain.
Definition AMReX_Geometry.H:218
static void Setup(const RealBox *rb=nullptr, int coord=-1, int const *isper=nullptr)
Read static values from ParmParse database.
Definition AMReX_Geometry.cpp:108
GpuArray< int, 3 > isPeriodicArray() const noexcept
Definition AMReX_Geometry.H:353
const RealBox & ProbDomain() const noexcept
Returns the problem domain.
Definition AMReX_Geometry.H:176
Periodicity periodicity() const noexcept
Definition AMReX_Geometry.H:361
Array< int, 3 > isPeriodic() const noexcept
Definition AMReX_Geometry.H:348
void GetDLogA(MultiFab &dloga, const BoxArray &grds, const DistributionMapping &dm, int dir, int grow) const
Compute d(log(A))/dr at cell centers in given region and stuff the results into the passed MultiFab.
Box growNonPeriodicDomain(IntVect const &ngrow) const noexcept
Return domain box with non-periodic directions grown by ngrow.
Definition AMReX_Geometry.cpp:485
void periodicShift(const Box &target, const Box &src, Vector< IntVect > &out) const noexcept
Compute Array of shifts which will translate src so that it will intersect target with non-zero inter...
Definition AMReX_Geometry.cpp:375
Real ProbLength(int dir) const noexcept
Returns length of problem domain in specified dimension.
Definition AMReX_Geometry.H:214
void computeRoundoffDomain()
Compute the roundoff domain.
Definition AMReX_Geometry.cpp:521
const Real * ProbLo() const noexcept
Returns the lo end of the problem domain in each dimension.
Definition AMReX_Geometry.H:184
GeometryData data() const noexcept
Returns non-static copy of geometry's stored data.
Definition AMReX_Geometry.H:125
Geometry() noexcept
The default constructor.
Definition AMReX_Geometry.cpp:45
GpuArray< Real, 3 > ProbHiArray() const noexcept
Definition AMReX_Geometry.H:196
Box growPeriodicDomain(IntVect const &ngrow) const noexcept
Return domain box with periodic directions grown by ngrow.
Definition AMReX_Geometry.cpp:497
Real ProbSize() const noexcept
Returns the overall size of the domain by multiplying the ProbLength's together.
Definition AMReX_Geometry.H:209
bool insideRoundoffDomain(ParticleReal x, ParticleReal y, ParticleReal z) const
Returns true if a point is inside the roundoff domain. All particles with positions inside the roundo...
Definition AMReX_Geometry.cpp:720
GpuArray< Real, 3 > ProbLoArray() const noexcept
Definition AMReX_Geometry.H:192
bool isAllPeriodic() const noexcept
Is domain periodic in all directions?
Definition AMReX_Geometry.H:344
bool isPeriodic(int dir) const noexcept
Is the domain periodic in the specified direction?
Definition AMReX_Geometry.H:337
static void ResetDefaultPeriodicity(const Array< int, 3 > &is_per) noexcept
Definition AMReX_Geometry.cpp:192
GpuArray< ParticleReal, 3 > const & RoundOffLo() const
Returns roundoff domain's lower end.
Definition AMReX_Geometry.H:458
static void ResetDefaultCoord(int coord) noexcept
Definition AMReX_Geometry.cpp:199
void define(const Box &dom, const RealBox *rb=nullptr, int coord=-1, int const *is_per=nullptr)
Definition AMReX_Geometry.cpp:70
GpuArray< ParticleReal, 3 > const & RoundOffHi() const
Returns roundoff domain's higher end.
Definition AMReX_Geometry.H:462
Periodicity periodicity(const Box &b) const noexcept
Definition AMReX_Geometry.H:367
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
This provides length of period for periodic domains. 0 means it is not periodic in that direction....
Definition AMReX_Periodicity.H:17
A Box with real dimensions.
Definition AMReX_RealBox.H:28
__host__ __device__ const Real * lo() const &noexcept
Returns lo side.
Definition AMReX_RealBox.H:53
__host__ __device__ const Real * hi() const &noexcept
Returns hide side.
Definition AMReX_RealBox.H:58
__host__ __device__ Real length(int dir) const noexcept
Returns length in specified direction.
Definition AMReX_RealBox.H:69
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
amrex_particle_real ParticleReal
Floating Point Type for Particles.
Definition AMReX_REAL.H:90
__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
std::array< T, N > Array
Definition AMReX_Array.H:26
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__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
const Geometry & DefaultGeometry()
Definition AMReX_Geometry.H:539
std::istream & operator>>(std::istream &is, BoxND< dim > &bx)
Read from istream.
Definition AMReX_Box.H:1825
Definition AMReX_Geometry.H:26
Box domain
Definition AMReX_Geometry.H:57
__host__ __device__ Real ProbLo(int dir) const noexcept
Definition AMReX_Geometry.H:38
__host__ __device__ Real CellSize(int dir) const noexcept
Definition AMReX_Geometry.H:32
__host__ __device__ int Coord() const noexcept
Coordinates type.
Definition AMReX_Geometry.H:53
__host__ __device__ const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:47
Real dx[3]
Definition AMReX_Geometry.H:58
__host__ __device__ Real ProbHi(int dir) const noexcept
Definition AMReX_Geometry.H:44
__host__ __device__ const Real * CellSize() const noexcept
Returns the cellsize for each coordinate direction.
Definition AMReX_Geometry.H:29
int coord
Definition AMReX_Geometry.H:60
__host__ __device__ int isPeriodic(const int i) const noexcept
Returns whether the domain is periodic in the given direction.
Definition AMReX_Geometry.H:50
int is_periodic[3]
Definition AMReX_Geometry.H:59
__host__ __device__ const Real * ProbLo() const noexcept
Returns the lo end of the problem domain in each dimension.
Definition AMReX_Geometry.H:35
RealBox prob_domain
Definition AMReX_Geometry.H:56
__host__ __device__ const Real * ProbHi() const noexcept
Returns the hi end of the problem domain in each dimension.
Definition AMReX_Geometry.H:41
Fixed-size array that can be used on GPU.
Definition AMReX_Array.H:43