Block-Structured AMR Software Framework
amrex::Geometry Class Reference

Rectangular problem domain geometry. More...

#include <AMReX_Geometry.H>

Inheritance diagram for amrex::Geometry:
amrex::CoordSys

Public Member Functions

 Geometry () noexcept
 The default constructor. More...
 
 Geometry (const Box &dom, const RealBox *rb=nullptr, int coord=-1, int const *is_per=nullptr) noexcept
 
 Geometry (const Box &dom, const RealBox &rb, int coord, Array< int, AMREX_SPACEDIM > const &is_per) noexcept
 
 ~Geometry ()=default
 
 Geometry (const Geometry &rhs)=default
 
 Geometry (Geometry &&rhs) noexcept=default
 
Geometryoperator= (const Geometry &rhs)=default
 
Geometryoperator= (Geometry &&rhs) noexcept=default
 
GeometryData data () const noexcept
 Returns non-static copy of geometry's stored data. More...
 
void define (const Box &dom, const RealBox *rb=nullptr, int coord=-1, int const *is_per=nullptr) noexcept
 
void define (const Box &dom, const RealBox &rb, int coord, Array< int, AMREX_SPACEDIM > const &is_per) noexcept
 
const RealBoxProbDomain () const noexcept
 Returns the problem domain. More...
 
void ProbDomain (const RealBox &rb) noexcept
 Sets the problem domain. More...
 
const Real * ProbLo () const noexcept
 Returns the lo end of the problem domain in each dimension. More...
 
const Real * ProbHi () const noexcept
 Returns the hi end of the problem domain in each dimension. More...
 
Real ProbLo (int dir) const noexcept
 Returns the lo end of the problem domain in specified direction. More...
 
Real ProbHi (int dir) const noexcept
 Returns the hi end of the problem domain in specified direction. More...
 
GpuArray< Real, AMREX_SPACEDIM > ProbLoArray () const noexcept
 
GpuArray< Real, AMREX_SPACEDIM > ProbHiArray () const noexcept
 
GpuArray< ParticleReal, AMREX_SPACEDIM > ProbLoArrayInParticleReal () const noexcept
 
GpuArray< ParticleReal, AMREX_SPACEDIM > ProbHiArrayInParticleReal () const noexcept
 
Real ProbSize () const noexcept
 Returns the overall size of the domain by multiplying the ProbLength's together. More...
 
Real ProbLength (int dir) const noexcept
 Returns length of problem domain in specified dimension. More...
 
const BoxDomain () const noexcept
 Returns our rectangular domain. More...
 
void Domain (const Box &bx) noexcept
 Sets our rectangular domain. More...
 
void GetVolume (MultiFab &vol, const BoxArray &grds, const DistributionMapping &dm, int grow) const
 Define a multifab of areas and volumes with given grow factor. More...
 
void GetVolume (MultiFab &vol) const
 Fill the pre-built multifab with volume. More...
 
void GetVolume (FArrayBox &vol, const BoxArray &grds, int idx, int grow) const
 
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. More...
 
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. More...
 
void GetFaceArea (MultiFab &area, int dir) const
 Fill the pre-built multifab with area. More...
 
void GetFaceArea (FArrayBox &area, const BoxArray &grds, int idx, int dir, int grow) const
 
bool isPeriodic (int dir) const noexcept
 Is the domain periodic in the specified direction? More...
 
bool isAnyPeriodic () const noexcept
 Is domain periodic in any direction? More...
 
bool isAllPeriodic () const noexcept
 Is domain periodic in all directions? More...
 
Array< int, AMREX_SPACEDIM > isPeriodic () const noexcept
 
GpuArray< int, AMREX_SPACEDIM > isPeriodicArray () const noexcept
 
int period (int dir) const noexcept
 What's period in specified direction? More...
 
Periodicity periodicity () const noexcept
 
Periodicity periodicity (const Box &b) const noexcept
 
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 intersection. the array will be resized internally, so anything previously there will be gone DO NOT return non-periodic shifts, even if the box's do intersect without shifting. The logic is that you will only do this as a special case if there is some periodicity. More...
 
Box growNonPeriodicDomain (IntVect const &ngrow) const noexcept
 Return domain box with non-periodic directions grown by ngrow. More...
 
Box growNonPeriodicDomain (int ngrow) const noexcept
 Return domain box with non-periodic directions grown by ngrow. More...
 
