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

A Rectangular Domain on an Integer Lattice. More...

#include <AMReX_Box.H>

Public Member Functions

constexpr AMREX_GPU_HOST_DEVICE Box () noexcept
 
constexpr AMREX_GPU_HOST_DEVICE Box (const IntVect &small, const IntVect &big) noexcept
 Construct cell-centered type Box. More...
 
AMREX_GPU_HOST_DEVICE Box (const IntVect &small, const int *vec_len) noexcept
 Construct box with specified lengths. More...
 
AMREX_GPU_HOST_DEVICE Box (const IntVect &small, const IntVect &big, const IntVect &typ) noexcept
 Construct Box with given type. small and big are expected to be consistent with given type. More...
 
AMREX_GPU_HOST_DEVICE Box (const IntVect &small, const IntVect &big, IndexType t) noexcept
 Construct dimension specific Boxes. More...
 
template<typename T >
AMREX_GPU_HOST_DEVICE Box (Array4< T > const &a) noexcept
 
AMREX_GPU_HOST_DEVICE const IntVectsmallEnd () const &noexcept
 Get the smallend of the box. More...
 
const IntVectsmallEnd () &&=delete
 Get the smallend of the box. More...
 
AMREX_GPU_HOST_DEVICE int smallEnd (int dir) const &noexcept
 Returns the coordinate of the low end in the given direction. More...
 
AMREX_GPU_HOST_DEVICE const IntVectbigEnd () const &noexcept
 Get the bigend. More...
 
const IntVectbigEnd () &&=delete
 Get the bigend. More...
 
AMREX_GPU_HOST_DEVICE int bigEnd (int dir) const noexcept
 Returns the coordinate of the high end in the given direction. More...
 
AMREX_GPU_HOST_DEVICE IndexType ixType () const noexcept
 Returns the indexing type. More...
 
AMREX_GPU_HOST_DEVICE IntVect type () const noexcept
 Returns the indexing type. More...
 
AMREX_GPU_HOST_DEVICE IndexType::CellIndex type (int dir) const noexcept
 Returns the indexing type in the specified direction. More...
 
AMREX_GPU_HOST_DEVICE IntVect size () const noexcept
 Return the length of the Box. More...
 
AMREX_GPU_HOST_DEVICE IntVect length () const noexcept
 Return the length of the Box. More...
 
AMREX_GPU_HOST_DEVICE int length (int dir) const noexcept
 Return the length of the Box in given direction. More...
 
AMREX_GPU_HOST_DEVICE GpuArray< int, 3 > length3d () const noexcept
 
AMREX_GPU_HOST_DEVICE GpuArray< int, 3 > loVect3d () const noexcept
 
AMREX_GPU_HOST_DEVICE GpuArray< int, 3 > hiVect3d () const noexcept
 
AMREX_GPU_HOST_DEVICE const intloVect () const &noexcept
 Returns a constant pointer the array of low end coordinates. Useful for calls to FORTRAN. More...
 
AMREX_GPU_HOST_DEVICE const intloVect () &&=delete
 
AMREX_GPU_HOST_DEVICE const inthiVect () const &noexcept
 Returns a constant pointer the array of high end coordinates. Useful for calls to FORTRAN. More...
 
AMREX_GPU_HOST_DEVICE const inthiVect () &&=delete
 
AMREX_GPU_HOST_DEVICE int operator[] (Orientation face) const noexcept
 Returns the coordinate normal to given face. More...
 
AMREX_GPU_HOST_DEVICE bool isEmpty () const noexcept
 Checks if it is an empty box. More...
 
AMREX_GPU_HOST_DEVICE bool ok () const noexcept
 Checks if it is a proper Box (including a valid type). More...
 
AMREX_GPU_HOST_DEVICE bool contains (const IntVect &p) const noexcept
 Returns true if argument is contained within Box. More...
 
AMREX_GPU_HOST_DEVICE bool contains (const Dim3 &p) const noexcept
 Returns true if argument is contained within Box. More...
 
AMREX_GPU_HOST_DEVICE bool contains (int i, int, int) const noexcept
 Returns true if argument is contained within Box. More...
 
AMREX_GPU_HOST_DEVICE bool contains (const Box &b) const noexcept
 Returns true if argument is contained within Box. It is an error if the Boxes have different types. More...
 
AMREX_GPU_HOST_DEVICE bool strictly_contains (const IntVect &p) const noexcept
 Returns true if argument is strictly contained within Box. More...
 
AMREX_GPU_HOST_DEVICE bool strictly_contains (const Box &b) const noexcept
 Returns true if argument is strictly contained within Box. It is an error if the Boxes have different types. More...
 
AMREX_GPU_HOST_DEVICE bool strictly_contains (int i, int, int) const noexcept
 Returns true if argument is strictly contained within Box. More...
 
AMREX_GPU_HOST_DEVICE bool intersects (const Box &b) const noexcept
 Returns true if Boxes have non-null intersections. It is an error if the Boxes have different types. More...
 
