Block-Structured AMR Software Framework
 
Loading...
Searching...
No Matches
amrex::BoxND< dim > Class Template Reference

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

#include <AMReX_Box.H>

Public Member Functions

__host__ __device__ constexpr BoxND () noexcept
 
__host__ __device__ constexpr BoxND (const IntVectND< dim > &small, const IntVectND< dim > &big) noexcept
 Construct cell-centered type BoxND.
 
__host__ __device__ BoxND (const IntVectND< dim > &small, const int *vec_len) noexcept
 Construct BoxND with specified lengths.
 
__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.
 
__host__ __device__ BoxND (const IntVectND< dim > &small, const IntVectND< dim > &big, IndexTypeND< dim > t) noexcept
 Construct dimension specific Boxes.
 
template<typename T , int Tdim = dim, std::enable_if_t<(1<=Tdim &&Tdim<=3), int > = 0>
__host__ __device__ BoxND (Array4< T > const &a) noexcept
 
__host__ __device__ const IntVectND< dim > & smallEnd () const &noexcept
 Get the smallend of the BoxND.
 
const IntVectND< dim > & smallEnd () &&=delete
 Get the smallend of the BoxND.
 
__host__ __device__ int smallEnd (int dir) const &noexcept
 Returns the coordinate of the low end in the given direction.
 
__host__ __device__ const IntVectND< dim > & bigEnd () const &noexcept
 Get the bigend.
 
const IntVectND< dim > & bigEnd () &&=delete
 Get the bigend.
 
__host__ __device__ int bigEnd (int dir) const noexcept
 Returns the coordinate of the high end in the given direction.
 
__host__ __device__ IndexTypeND< dim > ixType () const noexcept
 Returns the indexing type.
 
__host__ __device__ IntVectND< dim > type () const noexcept
 Returns the indexing type.
 
__host__ __device__ IndexType::CellIndex type (int dir) const noexcept
 Returns the indexing type in the specified direction.
 
__host__ __device__ IntVectND< dim > size () const noexcept
 Return the length of the BoxND.
 
__host__ __device__ IntVectND< dim > length () const noexcept
 Return the length of the BoxND.
 
__host__ __device__ int length (int dir) const noexcept
 Return the length of the BoxND in given direction.
 
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ GpuArray< int, 3 > length3d () const noexcept
 
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ GpuArray< int, 3 > loVect3d () const noexcept
 
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ GpuArray< int, 3 > hiVect3d () const noexcept
 
__host__ __device__ const int * loVect () const &noexcept
 Returns a constant pointer the array of low end coordinates. Useful for calls to FORTRAN.
 
__host__ __device__ const int * loVect () &&=delete
 
__host__ __device__ const int * hiVect () const &noexcept
 Returns a constant pointer the array of high end coordinates. Useful for calls to FORTRAN.
 
__host__ __device__ const int * hiVect () &&=delete
 
__host__ __device__ int operator[] (Orientation face) const noexcept
 Returns the coordinate normal to given face.
 
__host__ __device__ bool isEmpty () const noexcept
 Checks if it is an empty BoxND.
 
__host__ __device__ bool ok () const noexcept
 Checks if it is a proper BoxND (including a valid type).
 
__host__ __device__ bool contains (const IntVectND< dim > &p) const noexcept
 Returns true if argument is contained within BoxND.
 
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ bool contains (const Dim3 &p) const noexcept
 Returns true if argument is contained within BoxND.
 
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ bool contains (int i, int j, int k) const noexcept
 Returns true if argument is contained within BoxND.
 
__host__ __device__ bool contains (const BoxND &b) const noexcept
 Returns true if argument is contained within BoxND. It is an error if the Boxes have different types.
 
__host__ __device__ bool strictly_contains (const IntVectND< dim > &p) const noexcept
 Returns true if argument is strictly contained within BoxND.
 
__host__ __device__ bool strictly_contains (const BoxND &b) const noexcept
 Returns true if argument is strictly contained within BoxND. It is an error if the Boxes have different types.
 
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ bool strictly_contains (const Dim3 &p) const noexcept
 Returns true if argument is strictly contained within BoxND.
 
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ bool strictly_contains (int i, int j, int k) const noexcept
 Returns true if argument is strictly contained within BoxND.
 
__host__ __device__ bool intersects (const BoxND &b) const noexcept
 Returns true if Boxes have non-null intersections. It is an error if the Boxes have different types.
 