Box growPeriodicDomain (IntVect const &ngrow) const noexcept
 Return domain box with periodic directions grown by ngrow. More...
 
Box growPeriodicDomain (int ngrow) const noexcept
 Return domain box with periodic directions grown by ngrow. More...
 
Array< int, AMREX_SPACEDIM > setPeriodicity (Array< int, AMREX_SPACEDIM > const &period) noexcept
 
void coarsen (IntVect const &rr)
 
void refine (IntVect const &rr)
 
bool outsideRoundoffDomain (AMREX_D_DECL(ParticleReal x, ParticleReal y, ParticleReal z)) const
 Returns true if a point is outside the roundoff domain. All particles with positions inside the roundoff domain are sure to be mapped to cells inside the Domain() box. Note that the same need not be true for all points inside ProbDomain(). More...
 
bool insideRoundoffDomain (AMREX_D_DECL(ParticleReal x, ParticleReal y, ParticleReal z)) const
 Returns true if a point is inside the roundoff domain. All particles with positions inside the roundoff domain are sure to be mapped to cells inside the Domain() box. Note that the same need not be true for all points inside ProbDomain(). More...
 
void computeRoundoffDomain ()
 Compute the roundoff domain. Public because it contains an extended host / device lambda. More...
 
- Public Member Functions inherited from amrex::CoordSys
bool Ok () const noexcept
 Is ok? More...
 
void SetCoord (CoordType coord) noexcept
 Set the CoordType. More...
 
CoordType Coord () const noexcept
 Returns the CoordType. More...
 
int CoordInt () const noexcept
 Returns the CoordType as an int. More...
 
bool IsSPHERICAL () const noexcept
 Is CoordType == SPHERICAL? More...
 
bool IsRZ () const noexcept
 Is CoordType == RZ? More...
 
bool IsCartesian () const noexcept
 Is CoordType == cartesian? More...
 
void SetOffset (const Real *x_lo) noexcept
 Sets the offset for each coordinate direction. More...
 
const Real * Offset () const noexcept
 Returns the offset. More...
 
Real Offset (int dir) const noexcept
 Returns the offset for the specified coordinate direction. More...
 
const Real * CellSize () const noexcept
 Returns the cellsize for each coordinate direction. More...
 
Real CellSize (int dir) const noexcept
 Returns the cellsize for the specified coordinate direction. More...
 
GpuArray< Real, AMREX_SPACEDIM > CellSizeArray () const noexcept
 
const Real * InvCellSize () const noexcept
 Returns the inverse cellsize for each coordinate direction. More...
 
Real InvCellSize (int dir) const noexcept
 Returns the inverse cellsize for the specified coordinate direction. More...
 
GpuArray< Real, AMREX_SPACEDIM > InvCellSizeArray () const noexcept
 
Real CellCenter (int point, int dir) const noexcept
 Returns location of cell center in specified direction. More...
 
void CellCenter (const IntVect &point, Vector< Real > &loc) const noexcept
 Return location of cell center. More...
 
void CellCenter (const IntVect &point, Real *loc) const noexcept
 Return location of cell center. More...
 
Real LoEdge (int point, int dir) const noexcept
 Returns location of lo edge in specified direction. More...
 
Real LoEdge (const IntVect &point, int dir) const noexcept
 Equivalent to LoEdge(point[dir], dir). More...
 
Real HiEdge (int point, int dir) const noexcept
 Returns location of hi edge in specified direction. More...
 
Real HiEdge (const IntVect &point, int dir) const noexcept
 Equivalent to HiEdge(point[dir], dir). More...
 
void LoFace (const IntVect &point, int dir, Vector< Real > &loc) const noexcept
 Sets location of lo face into loc. More...
 
void LoFace (const IntVect &point, int dir, Real *loc) const noexcept
 Sets location of lo face into loc. More...
 
void HiFace (const IntVect &point, int dir, Vector< Real > &loc) const noexcept
 Sets location of hi face into loc. More...
 
void HiFace (const IntVect &point, int dir, Real *loc) const noexcept
 Sets location of hi face into loc. More...
 
void LoNode (const IntVect &point, Vector< Real > &loc) const noexcept
 Return location of lower left hand corner. More...
 
void LoNode (const IntVect &point, Real *loc) const noexcept
 Return location of lower left hand corner. More...
 