AMREX_GPU_HOST_DEVICE bool sameSize (const Box &b) const noexcept
 Returns true is Boxes same size, ie translates of each other,. It is an error if they have different types. More...
 
AMREX_GPU_HOST_DEVICE bool sameType (const Box &b) const noexcept
 Returns true if Boxes have same type. More...
 
AMREX_GPU_HOST_DEVICE bool operator== (const Box &b) const noexcept
 Returns true if Boxes are identical (including type). More...
 
AMREX_GPU_HOST_DEVICE bool operator!= (const Box &b) const noexcept
 Returns true if Boxes differ (including type). More...
 
AMREX_GPU_HOST_DEVICE bool operator< (const Box &rhs) const noexcept
 
AMREX_GPU_HOST_DEVICE bool operator<= (const Box &rhs) const noexcept
 
AMREX_GPU_HOST_DEVICE bool operator> (const Box &rhs) const noexcept
 
AMREX_GPU_HOST_DEVICE bool operator>= (const Box &rhs) const noexcept
 
AMREX_GPU_HOST_DEVICE bool cellCentered () const noexcept
 Returns true if Box is cell-centered in all indexing directions. More...
 
AMREX_GPU_HOST_DEVICE Long numPts () const noexcept
 Returns the number of points contained in the Box. More...
 
AMREX_GPU_HOST_DEVICE double d_numPts () const noexcept
 Returns the number of points contained in the Box. This is intended for use only in diagnostic messages. More...
 
AMREX_GPU_HOST_DEVICE Long volume () const noexcept
 Return the volume, in indexing space, of region enclosed by this Box. This is identical to numPts() for CELL centered Box; otherwise, numPts() > volume(). More...
 
AMREX_GPU_HOST_DEVICE int longside (int &dir) const noexcept
 Returns length of longest side. dir is modified to give direction with longest side: 0...SPACEDIM-1. Ignores type. More...
 
AMREX_GPU_HOST_DEVICE int longside () const noexcept
 Returns length of longest side. Ignores type. More...
 
AMREX_GPU_HOST_DEVICE int shortside (int &dir) const noexcept
 Returns length of shortest side. dir is modified to give direction with shortest side: 0...SPACEDIM-1. Ignores type. More...
 
AMREX_GPU_HOST_DEVICE int shortside () const noexcept
 Returns length of shortest side. Ignores type. More...
 
AMREX_GPU_HOST_DEVICE Long index (const IntVect &v) const noexcept
 Returns offset of point from smallend; i.e. index(smallend) -> 0, bigend would return numPts()-1. Is used in accessing FArrayBox. More...
 
AMREX_GPU_HOST_DEVICE IntVect atOffset (Long offset) const noexcept
 Given the offset, compute IntVect. More...
 
AMREX_GPU_HOST_DEVICE GpuArray< int, 3 > atOffset3d (Long offset) const noexcept
 
AMREX_GPU_HOST_DEVICE BoxsetSmall (const IntVect &sm) noexcept
 Redefine the small end of the Box. More...
 
AMREX_GPU_HOST_DEVICE BoxsetSmall (int dir, int sm_index) noexcept
 Redefine the small end of the Box. More...
 
AMREX_GPU_HOST_DEVICE BoxsetBig (const IntVect &bg) noexcept
 Redefine the big end of the Box. More...
 
AMREX_GPU_HOST_DEVICE BoxsetBig (int dir, int bg_index) noexcept
 Redefine the big end of the Box. More...
 