__host__ __device__ bool sameSize (const BoxND &b) const noexcept
 Returns true is Boxes same size, ie translates of each other,. It is an error if they have different types.
 
__host__ __device__ bool sameType (const BoxND &b) const noexcept
 Returns true if Boxes have same type.
 
__host__ __device__ bool operator== (const BoxND &b) const noexcept
 Returns true if Boxes are identical (including type).
 
__host__ __device__ bool operator!= (const BoxND &b) const noexcept
 Returns true if Boxes differ (including type).
 
__host__ __device__ bool operator< (const BoxND &rhs) const noexcept
 
__host__ __device__ bool operator<= (const BoxND &rhs) const noexcept
 
__host__ __device__ bool operator> (const BoxND &rhs) const noexcept
 
__host__ __device__ bool operator>= (const BoxND &rhs) const noexcept
 
__host__ __device__ bool cellCentered () const noexcept
 Returns true if BoxND is cell-centered in all indexing directions.
 
void checkOverflow () const noexcept
 Assert that there are no int/Long overflows when calling length or numPts.
 
__host__ __device__ Long numPts () const noexcept
 Returns the number of points contained in the BoxND.
 
__host__ __device__ double d_numPts () const noexcept
 Returns the number of points contained in the BoxND. This is intended for use only in diagnostic messages.
 
__host__ __device__ Long volume () const noexcept
 Return the volume, in indexing space, of region enclosed by this BoxND. This is identical to numPts() for CELL centered BoxND; otherwise, numPts() > volume().
 
__host__ __device__ int longside (int &dir) const noexcept
 Returns length of longest side. dir is modified to give direction with longest side: 0...dim-1. Ignores type.
 
__host__ __device__ int longside () const noexcept
 Returns length of longest side. Ignores type.
 
__host__ __device__ int shortside (int &dir) const noexcept
 Returns length of shortest side. dir is modified to give direction with shortest side: 0...dim-1. Ignores type.
 
__host__ __device__ int shortside () const noexcept
 Returns length of shortest side. Ignores type.
 
__host__ __device__ Long index (const IntVectND< dim > &v) const noexcept
 Returns offset of point from smallend; i.e. index(smallend) -> 0, bigend would return numPts()-1. Is used in accessing FArrayBox.
 
__host__ __device__ IntVectND< dim > atOffset (Long offset) const noexcept
 Given the offset, compute IntVectND<dim>
 
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ GpuArray< int, 3 > atOffset3d (Long offset) const noexcept
 
__host__ __device__ BoxNDsetSmall (const IntVectND< dim > &sm) noexcept
 Redefine the small end of the BoxND.
 
__host__ __device__ BoxNDsetSmall (int dir, int sm_index) noexcept
 Redefine the small end of the BoxND.
 
__host__ __device__ BoxNDsetBig (const IntVectND< dim > &bg) noexcept
 Redefine the big end of the BoxND.
 
__host__ __device__ BoxNDsetBig (int dir, int bg_index) noexcept
 Redefine the big end of the BoxND.
 
__host__ __device__ BoxNDsetRange (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 BoxND if n_cells <= 0.
 
__host__ __device__ BoxNDsetType (const IndexTypeND< dim > &t) noexcept
 Set indexing type.
 
__host__ __device__ BoxNDshift (int dir, int nzones) noexcept
 Shift this BoxND nzones indexing positions in coordinate direction dir.
 
__host__ __device__ BoxNDshift (const IntVectND< dim > &iv) noexcept
 Equivalent to b.shift(0,iv[0]).shift(1,iv[1]) ....
 
__host__ __device__ BoxNDshiftHalf (int dir, int num_halfs) noexcept
 This member shifts the BoxND by "half" indices, thereby converting the BoxND 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).
 
__host__ __device__ BoxNDshiftHalf (const IntVectND< dim > &iv) noexcept
 Equivalent to b.shiftHalf(0,iv[0]).shiftHalf(1,iv[1]) ...
 
__host__ __device__ BoxNDconvert (IndexTypeND< dim > typ) noexcept
 Convert the BoxND from the current type into the argument type. This may change the BoxND 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.
 
__host__ __device__ BoxNDconvert (const IntVectND< dim > &typ) noexcept
 Convert the BoxND from the current type into the argument type. This may change the BoxND 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.
 
__host__ __device__ BoxNDsurroundingNodes () noexcept
 Convert to NODE type in all directions.
 
__host__ __device__ BoxNDsurroundingNodes (int dir) noexcept
 Convert to NODE type in given direction.
 