void HiNode (const IntVect &point, Vector< Real > &loc) const noexcept
 Return location of upper right hand corner. More...
 
void HiNode (const IntVect &point, Real *loc) const noexcept
 Return location of upper right hand corner. More...
 
IntVect CellIndex (const Real *point) const noexcept
 Returns cell centered index of cell containing point. This may return undesired results if point is on a cell boundary. More...
 
IntVect LowerIndex (const Real *point) const noexcept
 Returns node centered index of lower left hand corner of cell containing this point. More...
 
IntVect UpperIndex (const Real *point) const noexcept
 Returns node centered index of upper right hand corner of cell containing this point. More...
 
void SetVolume (FArrayBox &a_volfab, const Box &region) const
 Compute cell volumes in given region and place them into input FAB. More...
 
void GetVolume (FArrayBox &vol, const Box &region) const
 Compute cell volumes in given region and place them into resize()d input FAB. More...
 
void SetDLogA (FArrayBox &a_dlogafab, const Box &region, int dir) const
 Compute d(log(A))/dr at cell centers in given region and place them into input FAB. More...
 
void GetDLogA (FArrayBox &dloga, const Box &region, int dir) const
 Compute d(log(A))/dr at cell centers in given region and return the results in the resize()d input FAB. More...
 
Real Volume (const IntVect &point) const
 Return the volume of the specified cell. More...
 
Real Volume (const Real xlo[AMREX_SPACEDIM], const Real xhi[AMREX_SPACEDIM]) const
 Return the volume of the specified cell. More...
 
void SetFaceArea (FArrayBox &a_areafab, const Box &region, int dir) const
 Compute area of cell faces in given region and given index direction and return the result in input FAB. More...
 
void GetFaceArea (FArrayBox &area, const Box &region, int dir) const
 Compute area of cell faces in given region and given index direction and return the result in resize()d input FAB. More...
 
Real AreaLo (const IntVect &point, int dir) const noexcept
 Returns lo face area of given cell in direction dir. More...
 
Real AreaHi (const IntVect &point, int dir) const noexcept
 Returns hi face area of given cell in direction dir. More...
 
void GetEdgeLoc (Vector< Real > &loc, const Box &region, int dir) const
 Return array of physical locations of cell edges in the resize()d input array. More...
 
void GetCellLoc (Vector< Real > &loc, const Box &region, int dir) const
 Return array of physical locations of cell centers in the resize()d input array. More...
 
void GetEdgeVolCoord (Vector< Real > &vc, const Box &region, int dir) const
 Return array of volume coordinates at cell edges for region in given direction. More...
 
void GetCellVolCoord (Vector< Real > &vc, const Box &region, int dir) const
 Return array of volume coordinates at cell centers for region in given direction. More...
 

Static Public Member Functions

static void Setup (const RealBox *rb=nullptr, int coord=-1, int const *isper=nullptr) noexcept
 Read static values from ParmParse database. More...
 
static void ResetDefaultProbDomain (const RealBox &rb) noexcept
 
static void ResetDefaultPeriodicity (const Array< int, AMREX_SPACEDIM > &is_per) noexcept
 
static void ResetDefaultCoord (int coord) noexcept
 
AMREX_GPU_HOST_DEVICE static AMREX_INLINE Real Volume (const IntVect &point, const GeometryData &geomdata)
 Return the volume of the specified cell. More...
 

Private Member Functions

void read_params ()
 

Private Attributes

bool is_periodic [AMREX_SPACEDIM] = {AMREX_D_DECL(false,false,false)}
 
RealBox prob_domain
 
GpuArray< ParticleReal, AMREX_SPACEDIM > roundoff_lo
 
GpuArray< ParticleReal, AMREX_SPACEDIM > roundoff_hi
 
Box domain
 

Additional Inherited Members

- Public Types inherited from amrex::CoordSys
enum  CoordType { undef = -1 , cartesian = 0 , RZ = 1 , SPHERICAL = 2 }
 
- Protected Attributes inherited from amrex::CoordSys
CoordType c_sys = undef
 
Real offset [AMREX_SPACEDIM]
 
Real dx [AMREX_SPACEDIM] = {AMREX_D_DECL(0.,0.,0.)}
 
Real inv_dx [AMREX_SPACEDIM]
 
bool ok = false
 

Detailed Description

