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

A collection of Boxes stored in an Array. More...

#include <AMReX_BoxArray.H>

Classes

struct  RefID
 

Public Member Functions

 BoxArray () noexcept
 
 BoxArray (const BoxArray &rhs)=default
 
 BoxArray (BoxArray &&rhs) noexcept=default
 
BoxArrayoperator= (BoxArray const &rhs)=default
 
BoxArrayoperator= (BoxArray &&rhs) noexcept=default
 
 ~BoxArray () noexcept=default
 
 BoxArray (const Box &bx)
 Make a boxarray out of a single box. More...
 
 BoxArray (size_t n)
 Construct a BoxArray of the specified size. More...
 
 BoxArray (const Box *bxvec, int nbox)
 Construct a BoxArray from an array of Boxes of size nbox. More...
 
 BoxArray (const BoxList &bl)
 Construct a BoxArray from a BoxList. More...
 
 BoxArray (BoxList &&bl) noexcept
 
 BoxArray (const BoxArray &rhs, const BATransformer &trans)
 
 BoxArray (BoxList &&bl, IntVect const &max_grid_size)
 
void define (const Box &bx)
 Initialize the BoxArray from a single box. It is an error if the BoxArray has already been initialized. More...
 
void define (const BoxList &bl)
 Initialize the BoxArray from the supplied BoxList. It is an error if the BoxArray has already been initialized. More...
 
void define (BoxList &&bl) noexcept
 
void clear ()
 Remove all Boxes from the BoxArray. More...
 
void resize (Long len)
 Resize the BoxArray. See Vector<T>::resize() for the gory details. More...
 
Long size () const noexcept
 Return the number of boxes in the BoxArray. More...
 
Long capacity () const noexcept
 Return the number of boxes that can be held in the current allocated storage. More...
 
bool empty () const noexcept
 Return whether the BoxArray is empty. More...
 
Long numPts () const noexcept
 Returns the total number of cells contained in all boxes in the BoxArray. More...
 
double d_numPts () const noexcept
 Returns the total number of cells (in double type) contained in all boxes in the BoxArray. More...
 
int readFrom (std::istream &is)
 Initialize the BoxArray from the supplied istream. It is an error if the BoxArray has already been initialized. Note that the BoxArray in the istream must have been written using writeOn(). More...
 
std::ostream & writeOn (std::ostream &) const
 Output this BoxArray to a checkpoint file. More...
 
bool operator== (const BoxArray &rhs) const noexcept
 Are the BoxArrays equal? More...
 
bool operator!= (const BoxArray &rhs) const noexcept
 Are the BoxArrays not equal? More...
 
bool operator== (const Vector< Box > &bv) const noexcept
 
bool operator!= (const Vector< Box > &bv) const noexcept
 
bool CellEqual (const BoxArray &rhs) const noexcept
 Are the BoxArrays equal after conversion to cell-centered. More...
 
BoxArraymaxSize (int block_size)
 Forces each Box in BoxArray to have sides <= block_size. More...
 
BoxArraymaxSize (const IntVect &block_size)
 
BoxArrayminmaxSize (const IntVect &min_size, const IntVect &max_size)
 
BoxArrayrefine (int refinement_ratio)
 Refine each Box in the BoxArray to the specified ratio. More...
 
BoxArrayrefine (const IntVect &iv)
 Refine each Box in the BoxArray to the specified ratio. More...
 
BoxArraycoarsen (int refinement_ratio)
 Coarsen each Box in the BoxArray to the specified ratio. More...
 
bool coarsenable (int refinement_ratio, int min_width=1) const
 Coarsen each Box in the BoxArray to the specified ratio. More...
 
bool coarsenable (const IntVect &refinement_ratio, int min_width=1) const
 
bool coarsenable (const IntVect &refinement_ratio, const IntVect &min_width) const
 
BoxArraycoarsen (const IntVect &iv)
 Coarsen each Box in the BoxArray to the specified ratio. More...
 
BoxArraygrowcoarsen (int n, const IntVect &iv)
 Grow and then coarsen each Box in the BoxArray. More...
 
BoxArraygrowcoarsen (IntVect const &ngrow, const IntVect &iv)
 
BoxArraygrow (int n)
 Grow each Box in the BoxArray by the specified amount. More...
 
BoxArraygrow (const IntVect &iv)
 Grow each Box in the BoxArray by the specified amount. More...
 