__host__ __device__ BoxNDsurroundingNodes (Direction d) noexcept
 
__host__ __device__ BoxNDenclosedCells () noexcept
 Convert to CELL type in all directions.
 
__host__ __device__ BoxNDenclosedCells (int dir) noexcept
 Convert to CELL type in given direction.
 
__host__ __device__ BoxNDenclosedCells (Direction d) noexcept
 
__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.
 
__host__ __device__ BoxNDoperator&= (const BoxND &rhs) noexcept
 Intersect this BoxND with its argument. The Boxes MUST be of the same type.
 
__host__ __device__ BoxNDminBox (const BoxND &b) noexcept
 Modify BoxND to that of the minimum BoxND containing both the original BoxND and the argument. Both Boxes must have identical type.
 
__host__ __device__ BoxNDoperator+= (const IntVectND< dim > &v) noexcept
 Shift BoxND (relative) by given IntVectND<dim>.
 
__host__ __device__ BoxND operator+ (const IntVectND< dim > &v) const noexcept
 Shift BoxND (relative) by given IntVectND<dim>.
 
__host__ __device__ BoxNDoperator-= (const IntVectND< dim > &v) noexcept
 Shift BoxND (relative) by given IntVectND<dim>.
 
__host__ __device__ BoxND operator- (const IntVectND< dim > &v) const noexcept
 Shift BoxND (relative) by given IntVectND<dim>.
 
__host__ __device__ BoxND chop (int dir, int chop_pnt) noexcept
 Chop the BoxND at the chop_pnt in the dir direction returns one BoxND, modifies the object BoxND. The union of the two is the original BoxND. The modified BoxND is the low end, the returned BoxND is the high end. If type(dir) = CELL, the Boxes are disjoint with the chop_pnt included in the high end (new BoxND). It is an ERROR if chop_pnt is the low end of the orig BoxND. 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 BoxND.
 
__host__ __device__ BoxNDgrow (int i) noexcept
 
__host__ __device__ BoxNDgrow (const IntVectND< dim > &v) noexcept
 Grow BoxND in each direction by specified amount.
 
__host__ __device__ BoxNDgrow (int idir, int n_cell) noexcept
 Grow the BoxND on the low and high end by n_cell cells in direction idir.
 
__host__ __device__ BoxNDgrow (Direction d, int n_cell) noexcept
 