Rectangular problem domain geometry.

This class describes problem domain and coordinate system for RECTANGULAR problem domains. Since the problem domain is RECTANGULAR, periodicity is meaningful.

Constructor & Destructor Documentation

◆ Geometry() [1/5]

amrex::Geometry::Geometry ( )
noexcept

The default constructor.

Leaves object in an unusable state. A "define" method must be called before use.

◆ Geometry() [2/5]

amrex::Geometry::Geometry ( const Box dom,
const RealBox rb = nullptr,
int  coord = -1,
int const *  is_per = nullptr 
)
explicitnoexcept

Constructs a fully-defined geometry object that describes the problem domain and coordinate system.

Parameters
domcell-centered Box specifying the index space bounds of the domain.
rbA pointer to a RealBox specifying the real space bounds of the domain. Optional. If not specified, will be ParmParsed from "geometry.prob_lo" and "geometry.prob_hi".
coordSpecifies the coordinate system type. Optional. Can be one of: CoordSys::cartesian CoordSys::RZ CoordSys::SPHERICAL Default is cartesian.
is_perpointer to memory space specifying periodicity in each coordinate direction. Optional. Default is non-periodic.

◆ Geometry() [3/5]

amrex::Geometry::Geometry ( const Box dom,
const RealBox rb,
int  coord,
Array< int, AMREX_SPACEDIM > const &  is_per 
)
noexcept

Constructs a fully-defined geometry object that describes the problem domain and coordinate system.

Parameters
domcell-centered Box specifying the index space bounds of the domain.
rbA pointer to a RealBox specifying the real space bounds of the domain. Optional. If not specified, will be ParmParsed from "geometry.prob_lo" and "geometry.prob_hi".
coordSpecifies the coordinate system type. Optional. Can be one of: CoordSys::cartesian CoordSys::RZ CoordSys::SPHERICAL Default is cartesian.
is_peramrex::Array specifying periodicity in each coordinate direction. Optional. Default is non-periodic.

◆ ~Geometry()

amrex::Geometry::~Geometry ( )
default

◆ Geometry() [4/5]

amrex::Geometry::Geometry ( const Geometry rhs)
default

◆ Geometry() [5/5]

amrex::Geometry::Geometry ( Geometry &&  rhs)
defaultnoexcept

Member Function Documentation

◆ coarsen()

void amrex::Geometry::coarsen ( IntVect const &  rr)
inline

◆ computeRoundoffDomain()

void amrex::Geometry::computeRoundoffDomain ( )

Compute the roundoff domain. Public because it contains an extended host / device lambda.

◆ data()

GeometryData amrex::Geometry::data ( ) const
inlinenoexcept

Returns non-static copy of geometry's stored data.

◆ define() [1/2]

void amrex::Geometry::define ( const Box dom,
const RealBox rb,
int  coord,
Array< int, AMREX_SPACEDIM > const &  is_per 
)
noexcept

Defines a geometry object that describes the problem domain and coordinate system for rectangular problem domains. Useful for when the Geometry object has been constructed with the default constructor.

Parameters
domcell-centered Box specifying the index space bounds of the domain.
rbA pointer to a RealBox specifying the real space bounds of the domain. Optional. If not specified, will be ParmParsed from "geometry.prob_lo" and "geometry.prob_hi".
coordSpecifies the coordinate system type. Optional. Can be one of: CoordSys::cartesian CoordSys::RZ CoordSys::SPHERICAL Default is cartesian.
is_peramrex::Array specifying periodicity in each coordinate direction. Optional. Default is non-periodic.

◆ define() [2/2]

void amrex::Geometry::define ( const Box dom,
const RealBox rb = nullptr,
int  coord = -1,
int const *  is_per = nullptr 
)
noexcept

Defines a geometry object that describes the problem domain and coordinate system for rectangular problem domains. Useful for when the Geometry object has been constructed with the default constructor.

Parameters
domcell-centered Box specifying the index space bounds of the domain.
rbA pointer to a RealBox specifying the real space bounds of the domain. Optional. If not specified, will be ParmParsed from "geometry.prob_lo" and "geometry.prob_hi".
coordSpecifies the coordinate system type. Optional. Can be one of: CoordSys::cartesian CoordSys::RZ CoordSys::SPHERICAL Default is cartesian.
is_perpointer to memory space specifying periodicity in each coordinate direction. Optional. Default is non-periodic.