BoxArraygrow (int idir, int n_cell)
 Grow each Box in the BoxArray on the low and high ends by n_cell cells in the idir direction. More...
 
BoxArraygrowLo (int idir, int n_cell)
 Grow each Box in the BoxArray on the low end by n_cell cells in the idir direction. More...
 
BoxArraygrowHi (int idir, int n_cell)
 Grow each Box in the BoxArray on the high end by n_cell cells in the idir direction. More...
 
BoxArraysurroundingNodes ()
 Apply surroundingNodes(Box) to each Box in BoxArray. See the documentation of Box for details. More...
 
BoxArraysurroundingNodes (int dir)
 Apply surroundingNodes(Box,int) to each Box in BoxArray. See the documentation of Box for details. More...
 
BoxArrayenclosedCells ()
 Apply Box::enclosedCells() to each Box in the BoxArray. More...
 
BoxArrayenclosedCells (int dir)
 Apply Box::enclosedCells(int) to each Box in the BoxArray. More...
 
BoxArrayconvert (IndexType typ)
 Apply Box::convert(IndexType) to each Box in the BoxArray. More...
 
BoxArrayconvert (const IntVect &iv)
 
BoxArrayconvert (Box(*fp)(const Box &))
 Apply function (*fp)(Box) to each Box in the BoxArray. More...
 
BoxArrayshift (int dir, int nzones)
 Apply Box::shift(int,int) to each Box in the BoxArray. More...
 
BoxArrayshift (const IntVect &iv)
 Apply Box::shift(const IntVect &iv) to each Box in the BoxArray. More...
 
void set (int i, const Box &ibox)
 Set element i in this BoxArray to Box ibox. More...
 
Box operator[] (int index) const noexcept
 Return element index of this BoxArray. More...
 
Box operator[] (const MFIter &mfi) const noexcept
 Return element index of this BoxArray. More...
 
Box get (int index) const noexcept
 Return element index of this BoxArray. More...
 
Box getCellCenteredBox (int index) const noexcept
 Return cell-centered box at element index of this BoxArray. More...
 
bool ok () const
 Return true if Box is valid and they all have the same IndexType. Is true by default if the BoxArray is empty. More...
 
bool isDisjoint () const
 Return true if set of intersecting Boxes in BoxArray is null. More...
 
BoxList boxList () const
 Create a BoxList from this BoxArray. More...
 
bool contains (const IntVect &v) const
 True if the IntVect is within any of the Boxes in this BoxArray. More...
 
bool contains (const Box &b, bool assume_disjoint_ba=false, const IntVect &ng=IntVect(0)) const
 True if the Box is contained in this BoxArray(+ng). The Box must also have the same IndexType as those in this BoxArray. More...
 
bool contains (const BoxArray &ba, bool assume_disjoint_ba=false, const IntVect &ng=IntVect(0)) const
 True if all Boxes in ba are contained in this BoxArray(+ng). More...
 
bool contains (const BoxArray &ba, Periodicity const &period) const
 True if all cells in ba are periodically contained in this BoxArray. More...
 
Box minimalBox () const
 Return smallest Box that contains all Boxes in this BoxArray. More...
 
Box minimalBox (Long &npts_avg_box) const
 
bool intersects (const Box &b, int ng=0) const
 True if the Box intersects with this BoxArray(+ghostcells). The Box must have the same IndexType as those in this BoxArray. More...
 
bool intersects (const Box &b, const IntVect &ng) const
 
std::vector< std::pair< int, Box > > intersections (const Box &bx) const
 Return intersections of Box and BoxArray. More...
 
std::vector< std::pair< int, Box > > intersections (const Box &bx, bool first_only, int ng) const
 Return intersections of Box and BoxArray(+ghostcells). More...
 
std::vector< std::pair< int, Box > > intersections (const Box &bx, bool first_only, const IntVect &ng) const
 
void intersections (const Box &bx, std::vector< std::pair< int, Box > > &isects) const
 intersect Box and BoxArray, then store the result in isects More...
 
void intersections (const Box &bx, std::vector< std::pair< int, Box > > &isects, bool first_only, int ng) const
 intersect Box and BoxArray(+ghostcells), then store the result in isects More...
 
void intersections (const Box &bx, std::vector< std::pair< int, Box > > &isects, bool first_only, const IntVect &ng) const
 
BoxList complementIn (const Box &b) const
 Return box - boxarray. More...
 
void complementIn (BoxList &bl, const Box &b) const
 