__host__ __device__ BoxNDgrowLo (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 BoxND by that number of cells.
 
__host__ __device__ BoxNDgrowLo (Direction d, int n_cell=1) noexcept
 
__host__ __device__ BoxNDgrowHi (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 BoxND by that number of cells.
 
__host__ __device__ BoxNDgrowHi (Direction d, int n_cell=1) noexcept
 
__host__ __device__ BoxNDgrow (Orientation face, int n_cell=1) noexcept
 Grow in the direction of the given face.
 
__host__ __device__ BoxNDrefine (int ref_ratio) noexcept
 Refine BoxND 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.
 
__host__ __device__ BoxNDrefine (const IntVectND< dim > &ref_ratio) noexcept
 
__host__ __device__ BoxNDcoarsen (int ref_ratio) noexcept
 Coarsen BoxND 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 BoxND must contain the original BoxND.
 
__host__ __device__ BoxNDcoarsen (const IntVectND< dim > &ref_ratio) noexcept
 Coarsen BoxND 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 BoxND must contain the original BoxND.
 
__host__ __device__ void next (IntVectND< dim > &) const noexcept
 Step through the rectangle. It is a runtime error to give a point not inside rectangle. Iteration may not be efficient.
 
__host__ __device__ bool isSquare () const noexcept
 
__host__ __device__ bool coarsenable (const IntVectND< dim > &refrat, const IntVectND< dim > &min_width) const noexcept
 
__host__ __device__ bool coarsenable (int refrat, int min_width=1) const noexcept
 
__host__ __device__ bool coarsenable (const IntVectND< dim > &refrat, int min_width=1) const noexcept
 
__host__ __device__ void normalize () noexcept
 
__host__ __device__ BoxNDmakeSlab (int direction, int slab_index) noexcept
 
BoxIteratorND< dim > iterator () const noexcept
 Returns a BoxIteratorND that can be used to loop over the IntVects contained by the BoxND.
 
template<int new_dim>
__host__ __device__ BoxND< new_dim > shrink () const noexcept
 Returns a new BoxND of dimension new_dim and assigns the first new_dim dimension of this BoxND to it.
 
template<int new_dim>
__host__ __device__ BoxND< new_dim > expand () const noexcept
 Returns a new BoxND of size new_dim and assigns all values of this BoxND to it and (small=0, big=0, typ=CELL) to the remaining elements.
 
template<int new_dim>
__host__ __device__ BoxND< new_dim > resize () const noexcept
 Returns a new BoxND of size new_dim by either shrinking or expanding this BoxND.
 

Static Public Member Functions

__host__ static __device__ BoxND TheUnitBox () noexcept
 This static member function returns a constant reference to an object of type BoxND representing the unit BoxND in dim-dimensional space.
 
__host__ static __device__ constexpr std::size_t ndims () noexcept
 
__host__ static __device__ constexpr int indims () noexcept
 

Private Attributes

IntVectND< dim > smallend
 
IntVectND< dim > bigend
 
IndexTypeND< dim > btype
 

Friends

class BoxCommHelper
 
MPI_Datatype ParallelDescriptor::Mpi_typemap ()
 

Detailed Description

template<int dim>
class amrex::BoxND< dim >

A Rectangular Domain on an Integer Lattice.

A BoxND is an abstraction for defining discrete regions of dim 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 BoxND. Boxes can exist in positive and negative indexing space.

Box is a dimension dependent alias to BoxND<AMREX_SPACEDIM>, so AMREX_SPACEDIM must be defined as either 1, 2, or 3 when compiling.

Constructor & Destructor Documentation

◆ BoxND() [1/6]

template<int dim>
__host__ __device__ constexpr amrex::BoxND< dim >::BoxND ( )
inlineconstexprnoexcept

◆ BoxND() [2/6]

template<int dim>
__host__ __device__ constexpr amrex::BoxND< dim >::BoxND ( const IntVectND< dim > &  small,
const IntVectND< dim > &  big 
)
inlineconstexprnoexcept

Construct cell-centered type BoxND.

◆ BoxND() [3/6]

template<int dim>
__host__ __device__ amrex::BoxND< dim >::BoxND ( const IntVectND< dim > &  small,
const int *  vec_len 
)
inlinenoexcept

Construct BoxND with specified lengths.

◆ BoxND() [4/6]

template<int dim>
__host__ __device__ amrex::BoxND< dim >::BoxND ( const IntVectND< dim > &  small,
const IntVectND< dim > &  big,
const IntVectND< dim > &  typ 
)
inlinenoexcept

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

◆ BoxND() [5/6]

template<int dim>
__host__ __device__ amrex::BoxND< dim >::BoxND ( const IntVectND< dim > &  small,
const IntVectND< dim > &  big,
IndexTypeND< dim >  t 
)
inlinenoexcept

Construct dimension specific Boxes.

◆ BoxND() [6/6]

template<int dim>
template<typename T , int Tdim = dim, std::enable_if_t<(1<=Tdim &&Tdim<=3), int > = 0>
__host__ __device__ amrex::BoxND< dim >::BoxND ( Array4< T > const &  a)
inlineexplicitnoexcept

Member Function Documentation

◆ atOffset()

template<int dim>
__host__ __device__ IntVectND< dim > amrex::BoxND< dim >::atOffset ( Long  offset) const
inlinenoexcept

Given the offset, compute IntVectND<dim>

◆ atOffset3d()

template<int dim>
template<int N, std::enable_if_t<(1<=N &&N<=3), int > >
__host__ __device__ GpuArray< int, 3 > amrex::BoxND< dim >::atOffset3d ( Long  offset) const
inlinenoexcept

◆ bigEnd() [1/3]

template<int dim>
const IntVectND< dim > & amrex::BoxND< dim >::bigEnd ( ) &&
delete

Get the bigend.

◆ bigEnd() [2/3]

template<int dim>
__host__ __device__ const IntVectND< dim > & amrex::BoxND< dim >::bigEnd ( ) const &
inlinenoexcept

Get the bigend.

◆ bigEnd() [3/3]

template<int dim>
__host__ __device__ int amrex::BoxND< dim >::bigEnd ( int  dir) const
inlinenoexcept

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

◆ cellCentered()

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::cellCentered ( ) const
inlinenoexcept

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

◆ checkOverflow()

template<int dim>
void amrex::BoxND< dim >::checkOverflow ( ) const
inlinenoexcept

Assert that there are no int/Long overflows when calling length or numPts.

◆ chop()

template<int dim>
__host__ __device__ BoxND< dim > amrex::BoxND< dim >::chop ( int  dir,
int  chop_pnt 
)
inlinenoexcept

Chop the BoxND at the chop_pnt in the dir direction returns one BoxND, modifies the object BoxND. The union of the two is the original BoxND. The modified BoxND is the low end, the returned BoxND is the high end. If type(dir) = CELL, the Boxes are disjoint with the chop_pnt included in the high end (new BoxND). It is an ERROR if chop_pnt is the low end of the orig BoxND. 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 BoxND.

◆ coarsen() [1/2]

template<int dim>
__host__ __device__ BoxND< dim > & amrex::BoxND< dim >::coarsen ( const IntVectND< dim > &  ref_ratio)
inlinenoexcept

Coarsen BoxND 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 BoxND must contain the original BoxND.

◆ coarsen() [2/2]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::coarsen ( int  ref_ratio)
inlinenoexcept

Coarsen BoxND 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 BoxND must contain the original BoxND.

◆ coarsenable() [1/3]

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::coarsenable ( const IntVectND< dim > &  refrat,
const IntVectND< dim > &  min_width 
) const
inlinenoexcept

◆ coarsenable() [2/3]

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::coarsenable ( const IntVectND< dim > &  refrat,
int  min_width = 1 
) const
inlinenoexcept

◆ coarsenable() [3/3]

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::coarsenable ( int  refrat,
int  min_width = 1 
) const
inlinenoexcept

◆ contains() [1/4]

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::contains ( const BoxND< dim > &  b) const
inlinenoexcept

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

◆ contains() [2/4]

template<int dim>
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ bool amrex::BoxND< dim >::contains ( const Dim3 p) const
inlinenoexcept

Returns true if argument is contained within BoxND.

◆ contains() [3/4]

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::contains ( const IntVectND< dim > &  p) const
inlinenoexcept

Returns true if argument is contained within BoxND.

◆ contains() [4/4]

template<int dim>
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ bool amrex::BoxND< dim >::contains ( int  i,
int  j,
int  k 
) const
inlinenoexcept

Returns true if argument is contained within BoxND.

◆ convert() [1/2]

template<int dim>
__host__ __device__ BoxND< dim > & amrex::BoxND< dim >::convert ( const IntVectND< dim > &  typ)
inlinenoexcept

Convert the BoxND from the current type into the argument type. This may change the BoxND 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]

template<int dim>
__host__ __device__ BoxND< dim > & amrex::BoxND< dim >::convert ( IndexTypeND< dim >  typ)
inlinenoexcept

Convert the BoxND from the current type into the argument type. This may change the BoxND 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()

template<int dim>
__host__ __device__ double amrex::BoxND< dim >::d_numPts ( ) const
inlinenoexcept

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

◆ enclosedCells() [1/3]

template<int dim>
__host__ __device__ BoxND< dim > & amrex::BoxND< dim >::enclosedCells ( )
inlinenoexcept

Convert to CELL type in all directions.

◆ enclosedCells() [2/3]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::enclosedCells ( Direction  d)
inlinenoexcept

◆ enclosedCells() [3/3]

template<int dim>
__host__ __device__ BoxND< dim > & amrex::BoxND< dim >::enclosedCells ( int  dir)
inlinenoexcept

Convert to CELL type in given direction.

◆ expand()

template<int dim>
template<int new_dim>
__host__ __device__ BoxND< new_dim > amrex::BoxND< dim >::expand ( ) const
inlinenoexcept

Returns a new BoxND of size new_dim and assigns all values of this BoxND to it and (small=0, big=0, typ=CELL) to the remaining elements.

◆ grow() [1/5]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::grow ( const IntVectND< dim > &  v)
inlinenoexcept

Grow BoxND in each direction by specified amount.

◆ grow() [2/5]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::grow ( Direction  d,
int  n_cell 
)
inlinenoexcept

◆ grow() [3/5]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::grow ( int  i)
inlinenoexcept

◆ grow() [4/5]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::grow ( int  idir,
int  n_cell 
)
inlinenoexcept

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

◆ grow() [5/5]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::grow ( Orientation  face,
int  n_cell = 1 
)
inlinenoexcept

Grow in the direction of the given face.

◆ growHi() [1/2]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::growHi ( Direction  d,
int  n_cell = 1 
)
inlinenoexcept

◆ growHi() [2/2]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::growHi ( int  idir,
int  n_cell = 1 
)
inlinenoexcept

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

◆ growLo() [1/2]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::growLo ( Direction  d,
int  n_cell = 1 
)
inlinenoexcept

◆ growLo() [2/2]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::growLo ( int  idir,
int  n_cell = 1 
)
inlinenoexcept

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

◆ hiVect() [1/2]

template<int dim>
__host__ __device__ const int * amrex::BoxND< dim >::hiVect ( ) &&
delete

◆ hiVect() [2/2]

template<int dim>
__host__ __device__ const int * amrex::BoxND< dim >::hiVect ( ) const &
inlinenoexcept

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

◆ hiVect3d()

template<int dim>
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ GpuArray< int, 3 > amrex::BoxND< dim >::hiVect3d ( ) const
inlinenoexcept

◆ index()

template<int dim>
__host__ __device__ Long amrex::BoxND< dim >::index ( const IntVectND< dim > &  v) const
inlinenoexcept

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

◆ indims()

template<int dim>
__host__ static __device__ constexpr int amrex::BoxND< dim >::indims ( )
inlinestaticconstexprnoexcept

◆ intersects()

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::intersects ( const BoxND< dim > &  b) const
inlinenoexcept

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

◆ isEmpty()

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::isEmpty ( ) const
inlinenoexcept

Checks if it is an empty BoxND.

◆ isSquare()

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::isSquare ( ) const
inlinenoexcept

◆ iterator()

template<int dim>
BoxIteratorND< dim > amrex::BoxND< dim >::iterator ( ) const
inlinenoexcept

Returns a BoxIteratorND that can be used to loop over the IntVects contained by the BoxND.

...
for (IntVect iv : b.iterator()) {
// do operations involving iv
}
BoxIteratorND< dim > iterator() const noexcept
Returns a BoxIteratorND that can be used to loop over the IntVects contained by the BoxND.
Definition AMReX_Box.H:800

◆ ixType()

template<int dim>
__host__ __device__ IndexTypeND< dim > amrex::BoxND< dim >::ixType ( ) const
inlinenoexcept

Returns the indexing type.

◆ length() [1/2]

template<int dim>
__host__ __device__ IntVectND< dim > amrex::BoxND< dim >::length ( ) const
inlinenoexcept

Return the length of the BoxND.

◆ length() [2/2]

template<int dim>
__host__ __device__ int amrex::BoxND< dim >::length ( int  dir) const
inlinenoexcept

Return the length of the BoxND in given direction.

◆ length3d()

template<int dim>
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ GpuArray< int, 3 > amrex::BoxND< dim >::length3d ( ) const
inlinenoexcept

◆ longside() [1/2]

template<int dim>
__host__ __device__ int amrex::BoxND< dim >::longside ( ) const
inlinenoexcept

Returns length of longest side. Ignores type.

◆ longside() [2/2]

template<int dim>
__host__ __device__ int amrex::BoxND< dim >::longside ( int &  dir) const
inlinenoexcept

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

◆ loVect() [1/2]

template<int dim>
__host__ __device__ const int * amrex::BoxND< dim >::loVect ( ) &&
delete

◆ loVect() [2/2]

template<int dim>
__host__ __device__ const int * amrex::BoxND< dim >::loVect ( ) const &
inlinenoexcept

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

◆ loVect3d()

template<int dim>
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ GpuArray< int, 3 > amrex::BoxND< dim >::loVect3d ( ) const
inlinenoexcept

◆ makeSlab()

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::makeSlab ( int  direction,
int  slab_index 
)
inlinenoexcept

◆ minBox()

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::minBox ( const BoxND< dim > &  b)
inlinenoexcept

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

◆ ndims()

template<int dim>
__host__ static __device__ constexpr std::size_t amrex::BoxND< dim >::ndims ( )
inlinestaticconstexprnoexcept

◆ next()

template<int dim>
__host__ __device__ void amrex::BoxND< dim >::next ( IntVectND< dim > &  p) const
inlinenoexcept

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

◆ normalize()

template<int dim>
__host__ __device__ void amrex::BoxND< dim >::normalize ( )
inlinenoexcept

◆ numPts()

template<int dim>
__host__ __device__ Long amrex::BoxND< dim >::numPts ( ) const
inlinenoexcept

Returns the number of points contained in the BoxND.

◆ ok()

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::ok ( ) const
inlinenoexcept

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

◆ operator!=()

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::operator!= ( const BoxND< dim > &  b) const
inlinenoexcept

Returns true if Boxes differ (including type).

◆ operator&()

template<int dim>
__host__ __device__ BoxND amrex::BoxND< dim >::operator& ( const BoxND< dim > &  rhs) const
inlinenoexcept

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

◆ operator&=()

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::operator&= ( const BoxND< dim > &  rhs)
inlinenoexcept

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

◆ operator+()

template<int dim>
__host__ __device__ BoxND amrex::BoxND< dim >::operator+ ( const IntVectND< dim > &  v) const
inlinenoexcept

Shift BoxND (relative) by given IntVectND<dim>.

◆ operator+=()

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::operator+= ( const IntVectND< dim > &  v)
inlinenoexcept

Shift BoxND (relative) by given IntVectND<dim>.

◆ operator-()

template<int dim>
__host__ __device__ BoxND amrex::BoxND< dim >::operator- ( const IntVectND< dim > &  v) const
inlinenoexcept

Shift BoxND (relative) by given IntVectND<dim>.

◆ operator-=()

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::operator-= ( const IntVectND< dim > &  v)
inlinenoexcept

Shift BoxND (relative) by given IntVectND<dim>.

◆ operator<()

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::operator< ( const BoxND< dim > &  rhs) const
inlinenoexcept

◆ operator<=()

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::operator<= ( const BoxND< dim > &  rhs) const
inlinenoexcept

◆ operator==()

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::operator== ( const BoxND< dim > &  b) const
inlinenoexcept

Returns true if Boxes are identical (including type).

◆ operator>()

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::operator> ( const BoxND< dim > &  rhs) const
inlinenoexcept

◆ operator>=()

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::operator>= ( const BoxND< dim > &  rhs) const
inlinenoexcept

◆ operator[]()

template<int dim>
__host__ __device__ int amrex::BoxND< dim >::operator[] ( Orientation  face) const
inlinenoexcept

Returns the coordinate normal to given face.

◆ refine() [1/2]

template<int dim>
__host__ __device__ BoxND< dim > & amrex::BoxND< dim >::refine ( const IntVectND< dim > &  ref_ratio)
inlinenoexcept

◆ refine() [2/2]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::refine ( int  ref_ratio)
inlinenoexcept

Refine BoxND 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.

◆ resize()

template<int dim>
template<int new_dim>
__host__ __device__ BoxND< new_dim > amrex::BoxND< dim >::resize ( ) const
inlinenoexcept

Returns a new BoxND of size new_dim by either shrinking or expanding this BoxND.

◆ sameSize()

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::sameSize ( const BoxND< dim > &  b) const
inlinenoexcept

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

◆ sameType()

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::sameType ( const BoxND< dim > &  b) const
inlinenoexcept

Returns true if Boxes have same type.

◆ setBig() [1/2]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::setBig ( const IntVectND< dim > &  bg)
inlinenoexcept

Redefine the big end of the BoxND.

◆ setBig() [2/2]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::setBig ( int  dir,
int  bg_index 
)
inlinenoexcept

Redefine the big end of the BoxND.

◆ setRange()

template<int dim>
__host__ __device__ BoxND< dim > & amrex::BoxND< dim >::setRange ( int  dir,
int  sm_index,
int  n_cells = 1 
)
inlinenoexcept

Set the entire range in a given direction, starting at sm_index with length n_cells. NOTE: This will yield an illegal BoxND if n_cells <= 0.

◆ setSmall() [1/2]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::setSmall ( const IntVectND< dim > &  sm)
inlinenoexcept

Redefine the small end of the BoxND.

◆ setSmall() [2/2]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::setSmall ( int  dir,
int  sm_index 
)
inlinenoexcept

Redefine the small end of the BoxND.

◆ setType()

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::setType ( const IndexTypeND< dim > &  t)
inlinenoexcept

Set indexing type.

◆ shift() [1/2]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::shift ( const IntVectND< dim > &  iv)
inlinenoexcept

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

◆ shift() [2/2]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::shift ( int  dir,
int  nzones 
)
inlinenoexcept

Shift this BoxND nzones indexing positions in coordinate direction dir.

◆ shiftHalf() [1/2]

template<int dim>
__host__ __device__ BoxND< dim > & amrex::BoxND< dim >::shiftHalf ( const IntVectND< dim > &  iv)
inlinenoexcept

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

◆ shiftHalf() [2/2]

template<int dim>
__host__ __device__ BoxND< dim > & amrex::BoxND< dim >::shiftHalf ( int  dir,
int  num_halfs 
)
inlinenoexcept

This member shifts the BoxND by "half" indices, thereby converting the BoxND 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]

