Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
amrex::AuxBoundaryData Class Reference

#include <AMReX_AuxBoundaryData.H>

Public Member Functions

 AuxBoundaryData () noexcept=default
 Construct an empty container; call initialize() before use.
 
 AuxBoundaryData (const BoxArray &ba, int n_grow, int n_comp, const Geometry &geom)
 Construct and allocate storage for boundary data.
 
 ~AuxBoundaryData ()=default
 Trivial destructor; MultiFab members handle their own cleanup.
 
 AuxBoundaryData (AuxBoundaryData &&rhs)=default
 Move constructor keeps the bookkeeping simple when containers relocate.
 
AuxBoundaryDataoperator= (AuxBoundaryData &&rhs)=default
 Move assignment likewise transfers ownership of the cached MultiFab.
 
 AuxBoundaryData (const AuxBoundaryData &rhs)
 Deep-copy constructor.
 
AuxBoundaryDataoperator= (const AuxBoundaryData &rhs)=delete
 Copy assignment is disabled; construct a new container instead.
 
void copyTo (MultiFab &destmf, int src_comp, int dst_comp, int num_comp) const
 Copy data from the auxiliary storage into destmf.
 
void copyFrom (const MultiFab &srcmf, int src_comp, int dst_comp, int num_comp, int src_ng=0)
 Populate the auxiliary storage from srcmf.
 
size_t size () const noexcept
 Number of FArrayBoxes currently cached (0 if empty/uninitialized).
 
void copy (const AuxBoundaryData &src, int src_comp, int dst_comp, int num_comp)
 Copy data from another AuxBoundaryData src (component-wise).
 
void initialize (const BoxArray &ba, int n_grow, int n_comp, const Geometry &geom)
 Allocate and configure storage.
 
const BoxArrayequivBoxArray () const noexcept
 
void setVal (Real r)
 Fill all stored data with scalar value r (if initialized).
 
const DistributionMappingDistributionMap () const noexcept
 Distribution mapping associated with the cached grids.
 
FArrayBoxoperator[] (const MFIter &mfi) noexcept
 Mutable access to the boundary data at given mfi.
 
const FArrayBoxoperator[] (const MFIter &mfi) const noexcept
 Read-only access to the boundary data at given mfi.
 
int nGrow () const noexcept
 Number of ghost cells captured when the container was initialized.
 
int nComp () const noexcept
 Number of components stored per entry.
 
bool isEmpty () const noexcept
 True if the container manages no data (either default constructed or after reset()).
 

Constructor & Destructor Documentation

◆ AuxBoundaryData() [1/4]

amrex::AuxBoundaryData::AuxBoundaryData ( )
defaultnoexcept

Construct an empty container; call initialize() before use.

◆ AuxBoundaryData() [2/4]

amrex::AuxBoundaryData::AuxBoundaryData ( const BoxArray ba,
int  n_grow,
int  n_comp,
const Geometry geom 
)

Construct and allocate storage for boundary data.

Parameters
baBoxArray describing the owning level.
n_growNumber of grow cells to capture around each box.
n_compNumber of components per FArrayBox.
geomGeometry describing the level (used for periodicity).

◆ ~AuxBoundaryData()

amrex::AuxBoundaryData::~AuxBoundaryData ( )
default

Trivial destructor; MultiFab members handle their own cleanup.

◆ AuxBoundaryData() [3/4]

amrex::AuxBoundaryData::AuxBoundaryData ( AuxBoundaryData &&  rhs)
default

Move constructor keeps the bookkeeping simple when containers relocate.

◆ AuxBoundaryData() [4/4]

amrex::AuxBoundaryData::AuxBoundaryData ( const AuxBoundaryData rhs)

Deep-copy constructor.

Parameters
rhsSource AuxBoundaryData whose storage is duplicated.

Member Function Documentation

◆ copy()

void amrex::AuxBoundaryData::copy ( const AuxBoundaryData src,
int  src_comp,
int  dst_comp,
int  num_comp 
)

Copy data from another AuxBoundaryData src (component-wise).

Parameters
srcSource container.
src_compStarting component inside src.
dst_compDestination component inside this container.
num_compNumber of components to copy.

◆ copyFrom()

void amrex::AuxBoundaryData::copyFrom ( const MultiFab srcmf,
int  src_comp,
int  dst_comp,
int  num_comp,
int  src_ng = 0 
)

Populate the auxiliary storage from srcmf.

Parameters
srcmfSource MultiFab supplying data.
src_compStarting component inside srcmf.
dst_compDestination component inside this container.
num_compNumber of components to copy.
src_ngNumber of valid grow cells present in srcmf.

◆ copyTo()

void amrex::AuxBoundaryData::copyTo ( MultiFab destmf,
int  src_comp,
int  dst_comp,
int  num_comp 
) const

Copy data from the auxiliary storage into destmf.

Parameters
destmfDestination MultiFab.
src_compStarting component inside this container.
dst_compDestination component index.
num_compNumber of components to copy.

◆ DistributionMap()

const DistributionMapping & amrex::AuxBoundaryData::DistributionMap ( ) const
inlinenoexcept

Distribution mapping associated with the cached grids.

◆ equivBoxArray()

const BoxArray & amrex::AuxBoundaryData::equivBoxArray ( ) const
inlinenoexcept
Returns
BoxArray describing the cached boundary regions.

◆ initialize()

void amrex::AuxBoundaryData::initialize ( const BoxArray ba,
int  n_grow,
int  n_comp,
const Geometry geom 
)

Allocate and configure storage.

Parameters
baBoxArray describing the owning level.
n_growNumber of grow cells to capture.
n_compNumber of components per face.
geomGeometry describing the level.

◆ isEmpty()

bool amrex::AuxBoundaryData::isEmpty ( ) const
inlinenoexcept

True if the container manages no data (either default constructed or after reset()).

◆ nComp()

int amrex::AuxBoundaryData::nComp ( ) const
inlinenoexcept

Number of components stored per entry.

◆ nGrow()

int amrex::AuxBoundaryData::nGrow ( ) const
inlinenoexcept

Number of ghost cells captured when the container was initialized.

◆ operator=() [1/2]

AuxBoundaryData & amrex::AuxBoundaryData::operator= ( AuxBoundaryData &&  rhs)
default

Move assignment likewise transfers ownership of the cached MultiFab.

◆ operator=() [2/2]

AuxBoundaryData & amrex::AuxBoundaryData::operator= ( const AuxBoundaryData rhs)
delete

Copy assignment is disabled; construct a new container instead.

◆ operator[]() [1/2]

const FArrayBox & amrex::AuxBoundaryData::operator[] ( const MFIter mfi) const
inlinenoexcept

Read-only access to the boundary data at given mfi.

◆ operator[]() [2/2]

FArrayBox & amrex::AuxBoundaryData::operator[] ( const MFIter mfi)
inlinenoexcept

Mutable access to the boundary data at given mfi.

◆ setVal()

void amrex::AuxBoundaryData::setVal ( Real  r)
inline

Fill all stored data with scalar value r (if initialized).

◆ size()

size_t amrex::AuxBoundaryData::size ( ) const
inlinenoexcept

Number of FArrayBoxes currently cached (0 if empty/uninitialized).


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