AMREX_GPU_HOST_DEVICE BoxsetRange (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. NOTE: This will yield an illegal Box if n_cells <= 0. More...
 
AMREX_GPU_HOST_DEVICE BoxsetType (const IndexType &t) noexcept
 Set indexing type. More...
 
AMREX_GPU_HOST_DEVICE Boxshift (int dir, int nzones) noexcept
 Shift this Box nzones indexing positions in coordinate direction dir. More...
 
AMREX_GPU_HOST_DEVICE Boxshift (const IntVect &iv) noexcept
 Equivalent to b.shift(0,iv[0]).shift(1,iv[1]) .... More...
 
AMREX_GPU_HOST_DEVICE BoxshiftHalf (int dir, int num_halfs) noexcept
 This member shifts the Box by "half" indices, thereby converting the Box from type CELL to NODE and visa-versa. b.shiftHalf(0,1) shifts b to the right by 1/2 cells. b.shiftHalf(1,-3) shifts b in the -j direction by 3/2 cells. NOTE: If num_halfs is EVEN the shift is num_halfs/2 full zones and hence will not change the type. This is: b.shifthalf(4) == b.shift(2). More...
 
AMREX_GPU_HOST_DEVICE BoxshiftHalf (const IntVect &iv) noexcept
 Equivalent to b.shiftHalf(0,iv[0]).shiftHalf(1,iv[1]) ... More...
 
AMREX_GPU_HOST_DEVICE Boxconvert (IndexType typ) noexcept
 Convert the Box from the current type into the argument type. This may change the Box coordinates: type CELL -> NODE : increase coordinate by one on high end type NODE -> CELL : reduce coordinate by one on high end other type mappings make no change. More...
 
AMREX_GPU_HOST_DEVICE Boxconvert (const IntVect &typ) noexcept
 Convert the Box from the current type into the argument type. This may change the Box coordinates: type CELL -> NODE : increase coordinate by one on high end type NODE -> CELL : reduce coordinate by one on high end other type mappings make no change. More...
 
AMREX_GPU_HOST_DEVICE BoxsurroundingNodes () noexcept
 Convert to NODE type in all directions. More...
 
AMREX_GPU_HOST_DEVICE BoxsurroundingNodes (int dir) noexcept
 Convert to NODE type in given direction. More...
 
AMREX_GPU_HOST_DEVICE BoxsurroundingNodes (Direction d) noexcept
 
AMREX_GPU_HOST_DEVICE BoxenclosedCells () noexcept
 Convert to CELL type in all directions. More...
 
AMREX_GPU_HOST_DEVICE BoxenclosedCells (int dir) noexcept
 Convert to CELL type in given direction. More...
 
AMREX_GPU_HOST_DEVICE BoxenclosedCells (Direction d) noexcept
 
AMREX_GPU_HOST_DEVICE Box operator& (const Box &rhs) const noexcept
 Return Box that is intersection of this Box and argument. The Boxes MUST be of same type. More...
 
AMREX_GPU_HOST_DEVICE Boxoperator&= (const Box &rhs) noexcept
 Intersect this Box with its argument. The Boxes MUST be of the same type. More...
 
AMREX_GPU_HOST_DEVICE BoxminBox (const Box &b) noexcept
 Modify Box to that of the minimum Box containing both the original Box and the argument. Both Boxes must have identical type. More...
 
AMREX_GPU_HOST_DEVICE Boxoperator+= (const IntVect &v) noexcept
 Shift Box (relative) by given IntVect. More...
 
AMREX_GPU_HOST_DEVICE Box operator+ (const IntVect &v) const noexcept
 Shift Box (relative) by given IntVect. More...
 
AMREX_GPU_HOST_DEVICE Boxoperator-= (const IntVect &v) noexcept
 Shift Box (relative) by given IntVect. More...
 
AMREX_GPU_HOST_DEVICE Box operator- (const IntVect &v) const noexcept
 Shift Box (relative) by given IntVect. More...
 
AMREX_GPU_HOST_DEVICE Box chop (int dir, int chop_pnt) noexcept
 Chop the Box at the chop_pnt in the dir direction returns one Box, modifies the object Box. The union of the two is the original Box. The modified Box is the low end, the returned Box is the high end. If type(dir) = CELL, the Boxes are disjoint with the chop_pnt included in the high end (new Box). It is an ERROR if chop_pnt is the low end of the orig Box. If type(dir) = NODE, the chop_pnt is included in both Boxes but is the only point in common. It is also an error if the chop_pnt is an end node of the Box. More...
 
AMREX_GPU_HOST_DEVICE Boxgrow (int i) noexcept
 
AMREX_GPU_HOST_DEVICE Boxgrow (const IntVect &v) noexcept
 Grow Box in each direction by specified amount. More...
 
AMREX_GPU_HOST_DEVICE Boxgrow (int idir, int n_cell) noexcept
 Grow the Box on the low and high end by n_cell cells in direction idir. More...
 
AMREX_GPU_HOST_DEVICE Boxgrow (Direction d, int n_cell) noexcept
 
AMREX_GPU_HOST_DEVICE BoxgrowLo (int idir, int n_cell=1) noexcept
 Grow the Box on the low end by n_cell cells in direction idir. NOTE: n_cell negative shrinks the Box by that number of cells. More...
 
AMREX_GPU_HOST_DEVICE BoxgrowLo (Direction d, int n_cell=1) noexcept
 
AMREX_GPU_HOST_DEVICE BoxgrowHi (int idir, int n_cell=1) noexcept
 Grow the Box on the high end by n_cell cells in direction idir. NOTE: n_cell negative shrinks the Box by that number of cells. More...
 
AMREX_GPU_HOST_DEVICE BoxgrowHi (Direction d, int n_cell=1) noexcept
 
AMREX_GPU_HOST_DEVICE Boxgrow (Orientation face, int n_cell=1) noexcept
 Grow in the direction of the given face. More...
 
AMREX_GPU_HOST_DEVICE Boxrefine (int ref_ratio) noexcept
 Refine Box by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo*ratio and hi <- (hi+1)*ratio - 1. NOTE: if type(dir) = NODE centered: lo <- lo*ratio and hi <- hi*ratio. More...
 
AMREX_GPU_HOST_DEVICE Boxrefine (const IntVect &ref_ratio) noexcept
 
AMREX_GPU_HOST_DEVICE Boxcoarsen (int ref_ratio) noexcept
 Coarsen Box by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo/ratio and hi <- hi/ratio. NOTE: if type(dir) = NODE centered: lo <- lo/ratio and hi <- hi/ratio + ((hiratio)==0 ? 0 : 1). That is, refinement of coarsened Box must contain the original Box. More...
 
AMREX_GPU_HOST_DEVICE Boxcoarsen (const IntVect &ref_ratio) noexcept
 Coarsen Box by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo/ratio and hi <- hi/ratio. NOTE: if type(dir) = NODE centered: lo <- lo/ratio and hi <- hi/ratio + ((hiratio)==0 ? 0 : 1). That is, refinement of coarsened Box must contain the original Box. More...
 
AMREX_GPU_HOST_DEVICE void next (IntVect &) const noexcept
 Step through the rectangle. It is a runtime error to give a point not inside rectangle. Iteration may not be efficient. More...
 
AMREX_GPU_HOST_DEVICE bool isSquare () const noexcept
 
AMREX_GPU_HOST_DEVICE bool coarsenable (const IntVect &refrat, const IntVect &min_width) const noexcept
 
AMREX_GPU_HOST_DEVICE bool coarsenable (int refrat, int min_width=1) const noexcept
 
AMREX_GPU_HOST_DEVICE bool coarsenable (const IntVect &refrat, int min_width=1) const noexcept
 
AMREX_GPU_HOST_DEVICE void normalize () noexcept
 
AMREX_GPU_HOST_DEVICE BoxmakeSlab (int direction, int slab_index) noexcept
 

Static Public Member Functions

static AMREX_GPU_HOST_DEVICE Box TheUnitBox () noexcept
 This static member function returns a constant reference to an object of type Box representing the unit box in AMREX_SPACEDIM-dimensional space. More...
 

Private Attributes

IntVect smallend
 
IntVect bigend
 
IndexType btype
 

Friends

class BoxCommHelper
 
MPI_Datatype ParallelDescriptor::Mpi_typemap ()
 
AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 lbound (Box const &box) noexcept
 
AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 ubound (Box const &box) noexcept
 
AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 begin (Box const &box) noexcept
 
AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 end (Box const &box) noexcept
 
AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 length (Box const &box) noexcept
 
AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 max_lbound (Box const &, Box const &) noexcept
 
AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 max_lbound (Box const &, Dim3 const &) noexcept
 
AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 min_ubound (Box const &, Box const &) noexcept
 
AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 min_ubound (Box const &, Dim3 const &) noexcept
 
AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Box minBox (Box const &, Box const &, IndexType) noexcept
 

Detailed Description

A Rectangular Domain on an Integer Lattice.

A Box is an abstraction for defining discrete regions of SPACEDIM indexing space. Boxes have an IndexType, which defines IndexType::CELL or IndexType::NODE based points for each direction and a low and high INTVECT which defines the lower and upper corners of the Box. Boxes can exist in positive and negative indexing space.

Box is a dimension dependent class, so SPACEDIM must be defined as either 1, 2, or 3 when compiling.

Constructor & Destructor Documentation

◆ Box() [1/6]

constexpr AMREX_GPU_HOST_DEVICE amrex::Box::Box ( )
inlineconstexprnoexcept

◆ Box() [2/6]

constexpr AMREX_GPU_HOST_DEVICE amrex::Box::Box ( const IntVect small,
const IntVect big 
)
inlineconstexprnoexcept

Construct cell-centered type Box.

◆ Box() [3/6]

AMREX_GPU_HOST_DEVICE amrex::Box::Box ( const IntVect small,
const int vec_len 
)
inlinenoexcept

Construct box with specified lengths.

◆ Box() [4/6]

AMREX_GPU_HOST_DEVICE amrex::Box::Box ( const IntVect small,
const IntVect big,
const IntVect typ 
)
inlinenoexcept

Construct Box with given type. small and big are expected to be consistent with given type.

◆ Box() [5/6]

AMREX_GPU_HOST_DEVICE amrex::Box::Box ( const IntVect small,
const IntVect big,
IndexType  t 
)
inlinenoexcept

Construct dimension specific Boxes.

◆ Box() [6/6]

template<typename T >
AMREX_GPU_HOST_DEVICE amrex::Box::Box ( Array4< T > const &  a)
inlineexplicitnoexcept

Member Function Documentation

◆ atOffset()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE IntVect amrex::Box::atOffset ( Long  offset) const
noexcept

Given the offset, compute IntVect.

◆ atOffset3d()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE GpuArray< int, 3 > amrex::Box::atOffset3d ( Long  offset) const
noexcept

◆ bigEnd() [1/3]

const IntVect& amrex::Box::bigEnd ( ) &&
delete

Get the bigend.

◆ bigEnd() [2/3]

AMREX_GPU_HOST_DEVICE const IntVect& amrex::Box::bigEnd ( ) const &
inlinenoexcept

Get the bigend.

◆ bigEnd() [3/3]

AMREX_GPU_HOST_DEVICE int amrex::Box::bigEnd ( int  dir) const
inlinenoexcept

Returns the coordinate of the high end in the given direction.

◆ cellCentered()

AMREX_GPU_HOST_DEVICE bool amrex::Box::cellCentered ( ) const
inlinenoexcept

Returns true if Box is cell-centered in all indexing directions.

◆ chop()

AMREX_GPU_HOST_DEVICE Box amrex::Box::chop ( int  dir,
int  chop_pnt 
)
inlinenoexcept

Chop the Box at the chop_pnt in the dir direction returns one Box, modifies the object Box. The union of the two is the original Box. The modified Box is the low end, the returned Box is the high end. If type(dir) = CELL, the Boxes are disjoint with the chop_pnt included in the high end (new Box). It is an ERROR if chop_pnt is the low end of the orig Box. If type(dir) = NODE, the chop_pnt is included in both Boxes but is the only point in common. It is also an error if the chop_pnt is an end node of the Box.

◆ coarsen() [1/2]

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Box & amrex::Box::coarsen ( const IntVect ref_ratio)
noexcept

Coarsen Box by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo/ratio and hi <- hi/ratio. NOTE: if type(dir) = NODE centered: lo <- lo/ratio and hi <- hi/ratio + ((hiratio)==0 ? 0 : 1). That is, refinement of coarsened Box must contain the original Box.

◆ coarsen() [2/2]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::coarsen ( int  ref_ratio)
inlinenoexcept

Coarsen Box by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo/ratio and hi <- hi/ratio. NOTE: if type(dir) = NODE centered: lo <- lo/ratio and hi <- hi/ratio + ((hiratio)==0 ? 0 : 1). That is, refinement of coarsened Box must contain the original Box.

◆ coarsenable() [1/3]

AMREX_GPU_HOST_DEVICE bool amrex::Box::coarsenable ( const IntVect refrat,
const IntVect min_width 
) const
inlinenoexcept

◆ coarsenable() [2/3]

AMREX_GPU_HOST_DEVICE bool amrex::Box::coarsenable ( const IntVect refrat,
int  min_width = 1 
) const
inlinenoexcept

◆ coarsenable() [3/3]

AMREX_GPU_HOST_DEVICE bool amrex::Box::coarsenable ( int  refrat,
int  min_width = 1 
) const
inlinenoexcept

◆ contains() [1/4]

AMREX_GPU_HOST_DEVICE bool amrex::Box::contains ( const Box b) const
inlinenoexcept

Returns true if argument is contained within Box. It is an error if the Boxes have different types.

◆ contains() [2/4]

AMREX_GPU_HOST_DEVICE bool amrex::Box::contains ( const Dim3 p) const
inlinenoexcept

Returns true if argument is contained within Box.

◆ contains() [3/4]

AMREX_GPU_HOST_DEVICE bool amrex::Box::contains ( const IntVect p) const
inlinenoexcept

Returns true if argument is contained within Box.

◆ contains() [4/4]

AMREX_GPU_HOST_DEVICE bool amrex::Box::contains ( int  i,
int  ,
int   
) const
inlinenoexcept

Returns true if argument is contained within Box.

◆ convert() [1/2]

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Box & amrex::Box::convert ( const IntVect typ)
noexcept

Convert the Box from the current type into the argument type. This may change the Box coordinates: type CELL -> NODE : increase coordinate by one on high end type NODE -> CELL : reduce coordinate by one on high end other type mappings make no change.

◆ convert() [2/2]

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Box & amrex::Box::convert ( IndexType  typ)
noexcept

Convert the Box from the current type into the argument type. This may change the Box coordinates: type CELL -> NODE : increase coordinate by one on high end type NODE -> CELL : reduce coordinate by one on high end other type mappings make no change.

◆ d_numPts()

AMREX_GPU_HOST_DEVICE double amrex::Box::d_numPts ( ) const
inlinenoexcept

Returns the number of points contained in the Box. This is intended for use only in diagnostic messages.

◆ enclosedCells() [1/3]

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Box & amrex::Box::enclosedCells ( )
noexcept

Convert to CELL type in all directions.

◆ enclosedCells() [2/3]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::enclosedCells ( Direction  d)
inlinenoexcept

◆ enclosedCells() [3/3]

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Box & amrex::Box::enclosedCells ( int  dir)
noexcept

Convert to CELL type in given direction.

◆ grow() [1/5]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::grow ( const IntVect v)
inlinenoexcept

Grow Box in each direction by specified amount.

◆ grow() [2/5]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::grow ( Direction  d,
int  n_cell 
)
inlinenoexcept

◆ grow() [3/5]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::grow ( int  i)
inlinenoexcept

◆ grow() [4/5]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::grow ( int  idir,
int  n_cell 
)
inlinenoexcept

Grow the Box on the low and high end by n_cell cells in direction idir.

◆ grow() [5/5]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::grow ( Orientation  face,
int  n_cell = 1 
)
inlinenoexcept

Grow in the direction of the given face.

◆ growHi() [1/2]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::growHi ( Direction  d,
int  n_cell = 1 
)
inlinenoexcept

◆ growHi() [2/2]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::growHi ( int  idir,
int  n_cell = 1 
)
inlinenoexcept

Grow the Box on the high end by n_cell cells in direction idir. NOTE: n_cell negative shrinks the Box by that number of cells.

◆ growLo() [1/2]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::growLo ( Direction  d,
int  n_cell = 1 
)
inlinenoexcept

◆ growLo() [2/2]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::growLo ( int  idir,
int  n_cell = 1 
)
inlinenoexcept

Grow the Box on the low end by n_cell cells in direction idir. NOTE: n_cell negative shrinks the Box by that number of cells.

◆ hiVect() [1/2]

AMREX_GPU_HOST_DEVICE const int* amrex::Box::hiVect ( ) &&
delete

◆ hiVect() [2/2]

AMREX_GPU_HOST_DEVICE const int* amrex::Box::hiVect ( ) const &
inlinenoexcept

Returns a constant pointer the array of high end coordinates. Useful for calls to FORTRAN.

◆ hiVect3d()

AMREX_GPU_HOST_DEVICE GpuArray<int,3> amrex::Box::hiVect3d ( ) const
inlinenoexcept

◆ index()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Long amrex::Box::index ( const IntVect v) const
noexcept

Returns offset of point from smallend; i.e. index(smallend) -> 0, bigend would return numPts()-1. Is used in accessing FArrayBox.

◆ intersects()

AMREX_GPU_HOST_DEVICE bool amrex::Box::intersects ( const Box b) const
inlinenoexcept

Returns true if Boxes have non-null intersections. It is an error if the Boxes have different types.

◆ isEmpty()

AMREX_GPU_HOST_DEVICE bool amrex::Box::isEmpty ( ) const
inlinenoexcept

Checks if it is an empty box.

◆ isSquare()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool amrex::Box::isSquare ( ) const
noexcept

◆ ixType()

AMREX_GPU_HOST_DEVICE IndexType amrex::Box::ixType ( ) const
inlinenoexcept

Returns the indexing type.

◆ length() [1/2]

AMREX_GPU_HOST_DEVICE IntVect amrex::Box::length ( ) const
inlinenoexcept

Return the length of the Box.

◆ length() [2/2]

AMREX_GPU_HOST_DEVICE int amrex::Box::length ( int  dir) const
inlinenoexcept

Return the length of the Box in given direction.

◆ length3d()

AMREX_GPU_HOST_DEVICE GpuArray<int,3> amrex::Box::length3d ( ) const
inlinenoexcept

◆ longside() [1/2]

AMREX_GPU_HOST_DEVICE int amrex::Box::longside ( ) const
inlinenoexcept

Returns length of longest side. Ignores type.

◆ longside() [2/2]

AMREX_GPU_HOST_DEVICE int amrex::Box::longside ( int dir) const
inlinenoexcept

Returns length of longest side. dir is modified to give direction with longest side: 0...SPACEDIM-1. Ignores type.

◆ loVect() [1/2]

AMREX_GPU_HOST_DEVICE const int* amrex::Box::loVect ( ) &&
delete

◆ loVect() [2/2]

AMREX_GPU_HOST_DEVICE const int* amrex::Box::loVect ( ) const &
inlinenoexcept

Returns a constant pointer the array of low end coordinates. Useful for calls to FORTRAN.

◆ loVect3d()

AMREX_GPU_HOST_DEVICE GpuArray<int,3> amrex::Box::loVect3d ( ) const
inlinenoexcept

◆ makeSlab()

AMREX_GPU_HOST_DEVICE Box& amrex::Box::makeSlab ( int  direction,
int  slab_index 
)
inlinenoexcept

◆ minBox()

AMREX_GPU_HOST_DEVICE Box& amrex::Box::minBox ( const Box b)
inlinenoexcept

Modify Box to that of the minimum Box containing both the original Box and the argument. Both Boxes must have identical type.

◆ next()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void amrex::Box::next ( IntVect p) const
noexcept

Step through the rectangle. It is a runtime error to give a point not inside rectangle. Iteration may not be efficient.

◆ normalize()

AMREX_GPU_HOST_DEVICE void amrex::Box::normalize ( )
inlinenoexcept

◆ numPts()

AMREX_GPU_HOST_DEVICE Long amrex::Box::numPts ( ) const
inlinenoexcept

Returns the number of points contained in the Box.

◆ ok()

AMREX_GPU_HOST_DEVICE bool amrex::Box::ok ( ) const
inlinenoexcept

Checks if it is a proper Box (including a valid type).

◆ operator!=()

AMREX_GPU_HOST_DEVICE bool amrex::Box::operator!= ( const Box b) const
inlinenoexcept

Returns true if Boxes differ (including type).

◆ operator&()

AMREX_GPU_HOST_DEVICE Box amrex::Box::operator& ( const Box rhs) const
inlinenoexcept

Return Box that is intersection of this Box and argument. The Boxes MUST be of same type.

◆ operator&=()

AMREX_GPU_HOST_DEVICE Box& amrex::Box::operator&= ( const Box rhs)
inlinenoexcept

Intersect this Box with its argument. The Boxes MUST be of the same type.

◆ operator+()

AMREX_GPU_HOST_DEVICE Box amrex::Box::operator+ ( const IntVect v) const
inlinenoexcept

Shift Box (relative) by given IntVect.

◆ operator+=()

AMREX_GPU_HOST_DEVICE Box& amrex::Box::operator+= ( const IntVect v)
inlinenoexcept

Shift Box (relative) by given IntVect.

◆ operator-()

AMREX_GPU_HOST_DEVICE Box amrex::Box::operator- ( const IntVect v) const
inlinenoexcept

Shift Box (relative) by given IntVect.

◆ operator-=()

AMREX_GPU_HOST_DEVICE Box& amrex::Box::operator-= ( const IntVect v)
inlinenoexcept

Shift Box (relative) by given IntVect.

◆ operator<()

AMREX_GPU_HOST_DEVICE bool amrex::Box::operator< ( const Box rhs) const
inlinenoexcept

◆ operator<=()

AMREX_GPU_HOST_DEVICE bool amrex::Box::operator<= ( const Box rhs) const
inlinenoexcept

◆ operator==()

AMREX_GPU_HOST_DEVICE bool amrex::Box::operator== ( const Box b) const
inlinenoexcept

Returns true if Boxes are identical (including type).

◆ operator>()

AMREX_GPU_HOST_DEVICE bool amrex::Box::operator> ( const Box rhs) const
inlinenoexcept

◆ operator>=()

AMREX_GPU_HOST_DEVICE bool amrex::Box::operator>= ( const Box rhs) const
inlinenoexcept

◆ operator[]()

AMREX_GPU_HOST_DEVICE int amrex::Box::operator[] ( Orientation  face) const
inlinenoexcept

Returns the coordinate normal to given face.

◆ refine() [1/2]

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Box & amrex::Box::refine ( const IntVect ref_ratio)
noexcept

◆ refine() [2/2]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::refine ( int  ref_ratio)
inlinenoexcept

Refine Box by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo*ratio and hi <- (hi+1)*ratio - 1. NOTE: if type(dir) = NODE centered: lo <- lo*ratio and hi <- hi*ratio.

◆ sameSize()

AMREX_GPU_HOST_DEVICE bool amrex::Box::sameSize ( const Box b) const
inlinenoexcept

Returns true is Boxes same size, ie translates of each other,. It is an error if they have different types.

◆ sameType()

AMREX_GPU_HOST_DEVICE bool amrex::Box::sameType ( const Box b) const
inlinenoexcept

Returns true if Boxes have same type.

◆ setBig() [1/2]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::setBig ( const IntVect bg)
inlinenoexcept

Redefine the big end of the Box.

◆ setBig() [2/2]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::setBig ( int  dir,
int  bg_index 
)
inlinenoexcept

Redefine the big end of the Box.

◆ setRange()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Box & amrex::Box::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. NOTE: This will yield an illegal Box if n_cells <= 0.

◆ setSmall() [1/2]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::setSmall ( const IntVect sm)
inlinenoexcept

Redefine the small end of the Box.

◆ setSmall() [2/2]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::setSmall ( int  dir,
int  sm_index 
)
inlinenoexcept

Redefine the small end of the Box.

◆ setType()

AMREX_GPU_HOST_DEVICE Box& amrex::Box::setType ( const IndexType t)
inlinenoexcept

Set indexing type.

◆ shift() [1/2]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::shift ( const IntVect iv)
inlinenoexcept

Equivalent to b.shift(0,iv[0]).shift(1,iv[1]) ....

◆ shift() [2/2]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::shift ( int  dir,
int  nzones 
)
inlinenoexcept

Shift this Box nzones indexing positions in coordinate direction dir.

◆ shiftHalf() [1/2]

AMREX_GPU_HOST_DEVICE Box & amrex::Box::shiftHalf ( const IntVect iv)
inlinenoexcept

Equivalent to b.shiftHalf(0,iv[0]).shiftHalf(1,iv[1]) ...

◆ shiftHalf() [2/2]

AMREX_GPU_HOST_DEVICE Box & amrex::Box::shiftHalf ( int  dir,
int  num_halfs 
)
inlinenoexcept

This member shifts the Box by "half" indices, thereby converting the Box from type CELL to NODE and visa-versa. b.shiftHalf(0,1) shifts b to the right by 1/2 cells. b.shiftHalf(1,-3) shifts b in the -j direction by 3/2 cells. NOTE: If num_halfs is EVEN the shift is num_halfs/2 full zones and hence will not change the type. This is: b.shifthalf(4) == b.shift(2).

◆ shortside() [1/2]

AMREX_GPU_HOST_DEVICE int amrex::Box::shortside ( ) const
inlinenoexcept

Returns length of shortest side. Ignores type.

◆ shortside() [2/2]

AMREX_GPU_HOST_DEVICE int amrex::Box::shortside ( int dir) const
inlinenoexcept

Returns length of shortest side. dir is modified to give direction with shortest side: 0...SPACEDIM-1. Ignores type.

◆ size()

AMREX_GPU_HOST_DEVICE IntVect amrex::Box::size ( ) const
inlinenoexcept

Return the length of the Box.

◆ smallEnd() [1/3]

const IntVect& amrex::Box::smallEnd ( ) &&
delete

Get the smallend of the box.

◆ smallEnd() [2/3]

AMREX_GPU_HOST_DEVICE const IntVect& amrex::Box::smallEnd ( ) const &
inlinenoexcept

Get the smallend of the box.

◆ smallEnd() [3/3]

AMREX_GPU_HOST_DEVICE int amrex::Box::smallEnd ( int  dir) const &
inlinenoexcept

Returns the coordinate of the low end in the given direction.

◆ strictly_contains() [1/3]

AMREX_GPU_HOST_DEVICE bool amrex::Box::strictly_contains ( const Box b) const
inlinenoexcept

Returns true if argument is strictly contained within Box. It is an error if the Boxes have different types.

◆ strictly_contains() [2/3]

AMREX_GPU_HOST_DEVICE bool amrex::Box::strictly_contains ( const IntVect p) const
inlinenoexcept

Returns true if argument is strictly contained within Box.

◆ strictly_contains() [3/3]

AMREX_GPU_HOST_DEVICE bool amrex::Box::strictly_contains ( int  i,
int  ,
int   
) const
inlinenoexcept

Returns true if argument is strictly contained within Box.

◆ surroundingNodes() [1/3]

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Box & amrex::Box::surroundingNodes ( )
noexcept

Convert to NODE type in all directions.

◆ surroundingNodes() [2/3]

AMREX_GPU_HOST_DEVICE Box& amrex::Box::surroundingNodes ( Direction  d)
inlinenoexcept

◆ surroundingNodes() [3/3]

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Box & amrex::Box::surroundingNodes ( int  dir)
noexcept

Convert to NODE type in given direction.

◆ TheUnitBox()

static AMREX_GPU_HOST_DEVICE Box amrex::Box::TheUnitBox ( )
inlinestaticnoexcept

This static member function returns a constant reference to an object of type Box representing the unit box in AMREX_SPACEDIM-dimensional space.

◆ type() [1/2]

AMREX_GPU_HOST_DEVICE IntVect amrex::Box::type ( ) const
inlinenoexcept

Returns the indexing type.

◆ type() [2/2]

AMREX_GPU_HOST_DEVICE IndexType::CellIndex amrex::Box::type ( int  dir) const
inlinenoexcept

Returns the indexing type in the specified direction.

◆ volume()

AMREX_GPU_HOST_DEVICE Long amrex::Box::volume ( ) const
inlinenoexcept

Return the volume, in indexing space, of region enclosed by this Box. This is identical to numPts() for CELL centered Box; otherwise, numPts() > volume().

Friends And Related Function Documentation

◆ begin

AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 begin ( Box const &  box)
friend

◆ BoxCommHelper

friend class BoxCommHelper
friend

◆ end

AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 end ( Box const &  box)
friend

◆ lbound

AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 lbound ( Box const &  box)
friend

◆ length

AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 length ( Box const &  box)
friend

◆ max_lbound [1/2]

AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 max_lbound ( Box const &  b1,
Box const &  b2 
)
friend

◆ max_lbound [2/2]

AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 max_lbound ( Box const &  b1,
Dim3 const &  lo 
)
friend

◆ min_ubound [1/2]

AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 min_ubound ( Box const &  b1,
Box const &  b2 
)
friend

◆ min_ubound [2/2]

AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 min_ubound ( Box const &  b1,
Dim3 const &  hi 
)
friend

◆ minBox

AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Box minBox ( Box const &  b1,
Box const &  b2,
IndexType  typ 
)
friend

◆ ParallelDescriptor::Mpi_typemap

◆ ubound

AMREX_GPU_HOST_DEVICE friend AMREX_FORCE_INLINE Dim3 ubound ( Box const &  box)
friend

Member Data Documentation

◆ bigend

IntVect amrex::Box::bigend
private

◆ btype

IndexType amrex::Box::btype
private

◆ smallend

IntVect amrex::Box::smallend
private

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