◆ Domain() [1/2]

const Box& amrex::Geometry::Domain ( ) const
inlinenoexcept

Returns our rectangular domain.

◆ Domain() [2/2]

void amrex::Geometry::Domain ( const Box bx)
inlinenoexcept

Sets our rectangular domain.

◆ GetDLogA()

void amrex::Geometry::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.

◆ GetFaceArea() [1/3]

void amrex::Geometry::GetFaceArea ( FArrayBox area,
const BoxArray grds,
int  idx,
int  dir,
int  grow 
) const

◆ GetFaceArea() [2/3]

void amrex::Geometry::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.

◆ GetFaceArea() [3/3]

void amrex::Geometry::GetFaceArea ( MultiFab area,
int  dir 
) const

Fill the pre-built multifab with area.

◆ GetVolume() [1/3]

void amrex::Geometry::GetVolume ( FArrayBox vol,
const BoxArray grds,
int  idx,
int  grow 
) const

◆ GetVolume() [2/3]

void amrex::Geometry::GetVolume ( MultiFab vol) const

Fill the pre-built multifab with volume.

◆ GetVolume() [3/3]

void amrex::Geometry::GetVolume ( MultiFab vol,
const BoxArray grds,
const DistributionMapping dm,
int  grow 
) const

Define a multifab of areas and volumes with given grow factor.

◆ growNonPeriodicDomain() [1/2]

Box amrex::Geometry::growNonPeriodicDomain ( int  ngrow) const
noexcept

Return domain box with non-periodic directions grown by ngrow.

◆ growNonPeriodicDomain() [2/2]

Box amrex::Geometry::growNonPeriodicDomain ( IntVect const &  ngrow) const
noexcept

Return domain box with non-periodic directions grown by ngrow.

◆ growPeriodicDomain() [1/2]

Box amrex::Geometry::growPeriodicDomain ( int  ngrow) const
noexcept

Return domain box with periodic directions grown by ngrow.

◆ growPeriodicDomain() [2/2]

Box amrex::Geometry::growPeriodicDomain ( IntVect const &  ngrow) const
noexcept

Return domain box with periodic directions grown by ngrow.

◆ insideRoundoffDomain()

bool amrex::Geometry::insideRoundoffDomain ( AMREX_D_DECL(ParticleReal x, ParticleReal y, ParticleReal z)  ) const

Returns true if a point is inside the roundoff domain. All particles with positions inside the roundoff domain are sure to be mapped to cells inside the Domain() box. Note that the same need not be true for all points inside ProbDomain().

◆ isAllPeriodic()

bool amrex::Geometry::isAllPeriodic ( ) const
inlinenoexcept

Is domain periodic in all directions?

◆ isAnyPeriodic()

bool amrex::Geometry::isAnyPeriodic ( ) const
inlinenoexcept

Is domain periodic in any direction?

◆ isPeriodic() [1/2]

Array<int,AMREX_SPACEDIM> amrex::Geometry::isPeriodic ( ) const
inlinenoexcept

◆ isPeriodic() [2/2]

bool amrex::Geometry::isPeriodic ( int  dir) const
inlinenoexcept

Is the domain periodic in the specified direction?

◆ isPeriodicArray()

GpuArray<int,AMREX_SPACEDIM> amrex::Geometry::isPeriodicArray ( ) const
inlinenoexcept

◆ operator=() [1/2]

Geometry& amrex::Geometry::operator= ( const Geometry rhs)
default

◆ operator=() [2/2]

Geometry& amrex::Geometry::operator= ( Geometry &&  rhs)
defaultnoexcept

◆ outsideRoundoffDomain()

bool amrex::Geometry::outsideRoundoffDomain ( AMREX_D_DECL(ParticleReal x, ParticleReal y, ParticleReal z)  ) const

Returns true if a point is outside the roundoff domain. All particles with positions inside the roundoff domain are sure to be mapped to cells inside the Domain() box. Note that the same need not be true for all points inside ProbDomain().

◆ period()

int amrex::Geometry::period ( int  dir) const
inlinenoexcept

What's period in specified direction?

◆ periodicity() [1/2]

Periodicity amrex::Geometry::periodicity ( ) const
inlinenoexcept

◆ periodicity() [2/2]

Periodicity amrex::Geometry::periodicity ( const Box b) const
inlinenoexcept