void clear_hash_bin () const
 Clear out the internal hash table used by intersections. More...
 
void removeOverlap (bool simplify=true)
 Change the BoxArray to one with no overlap and then simplify it (see the simplify function in BoxList). More...
 
RefID getRefID () const noexcept
 Return a unique ID of the reference. More...
 
IndexType ixType () const noexcept
 Return index type of this BoxArray. More...
 
IntVect crseRatio () const noexcept
 Return crse ratio of this BoxArray. More...
 
void uniqify ()
 Make ourselves unique. More...
 
BoxList const & simplified_list () const
 
BoxArray simplified () const
 
BATransformer const & transformer () const
 
std::weak_ptr< BARefgetWeakRef () const
 
std::shared_ptr< BARef > const & getSharedRef () const
 
std::shared_ptr< BARef > & getSharedRef ()
 

Static Public Member Functions

static bool SameRefs (const BoxArray &lhs, const BoxArray &rhs)
 whether two BoxArrays share the same data More...
 
static void Initialize ()
 
static void Finalize ()
 

Static Public Attributes

static bool initialized
 

Private Member Functions

void type_update ()
 Update BoxArray index type according the box type, and then convert boxes to cell-centered. More...
 
BARef::HashTypegetHashMap () const
 
IntVect getDoiLo () const noexcept
 
IntVect getDoiHi () const noexcept
 

Private Attributes

BATransformer m_bat
 
std::shared_ptr< BARefm_ref
 The data – a reference-counted pointer to a Ref. More...
 
std::shared_ptr< BoxListm_simplified_list
 

Friends

class AmrMesh
 
class FabArrayBase
 

Detailed Description

A collection of Boxes stored in an Array.

It is a reference-counted concrete class, not a polymorphic one; i.e. you cannot use any of the List member functions with a BoxList.

Constructor & Destructor Documentation

◆ BoxArray() [1/10]

amrex::BoxArray::BoxArray ( )
noexcept

◆ BoxArray() [2/10]

amrex::BoxArray::BoxArray ( const BoxArray rhs)
default

◆ BoxArray() [3/10]

amrex::BoxArray::BoxArray ( BoxArray &&  rhs)
defaultnoexcept

◆ ~BoxArray()

amrex::BoxArray::~BoxArray ( )
defaultnoexcept

◆ BoxArray() [4/10]

amrex::BoxArray::BoxArray ( const Box bx)
explicit

Make a boxarray out of a single box.

◆ BoxArray() [5/10]

amrex::BoxArray::BoxArray ( size_t  n)
explicit

Construct a BoxArray of the specified size.

◆ BoxArray() [6/10]

amrex::BoxArray::BoxArray ( const Box bxvec,
int  nbox 
)

Construct a BoxArray from an array of Boxes of size nbox.

◆ BoxArray() [7/10]

amrex::BoxArray::BoxArray ( const BoxList bl)
explicit

Construct a BoxArray from a BoxList.

◆ BoxArray() [8/10]

amrex::BoxArray::BoxArray ( BoxList &&  bl)
explicitnoexcept

◆ BoxArray() [9/10]

amrex::BoxArray::BoxArray ( const BoxArray rhs,
const BATransformer trans 
)

◆ BoxArray() [10/10]

amrex::BoxArray::BoxArray ( BoxList &&  bl,
IntVect const &  max_grid_size 
)

Member Function Documentation

◆ boxList()

BoxList amrex::BoxArray::boxList ( ) const

Create a BoxList from this BoxArray.

◆ capacity()

Long amrex::BoxArray::capacity ( ) const
inlinenoexcept

Return the number of boxes that can be held in the current allocated storage.

◆ CellEqual()

bool amrex::BoxArray::CellEqual ( const BoxArray rhs) const
noexcept

Are the BoxArrays equal after conversion to cell-centered.

◆ clear()

void amrex::BoxArray::clear ( )

Remove all Boxes from the BoxArray.

◆ clear_hash_bin()

void amrex::BoxArray::clear_hash_bin ( ) const

Clear out the internal hash table used by intersections.

◆ coarsen() [1/2]

BoxArray& amrex::BoxArray::coarsen ( const IntVect iv)

Coarsen each Box in the BoxArray to the specified ratio.

◆ coarsen() [2/2]

BoxArray& amrex::BoxArray::coarsen ( int  refinement_ratio)

Coarsen each Box in the BoxArray to the specified ratio.

◆ coarsenable() [1/3]

bool amrex::BoxArray::coarsenable ( const IntVect refinement_ratio,
const IntVect min_width 
) const

◆ coarsenable() [2/3]

bool amrex::BoxArray::coarsenable ( const IntVect refinement_ratio,
int  min_width = 1 
) const

◆ coarsenable() [3/3]

bool amrex::BoxArray::coarsenable ( int  refinement_ratio,
int  min_width = 1 
) const

Coarsen each Box in the BoxArray to the specified ratio.

◆ complementIn() [1/2]

void amrex::BoxArray::complementIn ( BoxList bl,
const Box b 
) const

◆ complementIn() [2/2]

BoxList amrex::BoxArray::complementIn ( const Box b) const

Return box - boxarray.

◆ contains() [1/4]

bool amrex::BoxArray::contains ( const Box b,
bool  assume_disjoint_ba = false,
const IntVect ng = IntVect(0) 
) const

True if the Box is contained in this BoxArray(+ng). The Box must also have the same IndexType as those in this BoxArray.

◆ contains() [2/4]

bool amrex::BoxArray::contains ( const BoxArray ba,
bool  assume_disjoint_ba = false,
const IntVect ng = IntVect(0) 
) const

True if all Boxes in ba are contained in this BoxArray(+ng).

◆ contains() [3/4]

bool amrex::BoxArray::contains ( const BoxArray ba,
Periodicity const &  period 
) const

True if all cells in ba are periodically contained in this BoxArray.

If a cell after being periodically shifted is contained in this BoxArray, it's considered being periodically contained.

◆ contains() [4/4]

bool amrex::BoxArray::contains ( const IntVect v) const

True if the IntVect is within any of the Boxes in this BoxArray.

◆ convert() [1/3]

BoxArray& amrex::BoxArray::convert ( Box(*)(const Box &)  fp)

Apply function (*fp)(Box) to each Box in the BoxArray.

◆ convert() [2/3]

BoxArray& amrex::BoxArray::convert ( const IntVect iv)

◆ convert() [3/3]

BoxArray& amrex::BoxArray::convert ( IndexType  typ)

Apply Box::convert(IndexType) to each Box in the BoxArray.

◆ crseRatio()

IntVect amrex::BoxArray::crseRatio ( ) const
inlinenoexcept

Return crse ratio of this BoxArray.

◆ d_numPts()

double amrex::BoxArray::d_numPts ( ) const
noexcept

Returns the total number of cells (in double type) contained in all boxes in the BoxArray.

◆ define() [1/3]

void amrex::BoxArray::define ( BoxList &&  bl)
noexcept

◆ define() [2/3]

void amrex::BoxArray::define ( const Box bx)

Initialize the BoxArray from a single box. It is an error if the BoxArray has already been initialized.

◆ define() [3/3]

void amrex::BoxArray::define ( const BoxList bl)

Initialize the BoxArray from the supplied BoxList. It is an error if the BoxArray has already been initialized.

◆ empty()

bool amrex::BoxArray::empty ( ) const
inlinenoexcept

Return whether the BoxArray is empty.

◆ enclosedCells() [1/2]

BoxArray& amrex::BoxArray::enclosedCells ( )

Apply Box::enclosedCells() to each Box in the BoxArray.

◆ enclosedCells() [2/2]

BoxArray& amrex::BoxArray::enclosedCells ( int  dir)

Apply Box::enclosedCells(int) to each Box in the BoxArray.

◆ Finalize()

static void amrex::BoxArray::Finalize ( )
static

◆ get()

Box amrex::BoxArray::get ( int  index) const
inlinenoexcept

Return element index of this BoxArray.

◆ getCellCenteredBox()

Box amrex::BoxArray::getCellCenteredBox ( int  index) const
inlinenoexcept

Return cell-centered box at element index of this BoxArray.

◆ getDoiHi()

IntVect amrex::BoxArray::getDoiHi ( ) const
privatenoexcept

◆ getDoiLo()

IntVect amrex::BoxArray::getDoiLo ( ) const
privatenoexcept

◆ getHashMap()

BARef::HashType& amrex::BoxArray::getHashMap ( ) const
private

◆ getRefID()

RefID amrex::BoxArray::getRefID ( ) const
inlinenoexcept

Return a unique ID of the reference.

◆ getSharedRef() [1/2]

std::shared_ptr<BARef>& amrex::BoxArray::getSharedRef ( )

◆ getSharedRef() [2/2]

std::shared_ptr<BARef> const& amrex::BoxArray::getSharedRef ( ) const

◆ getWeakRef()

std::weak_ptr<BARef> amrex::BoxArray::getWeakRef ( ) const

◆ grow() [1/3]

BoxArray& amrex::BoxArray::grow ( const IntVect iv)

Grow each Box in the BoxArray by the specified amount.

◆ grow() [2/3]

BoxArray& amrex::BoxArray::grow ( int  idir,
int  n_cell 
)

Grow each Box in the BoxArray on the low and high ends by n_cell cells in the idir direction.

◆ grow() [3/3]

BoxArray& amrex::BoxArray::grow ( int  n)

Grow each Box in the BoxArray by the specified amount.

◆ growcoarsen() [1/2]

BoxArray& amrex::BoxArray::growcoarsen ( int  n,
const IntVect iv 
)

Grow and then coarsen each Box in the BoxArray.

◆ growcoarsen() [2/2]

BoxArray& amrex::BoxArray::growcoarsen ( IntVect const &  ngrow,
const IntVect iv 
)

◆ growHi()

BoxArray& amrex::BoxArray::growHi ( int  idir,
int  n_cell 
)

Grow each Box in the BoxArray on the high end by n_cell cells in the idir direction.

◆ growLo()

BoxArray& amrex::BoxArray::growLo ( int  idir,
int  n_cell 
)

Grow each Box in the BoxArray on the low end by n_cell cells in the idir direction.

◆ Initialize()

static void amrex::BoxArray::Initialize ( )
static

◆ intersections() [1/6]

std::vector< std::pair<int,Box> > amrex::BoxArray::intersections ( const Box bx) const

Return intersections of Box and BoxArray.

◆ intersections() [2/6]

std::vector< std::pair<int,Box> > amrex::BoxArray::intersections ( const Box bx,
bool  first_only,
const IntVect ng 
) const

◆ intersections() [3/6]

std::vector< std::pair<int,Box> > amrex::BoxArray::intersections ( const Box bx,
bool  first_only,
int  ng 
) const

Return intersections of Box and BoxArray(+ghostcells).

◆ intersections() [4/6]

void amrex::BoxArray::intersections ( const Box bx,
std::vector< std::pair< int, Box > > &  isects 
) const

intersect Box and BoxArray, then store the result in isects

◆ intersections() [5/6]

void amrex::BoxArray::intersections ( const Box bx,
std::vector< std::pair< int, Box > > &  isects,
bool  first_only,
const IntVect ng 
) const

◆ intersections() [6/6]

void amrex::BoxArray::intersections ( const Box bx,
std::vector< std::pair< int, Box > > &  isects,
bool  first_only,
int  ng 
) const

intersect Box and BoxArray(+ghostcells), then store the result in isects

◆ intersects() [1/2]

bool amrex::BoxArray::intersects ( const Box b,
const IntVect ng 
) const

◆ intersects() [2/2]

bool amrex::BoxArray::intersects ( const Box b,
int  ng = 0 
) const

True if the Box intersects with this BoxArray(+ghostcells). The Box must have the same IndexType as those in this BoxArray.

◆ isDisjoint()

bool amrex::BoxArray::isDisjoint ( ) const

Return true if set of intersecting Boxes in BoxArray is null.

◆ ixType()

IndexType amrex::BoxArray::ixType ( ) const
inlinenoexcept

Return index type of this BoxArray.

◆ maxSize() [1/2]

BoxArray& amrex::BoxArray::maxSize ( const IntVect block_size)

◆ maxSize() [2/2]

BoxArray& amrex::BoxArray::maxSize ( int  block_size)

Forces each Box in BoxArray to have sides <= block_size.

◆ minimalBox() [1/2]

Box amrex::BoxArray::minimalBox ( ) const

Return smallest Box that contains all Boxes in this BoxArray.

◆ minimalBox() [2/2]

Box amrex::BoxArray::minimalBox ( Long &  npts_avg_box) const

◆ minmaxSize()

BoxArray& amrex::BoxArray::minmaxSize ( const IntVect min_size,
const IntVect max_size 
)

Forces each Box in BoxArray to have sizes >= min_size and <= max_size. It's the caller's responsibility to make sure both the BoxArray and max_size are coarsenable by min_size.

◆ numPts()

Long amrex::BoxArray::numPts ( ) const
noexcept

Returns the total number of cells contained in all boxes in the BoxArray.

◆ ok()

bool amrex::BoxArray::ok ( ) const

Return true if Box is valid and they all have the same IndexType. Is true by default if the BoxArray is empty.

◆ operator!=() [1/2]

bool amrex::BoxArray::operator!= ( const BoxArray rhs) const
noexcept

Are the BoxArrays not equal?

◆ operator!=() [2/2]

bool amrex::BoxArray::operator!= ( const Vector< Box > &  bv) const
noexcept

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ operator==() [1/2]

bool amrex::BoxArray::operator== ( const BoxArray rhs) const
noexcept

Are the BoxArrays equal?

◆ operator==() [2/2]

bool amrex::BoxArray::operator== ( const Vector< Box > &  bv) const
noexcept

◆ operator[]() [1/2]

Box amrex::BoxArray::operator[] ( const MFIter mfi) const
noexcept

Return element index of this BoxArray.

◆ operator[]() [2/2]

Box amrex::BoxArray::operator[] ( int  index) const
inlinenoexcept

Return element index of this BoxArray.

◆ readFrom()

int amrex::BoxArray::readFrom ( std::istream &  is)

Initialize the BoxArray from the supplied istream. It is an error if the BoxArray has already been initialized. Note that the BoxArray in the istream must have been written using writeOn().

◆ refine() [1/2]

BoxArray& amrex::BoxArray::refine ( const IntVect iv)

Refine each Box in the BoxArray to the specified ratio.

◆ refine() [2/2]

BoxArray& amrex::BoxArray::refine ( int  refinement_ratio)

Refine each Box in the BoxArray to the specified ratio.

◆ removeOverlap()

void amrex::BoxArray::removeOverlap ( bool  simplify = true)

Change the BoxArray to one with no overlap and then simplify it (see the simplify function in BoxList).

◆ resize()

void amrex::BoxArray::resize ( Long  len)

Resize the BoxArray. See Vector<T>::resize() for the gory details.

◆ SameRefs()

static bool amrex::BoxArray::SameRefs ( const BoxArray lhs,
const BoxArray rhs 
)
inlinestatic

whether two BoxArrays share the same data

◆ set()

void amrex::BoxArray::set ( int  i,
const Box ibox 
)

Set element i in this BoxArray to Box ibox.

◆ shift() [1/2]

BoxArray& amrex::BoxArray::shift ( const IntVect iv)

Apply Box::shift(const IntVect &iv) to each Box in the BoxArray.

◆ shift() [2/2]

BoxArray& amrex::BoxArray::shift ( int  dir,
int  nzones 
)

Apply Box::shift(int,int) to each Box in the BoxArray.

◆ simplified()

BoxArray amrex::BoxArray::simplified ( ) const

◆ simplified_list()

BoxList const& amrex::BoxArray::simplified_list ( ) const

◆ size()

Long amrex::BoxArray::size ( ) const
inlinenoexcept

Return the number of boxes in the BoxArray.

◆ surroundingNodes() [1/2]

BoxArray& amrex::BoxArray::surroundingNodes ( )

Apply surroundingNodes(Box) to each Box in BoxArray. See the documentation of Box for details.

◆ surroundingNodes() [2/2]

BoxArray& amrex::BoxArray::surroundingNodes ( int  dir)

Apply surroundingNodes(Box,int) to each Box in BoxArray. See the documentation of Box for details.

◆ transformer()

BATransformer const& amrex::BoxArray::transformer ( ) const

◆ type_update()

void amrex::BoxArray::type_update ( )
private

Update BoxArray index type according the box type, and then convert boxes to cell-centered.

◆ uniqify()

void amrex::BoxArray::uniqify ( )

Make ourselves unique.

◆ writeOn()

std::ostream& amrex::BoxArray::writeOn ( std::ostream &  ) const

Output this BoxArray to a checkpoint file.

Friends And Related Function Documentation

◆ AmrMesh

friend class AmrMesh
friend

◆ FabArrayBase

friend class FabArrayBase
friend

Member Data Documentation

◆ initialized

bool amrex::BoxArray::initialized
static

◆ m_bat

BATransformer amrex::BoxArray::m_bat
private

◆ m_ref

std::shared_ptr<BARef> amrex::BoxArray::m_ref
private

The data – a reference-counted pointer to a Ref.

◆ m_simplified_list

std::shared_ptr<BoxList> amrex::BoxArray::m_simplified_list
mutableprivate

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