template<int dim>
__host__ __device__ int amrex::BoxND< dim >::shortside ( ) const
inlinenoexcept

Returns length of shortest side. Ignores type.

◆ shortside() [2/2]

template<int dim>
__host__ __device__ int amrex::BoxND< dim >::shortside ( int &  dir) const
inlinenoexcept

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

◆ shrink()

template<int dim>
template<int new_dim>
__host__ __device__ BoxND< new_dim > amrex::BoxND< dim >::shrink ( ) const
inlinenoexcept

Returns a new BoxND of dimension new_dim and assigns the first new_dim dimension of this BoxND to it.

◆ size()

template<int dim>
__host__ __device__ IntVectND< dim > amrex::BoxND< dim >::size ( ) const
inlinenoexcept

Return the length of the BoxND.

◆ smallEnd() [1/3]

template<int dim>
const IntVectND< dim > & amrex::BoxND< dim >::smallEnd ( ) &&
delete

Get the smallend of the BoxND.

◆ smallEnd() [2/3]

template<int dim>
__host__ __device__ const IntVectND< dim > & amrex::BoxND< dim >::smallEnd ( ) const &
inlinenoexcept

Get the smallend of the BoxND.

◆ smallEnd() [3/3]

template<int dim>
__host__ __device__ int amrex::BoxND< dim >::smallEnd ( int  dir) const &
inlinenoexcept

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

◆ strictly_contains() [1/4]

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::strictly_contains ( const BoxND< dim > &  b) const
inlinenoexcept

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

◆ strictly_contains() [2/4]

template<int dim>
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ bool amrex::BoxND< dim >::strictly_contains ( const Dim3 p) const
inlinenoexcept

Returns true if argument is strictly contained within BoxND.

◆ strictly_contains() [3/4]

template<int dim>
__host__ __device__ bool amrex::BoxND< dim >::strictly_contains ( const IntVectND< dim > &  p) const
inlinenoexcept

Returns true if argument is strictly contained within BoxND.

◆ strictly_contains() [4/4]

template<int dim>
template<int N = dim, std::enable_if_t<(1<=N &&N<=3), int > = 0>
__host__ __device__ bool amrex::BoxND< dim >::strictly_contains ( int  i,
int  j,
int  k 
) const
inlinenoexcept

Returns true if argument is strictly contained within BoxND.

◆ surroundingNodes() [1/3]

template<int dim>
__host__ __device__ BoxND< dim > & amrex::BoxND< dim >::surroundingNodes ( )
inlinenoexcept

Convert to NODE type in all directions.

◆ surroundingNodes() [2/3]

template<int dim>
__host__ __device__ BoxND & amrex::BoxND< dim >::surroundingNodes ( Direction  d)
inlinenoexcept

◆ surroundingNodes() [3/3]

template<int dim>
__host__ __device__ BoxND< dim > & amrex::BoxND< dim >::surroundingNodes ( int  dir)
inlinenoexcept

Convert to NODE type in given direction.

◆ TheUnitBox()

template<int dim>
__host__ static __device__ BoxND amrex::BoxND< dim >::TheUnitBox ( )
inlinestaticnoexcept

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

◆ type() [1/2]

template<int dim>
__host__ __device__ IntVectND< dim > amrex::BoxND< dim >::type ( ) const
inlinenoexcept

Returns the indexing type.

◆ type() [2/2]

template<int dim>
__host__ __device__ IndexType::CellIndex amrex::BoxND< dim >::type ( int  dir) const
inlinenoexcept

Returns the indexing type in the specified direction.

◆ volume()

template<int dim>
__host__ __device__ Long amrex::BoxND< dim >::volume ( ) const
inlinenoexcept

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

Friends And Related Symbol Documentation

◆ BoxCommHelper

template<int dim>
friend class BoxCommHelper
friend

◆ ParallelDescriptor::Mpi_typemap

template<int dim>
MPI_Datatype ParallelDescriptor::Mpi_typemap ( )
friend

Member Data Documentation

◆ bigend

template<int dim>
IntVectND<dim> amrex::BoxND< dim >::bigend
private

◆ btype

template<int dim>
IndexTypeND<dim> amrex::BoxND< dim >::btype
private

◆ smallend

template<int dim>
IntVectND<dim> amrex::BoxND< dim >::smallend
private

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