◆ periodicShift()

void amrex::Geometry::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 intersection. the array will be resized internally, so anything previously there will be gone DO NOT return non-periodic shifts, even if the box's do intersect without shifting. The logic is that you will only do this as a special case if there is some periodicity.

◆ ProbDomain() [1/2]

const RealBox& amrex::Geometry::ProbDomain ( ) const
inlinenoexcept

Returns the problem domain.

◆ ProbDomain() [2/2]

void amrex::Geometry::ProbDomain ( const RealBox rb)
inlinenoexcept

Sets the problem domain.

◆ ProbHi() [1/2]

const Real* amrex::Geometry::ProbHi ( ) const
inlinenoexcept

Returns the hi end of the problem domain in each dimension.

◆ ProbHi() [2/2]

Real amrex::Geometry::ProbHi ( int  dir) const
inlinenoexcept

Returns the hi end of the problem domain in specified direction.

◆ ProbHiArray()

GpuArray<Real,AMREX_SPACEDIM> amrex::Geometry::ProbHiArray ( ) const
inlinenoexcept

◆ ProbHiArrayInParticleReal()

GpuArray<ParticleReal,AMREX_SPACEDIM> amrex::Geometry::ProbHiArrayInParticleReal ( ) const
inlinenoexcept

◆ ProbLength()

Real amrex::Geometry::ProbLength ( int  dir) const
inlinenoexcept

Returns length of problem domain in specified dimension.

◆ ProbLo() [1/2]

const Real* amrex::Geometry::ProbLo ( ) const
inlinenoexcept

Returns the lo end of the problem domain in each dimension.

◆ ProbLo() [2/2]

Real amrex::Geometry::ProbLo ( int  dir) const
inlinenoexcept

Returns the lo end of the problem domain in specified direction.

◆ ProbLoArray()

GpuArray<Real,AMREX_SPACEDIM> amrex::Geometry::ProbLoArray ( ) const
inlinenoexcept

◆ ProbLoArrayInParticleReal()

GpuArray<ParticleReal,AMREX_SPACEDIM> amrex::Geometry::ProbLoArrayInParticleReal ( ) const
inlinenoexcept

◆ ProbSize()

Real amrex::Geometry::ProbSize ( ) const
inlinenoexcept

Returns the overall size of the domain by multiplying the ProbLength's together.

◆ read_params()

void amrex::Geometry::read_params ( )
private

◆ refine()

void amrex::Geometry::refine ( IntVect const &  rr)
inline

◆ ResetDefaultCoord()

void amrex::Geometry::ResetDefaultCoord ( int  coord)
staticnoexcept

◆ ResetDefaultPeriodicity()

void amrex::Geometry::ResetDefaultPeriodicity ( const Array< int, AMREX_SPACEDIM > &  is_per)
staticnoexcept

◆ ResetDefaultProbDomain()

void amrex::Geometry::ResetDefaultProbDomain ( const RealBox rb)
staticnoexcept

◆ setPeriodicity()

Array<int,AMREX_SPACEDIM> amrex::Geometry::setPeriodicity ( Array< int, AMREX_SPACEDIM > const &  period)
inlinenoexcept

Set periodicity flags and return the old flags. Note that, unlike Periodicity class, the flags are just boolean.

◆ Setup()

void amrex::Geometry::Setup ( const RealBox rb = nullptr,
int  coord = -1,
int const *  isper = nullptr 
)
staticnoexcept

Read static values from ParmParse database.

◆ Volume()

AMREX_GPU_HOST_DEVICE static AMREX_INLINE Real amrex::Geometry::Volume ( const IntVect point,
const GeometryData geomdata 
)
inlinestatic

Return the volume of the specified cell.

Member Data Documentation

◆ domain

Box amrex::Geometry::domain
private

◆ is_periodic

bool amrex::Geometry::is_periodic[AMREX_SPACEDIM] = {AMREX_D_DECL(false,false,false)}
private

◆ prob_domain

RealBox amrex::Geometry::prob_domain
private

◆ roundoff_hi

GpuArray<ParticleReal , AMREX_SPACEDIM> amrex::Geometry::roundoff_hi
private

◆ roundoff_lo

GpuArray<ParticleReal , AMREX_SPACEDIM> amrex::Geometry::roundoff_lo
private

The documentation for this class was generated from the following files: