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

A collection (stored as an array) of FArrayBox objects. More...

#include <AMReX_MultiFab.H>

Inheritance diagram for amrex::MultiFab:
amrex::FabArray< FArrayBox > amrex::FabArrayBase

Public Member Functions

 MultiFab () noexcept
 Constructs an empty MultiFab. More...
 
 MultiFab (Arena *a) noexcept
 Constructs an empty MultiFab. More...
 
 MultiFab (const BoxArray &bxs, const DistributionMapping &dm, int ncomp, int ngrow, const MFInfo &info=MFInfo(), const FabFactory< FArrayBox > &factory=FArrayBoxFactory())
 Constructs a MultiFab. More...
 
 MultiFab (const BoxArray &bxs, const DistributionMapping &dm, int ncomp, const IntVect &ngrow, const MFInfo &info=MFInfo(), const FabFactory< FArrayBox > &factory=FArrayBoxFactory())
 
 MultiFab (const MultiFab &rhs, MakeType maketype, int scomp, int ncomp)
 Make an alias MultiFab. More...
 
 ~MultiFab ()
 
 MultiFab (MultiFab &&rhs) noexcept
 
MultiFaboperator= (MultiFab &&rhs) noexcept=default
 
 MultiFab (const MultiFab &rhs)=delete
 
MultiFaboperator= (const MultiFab &rhs)=delete
 
void define (const BoxArray &bxs, const DistributionMapping &dm, int nvar, int ngrow, const MFInfo &info=MFInfo(), const FabFactory< FArrayBox > &factory=FArrayBoxFactory())
 
void define (const BoxArray &bxs, const DistributionMapping &dm, int nvar, const IntVect &ngrow, const MFInfo &info=MFInfo(), const FabFactory< FArrayBox > &factory=FArrayBoxFactory())
 
MultiFaboperator= (Real r)
 
Real min (int comp, int nghost=0, bool local=false) const
 Returns the minimum value contained in component comp of the MultiFab. More...
 
Real min (const Box &region, int comp, int nghost=0, bool local=false) const
 Identical to min() function, but confines its search to intersection of Box b and the MultiFab. More...
 
Real max (int comp, int nghost=0, bool local=false) const
 Returns the maximum value contained in component comp of the MultiFab. More...
 
Real max (const Box &region, int comp, int nghost=0, bool local=false) const
 Identical to the previous max() function, but confines its search to intersection of Box b and the MultiFab. More...
 
Real norm0 (int comp=0, int nghost=0, bool local=false, bool ignore_covered=false) const
 Returns the maximum absolute value contained in component comp of the MultiFab. More...
 
Real norminf (int comp=0, int nghost=0, bool local=false, bool ignore_covered=false) const
 
Real norm0 (const iMultiFab &mask, int comp=0, int nghost=0, bool local=false) const
 
Real norminf (const iMultiFab &mask, int comp=0, int nghost=0, bool local=false) const
 
Real norm0 (int comp, int ncomp, IntVect const &nghost, bool local=false, bool ignore_covered=false) const
 
Vector< Real > norm0 (const Vector< int > &comps, int nghost=0, bool local=false, bool ignore_covered=false) const
 Returns the maximum absolute values contained in each component of comps of the MultiFab. nghost ghost cells are used. More...
 
Vector< Real > norminf (const Vector< int > &comps, int nghost=0, bool local=false, bool ignore_covered=false) const
 
Real norm1 (int comp, const Periodicity &period, bool ignore_covered=false) const
 Returns the L1 norm of component comp over the MultiFab. More...
 
Real norm1 (int comp=0, int ngrow=0, bool local=false) const
 Returns the L1 norm of component comp over the MultiFab. ngrow ghost cells are used. More...
 
Vector< Real > norm1 (const Vector< int > &comps, int ngrow=0, bool local=false) const
 Returns the L1 norm of each component of "comps" over the MultiFab. ngrow ghost cells are used. More...
 
Real norm2 (int comp=0) const
 Returns the L2 norm of component comp over the MultiFab. No ghost cells are used. More...
 
Real norm2 (int comp, const Periodicity &period) const
 Returns the L2 norm of component comp over the MultiFab. No ghost cells are used. This version has no double counting for nodal data. More...
 
Vector< Real > norm2 (const Vector< int > &comps) const
 Returns the L2 norm of each component of "comps" over the MultiFab. No ghost cells are used. More...
 
Real sum (int comp=0, bool local=false) const
 Returns the sum of component "comp" over the MultiFab – no ghost cells are included. More...
 
Real sum (Box const &region, int comp=0, bool local=false) const
 Returns the sum of component "comp" in the given "region". – no ghost cells are included. More...
 
Real sum_unique (int comp=0, bool local=false, const Periodicity &period=Periodicity::NonPeriodic()) const
 Same as sum with local =false, but for non-cell-centered data, this only adds non-unique points that are owned by multiple boxes once. More...
 
Real sum_unique (Box const &region, int comp=0, bool local=false) const
 Returns the unique sum of component "comp" in the given region. Non-unique points owned by multiple boxes in the MultiFab are only added once. No ghost cells are included. This function does not take periodicity into account in the determination of uniqueness of points. More...
 
void plus (Real val, int comp, int num_comp, int nghost=0)
 Adds the scalar value val to the value of each cell in the specified subregion of the MultiFab. More...
 
void plus (Real val, const Box &region, int comp, int num_comp, int nghost=0)
 Identical to the previous version of plus(), with the restriction that the subregion is further constrained to the intersection with Box region. More...
 
void plus (Real val, int nghost)
 Adds the scalar value val to the value of each cell in the valid region of each component of the MultiFab. The value of nghost specifies the number of cells in the boundary region that should be modified. More...
 
void plus (Real val, const Box &region, int nghost)
 Adds the scalar value val to the value of each cell in the valid region of each component of the MultiFab, that also intersects the Box region. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified. More...
 
void mult (Real val, int comp, int num_comp, int nghost=0)
 Scales the value of each cell in the specified subregion of the MultiFab by the scalar val (a[i] <- a[i]*val). The subregion consists of the num_comp components starting at component comp. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified. More...
 
void mult (Real val, const Box &region, int comp, int num_comp, int nghost=0)
 Identical to the previous version of mult(), with the restriction that the subregion is further constrained to the intersection with Box region. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified. More...
 
void mult (Real val, int nghost=0)
 Scales the value of each cell in the valid region of each component of the MultiFab by the scalar val (a[i] <- a[i]*val). The value of nghost specifies the number of cells in the boundary region that should be modified. More...
 
void mult (Real val, const Box &region, int nghost=0)
 Scales the value of each cell in the valid region of each component of the MultiFab by the scalar val (a[i] <- a[i]*val), that also intersects the Box region. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified. More...
 
void invert (Real numerator, int comp, int num_comp, int nghost=0)
 Replaces the value of each cell in the specified subregion of the MultiFab with its reciprocal multiplied by the value of numerator. The subregion consists of the num_comp components starting at component comp. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified. More...
 
void invert (Real numerator, const Box &region, int comp, int num_comp, int nghost=0)
 Identical to the previous version of invert(), with the restriction that the subregion is further constrained to the intersection with Box region. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified. More...
 
void invert (Real numerator, int nghost)
 Replaces the value of each cell in the specified subregion of the MultiFab with its reciprocal multiplied by the value of numerator. The value of nghost specifies the number of cells in the boundary region that should be modified. More...
 
void invert (Real numerator, const Box &region, int nghost)
 Replaces the value of each cell in the specified subregion of the MultiFab, that also intersects the Box region, with its reciprocal multiplied by the value of numerator. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified. More...
 
void negate (int comp, int num_comp, int nghost=0)
 Negates the value of each cell in the specified subregion of the MultiFab. The subregion consists of the num_comp components starting at component comp. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified. More...
 
void negate (const Box &region, int comp, int num_comp, int nghost=0)
 Identical to the previous version of negate(), with the restriction that the subregion is further constrained to the intersection with Box region. More...
 
void negate (int nghost=0)
 Negates the value of each cell in the valid region of the MultiFab. The value of nghost specifies the number of cells in the boundary region that should be modified. More...
 
void negate (const Box &region, int nghost=0)
 Negates the value of each cell in the valid region of the MultiFab that also intersects the Box region. The value of nghost specifies the number of cells in the boundary region that should be modified. More...
 
IntVect minIndex (int comp, int nghost=0) const
 
IntVect maxIndex (int comp, int nghost=0) const
 
void plus (const MultiFab &mf, int strt_comp, int num_comp, int nghost)
 This function adds the values of the cells in mf to the corresponding cells of this MultiFab. mf is required to have the same BoxArray or "valid region" as this MultiFab. The addition is done only to num_comp components, starting with component number strt_comp. The parameter nghost specifies the number of boundary cells that will be modified. If nghost == 0, only the valid region of each FArrayBox will be modified. More...
 
void minus (const MultiFab &mf, int strt_comp, int num_comp, int nghost)
 This function subtracts the values of the cells in mf from the corresponding cells of this MultiFab. mf is required to have the same BoxArray or "valid region" as this MultiFab. The subtraction is done only to num_comp components, starting with component number strt_comp. The parameter nghost specifies the number of boundary cells that will be modified. If nghost == 0, only the valid region of each FArrayBox will be modified. More...
 
void divide (const MultiFab &mf, int strt_comp, int num_comp, int nghost)
 This function divides the values of the cells in mf from the corresponding cells of this MultiFab. mf is required to have the same BoxArray or "valid region" as this MultiFab. The division is done only to num_comp components, starting with component number strt_comp. The parameter nghost specifies the number of boundary cells that will be modified. If nghost == 0, only the valid region of each FArrayBox will be modified. Note, nothing is done to protect against divide by zero. More...
 
MultiFab deepCopy () const
 
bool is_finite (bool local=false) const
 Are the numbers in the MF finite (i.e., neither nan nor inf)? This may return true, even if the MF contains NaNs or Infs, if the machine doesn't support the appropriate NaN and Inf testing functions. More...
 
bool is_finite (int scomp, int ncomp, int ngrow=0, bool local=false) const
 
bool is_finite (int scomp, int ncomp, const IntVect &ngrow, bool local=false) const
 
bool contains_nan (bool local=false) const
 Are there any NaNs in the MF? This may return false, even if the MF contains NaNs, if the machine doesn't support the appropriate NaN testing functions. More...
 
bool contains_nan (int scomp, int ncomp, int ngrow=0, bool local=false) const
 
bool contains_nan (int scomp, int ncomp, const IntVect &ngrow, bool local=false) const
 
bool contains_inf (bool local=false) const
 Are there any Infs in the MF? This may return false, even if the MF contains Infs, if the machine doesn't support the appropriate Inf testing functions. This version checks all components and grow cells. More...
 
bool contains_inf (int scomp, int ncomp, int ngrow=0, bool local=false) const
 
bool contains_inf (int scomp, int ncomp, const IntVect &ngrow, bool local=false) const
 
std::unique_ptr< MultiFabOverlapMask (const Periodicity &period=Periodicity::NonPeriodic()) const
 Return a mask indicating how many duplicates are in each point. More...
 
std::unique_ptr< iMultiFabOwnerMask (const Periodicity &period=Periodicity::NonPeriodic()) const
 Owner is the grid with the lowest grid number containing the data. More...
 
void AverageSync (const Periodicity &period=Periodicity::NonPeriodic())
 Sync up nodal data via averaging. More...
 
void WeightedSync (const MultiFab &wgt, const Periodicity &period=Periodicity::NonPeriodic())
 Sync up nodal data with weights. More...
 
void OverrideSync (const iMultiFab &msk, const Periodicity &period=Periodicity::NonPeriodic())
 Sync up nodal data with owners overriding non-owners. More...
 
- Public Member Functions inherited from amrex::FabArray< FArrayBox >
 FabArray () noexcept
 Constructs an empty FabArray<FAB>. More...
 
 FabArray (Arena *a) noexcept
 Construct an empty FabArray<FAB> that has a default Arena. More...
 
 FabArray (const BoxArray &bxs, const DistributionMapping &dm, int nvar, int ngrow, const MFInfo &info=MFInfo(), const FabFactory< FArrayBox > &factory=DefaultFabFactory< FArrayBox >())
 Construct a FabArray<FAB> with a valid region defined by bxs and a region of definition defined by the grow factor ngrow and the number of components nvar. More...
 
 FabArray (const BoxArray &bxs, const DistributionMapping &dm, int nvar, const IntVect &ngrow, const MFInfo &info=MFInfo(), const FabFactory< FArrayBox > &factory=DefaultFabFactory< FArrayBox >())
 
 FabArray (const FabArray< FArrayBox > &rhs, MakeType maketype, int scomp, int ncomp)
 
 FabArray (FabArray< FArrayBox > &&rhs) noexcept
 
 FabArray (const FabArray< FArrayBox > &rhs)=delete
 
 ~FabArray ()
 The destructor – deletes all FABs in the array. More...
 
FabArray< FArrayBox > & operator= (FabArray< FArrayBox > &&rhs) noexcept
 
FabArray< FArrayBox > & operator= (const FabArray< FArrayBox > &rhs)=delete
 
FabArray< FArrayBox > & operator= (value_type val)
 Set all components in the entire region of each FAB to val. More...
 
void define (const BoxArray &bxs, const DistributionMapping &dm, int nvar, int ngrow, const MFInfo &info=MFInfo(), const FabFactory< FArrayBox > &factory=DefaultFabFactory< FArrayBox >())
 Define this FabArray identically to that performed by the constructor having an analogous function signature. This is only valid if this FabArray was defined using the default constructor. More...
 
void define (const BoxArray &bxs, const DistributionMapping &dm, int nvar, const IntVect &ngrow, const MFInfo &info=MFInfo(), const FabFactory< FArrayBox > &factory=DefaultFabFactory< FArrayBox >())
 
const FabFactory< FArrayBox > & Factory () const noexcept
 
Arenaarena () const noexcept
 
const Vector< std::string > & tags () const noexcept
 
bool hasEBFabFactory () const noexcept
 
value_typesingleChunkPtr () noexcept
 
value_type const * singleChunkPtr () const noexcept
 
std::size_t singleChunkSize () const noexcept
 
bool isAllRegular () const noexcept
 
bool ok () const
 Return true if the FabArray is well-defined. That is, the FabArray has a BoxArray and DistributionMapping, the FABs are allocated for each Box in the BoxArray and the sizes of the FABs and the number of components are consistent with the definition of the FabArray. More...
 
bool isDefined () const
 
const FArrayBoxoperator[] (const MFIter &mfi) const noexcept
 Return a constant reference to the FAB associated with mfi. More...
 
FArrayBoxoperator[] (const MFIter &mfi) noexcept
 Returns a reference to the FAB associated mfi. More...
 
const FArrayBoxoperator[] (int K) const noexcept
 Return a constant reference to the FAB associated with the Kth element. More...
 
FArrayBoxoperator[] (int K) noexcept
 Return a reference to the FAB associated with the Kth element. More...
 
const FArrayBoxget (const MFIter &mfi) const noexcept
 Return a constant reference to the FAB associated with mfi. More...
 
FArrayBoxget (const MFIter &mfi) noexcept
 Returns a reference to the FAB associated mfi. More...
 
const FArrayBoxget (int K) const noexcept
 Return a constant reference to the FAB associated with the Kth element. More...
 
FArrayBoxget (int K) noexcept
 Return a reference to the FAB associated with the Kth element. More...
 
FArrayBoxatLocalIdx (int L) noexcept
 Return a reference to the FAB associated with local index L. More...
 
const FArrayBoxatLocalIdx (int L) const noexcept
 
FArrayBoxfabPtr (const MFIter &mfi) noexcept
 Return pointer to FAB. More...
 
FArrayBox const * fabPtr (const MFIter &mfi) const noexcept
 
FArrayBoxfabPtr (int K) noexcept
 
FArrayBox const * fabPtr (int K) const noexcept
 
void prefetchToHost (const MFIter &mfi) const noexcept
 
void prefetchToDevice (const MFIter &mfi) const noexcept
 
Array4< typename FabArray< FArrayBox >::value_type const > array (const MFIter &mfi) const noexcept
 
Array4< typename FabArray< FArrayBox >::value_typearray (const MFIter &mfi) noexcept
 
Array4< typename FabArray< FArrayBox >::value_type const > array (int K) const noexcept
 
Array4< typename FabArray< FArrayBox >::value_typearray (int K) noexcept
 
Array4< typename FabArray< FArrayBox >::value_type const > array (const MFIter &mfi, int start_comp) const noexcept
 
Array4< typename FabArray< FArrayBox >::value_typearray (const MFIter &mfi, int start_comp) noexcept
 
Array4< typename FabArray< FArrayBox >::value_type const > array (int K, int start_comp) const noexcept
 
Array4< typename FabArray< FArrayBox >::value_typearray (int K, int start_comp) noexcept
 
Array4< typename FabArray< FArrayBox >::value_type const > const_array (const MFIter &mfi) const noexcept
 
Array4< typename FabArray< FArrayBox >::value_type const > const_array (int K) const noexcept
 
Array4< typename FabArray< FArrayBox >::value_type const > const_array (const MFIter &mfi, int start_comp) const noexcept
 
Array4< typename FabArray< FArrayBox >::value_type const > const_array (int K, int start_comp) const noexcept
 
MultiArray4< typename FabArray< FArrayBox >::value_typearrays () noexcept
 
MultiArray4< typename FabArray< FArrayBox >::value_type const > arrays () const noexcept
 
MultiArray4< typename FabArray< FArrayBox >::value_type const > const_arrays () const noexcept
 
void setFab (int boxno, std::unique_ptr< FArrayBox > elem)
 Explicitly set the Kth FAB in the FabArray to point to elem. More...
 
void setFab (int boxno, FArrayBox &&elem)
 Explicitly set the Kth FAB in the FabArray to point to elem. More...
 
void setFab (const MFIter &mfi, std::unique_ptr< FArrayBox > elem)
 Explicitly set the FAB associated with mfi in the FabArray to point to elem. More...
 
void setFab (const MFIter &mfi, FArrayBox &&elem)
 Explicitly set the FAB associated with mfi in the FabArray to point to elem. More...
 
AMREX_NODISCARD FArrayBoxrelease (int K)
 Release ownership of the FAB. This function is not thread safe. More...
 
AMREX_NODISCARD FArrayBoxrelease (const MFIter &mfi)
 Release ownership of the FAB. This function is not thread safe. More...
 
void clear ()
 Releases FAB memory in the FabArray. More...
 
void LocalCopy (FabArray< SFAB > const &src, int scomp, int dcomp, int ncomp, IntVect const &nghost)
 Perform local copy of FabArray data. More...
 
void LocalAdd (FabArray< FArrayBox > const &src, int scomp, int dcomp, int ncomp, IntVect const &nghost)
 Perform local addition of FabArray data. More...
 
void setVal (value_type val)
 Set all components in the entire region of each FAB to val. More...
 
void setVal (value_type val, int comp, int ncomp, int nghost=0)
 Set the value of num_comp components in the valid region of each FAB in the FabArray, starting at component comp to val. Also set the value of nghost boundary cells. More...
 
void setVal (value_type val, int comp, int ncomp, const IntVect &nghost)
 
void setVal (value_type val, const Box &region, int comp, int ncomp, int nghost=0)
 Set the value of num_comp components in the valid region of each FAB in the FabArray, starting at component comp, as well as nghost boundary cells, to val, provided they also intersect with the Box region. More...
 
void setVal (value_type val, const Box &region, int comp, int ncomp, const IntVect &nghost)
 
void setVal (value_type val, int nghost)
 Set all components in the valid region of each FAB in the FabArray to val, including nghost boundary cells. More...
 
void setVal (value_type val, const IntVect &nghost)
 
void setVal (value_type val, const Box &region, int nghost)
 Set all components in the valid region of each FAB in the FabArray to val, including nghost boundary cells, that also intersect the Box region. More...
 
void setVal (value_type val, const Box &region, const IntVect &nghost)
 
void setVal (value_type val, const CommMetaData &thecmd, int scomp, int ncomp)
 
void abs (int comp, int ncomp, int nghost=0)
 
void abs (int comp, int ncomp, const IntVect &nghost)
 
void plus (value_type val, int comp, int num_comp, int nghost=0)
 
void plus (value_type val, const Box &region, int comp, int num_comp, int nghost=0)
 
void mult (value_type val, int comp, int num_comp, int nghost=0)
 
void mult (value_type val, const Box &region, int comp, int num_comp, int nghost=0)
 
void invert (value_type numerator, int comp, int num_comp, int nghost=0)
 
void invert (value_type numerator, const Box &region, int comp, int num_comp, int nghost=0)
 
void setBndry (value_type val)
 Set all values in the boundary region to val. More...
 
void setBndry (value_type val, int strt_comp, int ncomp)
 Set ncomp values in the boundary region, starting at start_comp to val. More...
 
void setDomainBndry (value_type val, const Geometry &geom)
 Set all values outside the Geometry domain to val. More...
 
void setDomainBndry (value_type val, int strt_comp, int ncomp, const Geometry &geom)
 Set ncomp values outside the Geometry domain to val, starting at start_comp. More...
 
F::value_type sum (int comp, IntVect const &nghost, bool local=false) const
 Returns the sum of component "comp". More...
 
void ParallelAdd (const FabArray< FArrayBox > &src, const Periodicity &period=Periodicity::NonPeriodic())
 This function copies data from src to this FabArray. Each FAB in fa is intersected with all FABs in this FabArray and a copy is performed on the region of intersection. The intersection is restricted to the valid regions. More...
 
void ParallelAdd (const FabArray< FArrayBox > &src, int src_comp, int dest_comp, int num_comp, const Periodicity &period=Periodicity::NonPeriodic())
 This function copies data from src to this FabArray. Each FAB in src is intersected with all FABs in this FabArray and a copy is performed on the region of intersection. The intersection is restricted to the num_comp components starting at src_comp in the FabArray src, with the destination components in this FabArray starting at dest_comp. More...
 
void ParallelAdd (const FabArray< FArrayBox > &src, int src_comp, int dest_comp, int num_comp, int src_nghost, int dst_nghost, const Periodicity &period=Periodicity::NonPeriodic())
 Similar to the above function, except that source and destination are grown by src_nghost and dst_nghost, respectively. More...
 
void ParallelAdd (const FabArray< FArrayBox > &src, int src_comp, int dest_comp, int num_comp, const IntVect &src_nghost, const IntVect &dst_nghost, const Periodicity &period=Periodicity::NonPeriodic())
 
void ParallelCopy (const FabArray< FArrayBox > &src, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
 
void ParallelCopy (const FabArray< FArrayBox > &src, int src_comp, int dest_comp, int num_comp, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
 
void ParallelCopy (const FabArray< FArrayBox > &src, int src_comp, int dest_comp, int num_comp, int src_nghost, int dst_nghost, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
 
void ParallelCopy (const FabArray< FArrayBox > &src, int scomp, int dcomp, int ncomp, const IntVect &snghost, const IntVect &dnghost, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY, const FabArrayBase::CPC *a_cpc=nullptr)
 
void copy (const FabArray< FArrayBox > &src, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
 
void copy (const FabArray< FArrayBox > &src, int src_comp, int dest_comp, int num_comp, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
 
void copy (const FabArray< FArrayBox > &src, int src_comp, int dest_comp, int num_comp, int src_nghost, int dst_nghost, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
 
void copy (const FabArray< FArrayBox > &src, int src_comp, int dest_comp, int num_comp, const IntVect &src_nghost, const IntVect &dst_nghost, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
 
void ParallelAdd_nowait (const FabArray< FArrayBox > &src, const Periodicity &period=Periodicity::NonPeriodic())
 
void ParallelAdd_nowait (const FabArray< FArrayBox > &src, int src_comp, int dest_comp, int num_comp, const Periodicity &period=Periodicity::NonPeriodic())
 
void ParallelAdd_nowait (const FabArray< FArrayBox > &src, int src_comp, int dest_comp, int num_comp, int src_nghost, int dst_nghost, const Periodicity &period=Periodicity::NonPeriodic())
 
void ParallelAdd_nowait (const FabArray< FArrayBox > &src, int src_comp, int dest_comp, int num_comp, const IntVect &src_nghost, const IntVect &dst_nghost, const Periodicity &period=Periodicity::NonPeriodic())
 
void ParallelCopy_nowait (const FabArray< FArrayBox > &src, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
 
void ParallelCopy_nowait (const FabArray< FArrayBox > &src, int src_comp, int dest_comp, int num_comp, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
 
void ParallelCopy_nowait (const FabArray< FArrayBox > &src, int src_comp, int dest_comp, int num_comp, int src_nghost, int dst_nghost, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
 
void ParallelCopy_nowait (const FabArray< FArrayBox > &src, int scomp, int dcomp, int ncomp, const IntVect &snghost, const IntVect &dnghost, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY, const FabArrayBase::CPC *a_cpc=nullptr, bool to_ghost_cells_only=false)
 
void ParallelCopy_finish ()
 
void ParallelCopyToGhost (const FabArray< FArrayBox > &src, int scomp, int dcomp, int ncomp, const IntVect &snghost, const IntVect &dnghost, const Periodicity &period=Periodicity::NonPeriodic())
 
void ParallelCopyToGhost_nowait (const FabArray< FArrayBox > &src, int scomp, int dcomp, int ncomp, const IntVect &snghost, const IntVect &dnghost, const Periodicity &period=Periodicity::NonPeriodic())
 
void ParallelCopyToGhost_finish ()
 
void Redistribute (const FabArray< FArrayBox > &src, int scomp, int dcomp, int ncomp, const IntVect &nghost)
 Copy from src to this. this and src have the same BoxArray, but different DistributionMapping. More...
 
void copyTo (FArrayBox &dest, int nghost=0) const
 Copy the values contained in the intersection of the valid + nghost region of this FabArray with the FAB dest into dest. Note that FAB dest is assumed to be identical on each process. More...
 
void copyTo (FArrayBox &dest, int scomp, int dcomp, int ncomp, int nghost=0) const
 Copy the values contained in the intersection of the num_comp component valid + nghost region of this FabArray, starting at component src_comp, with the FAB dest into dest, starting at component dest_comp in dest. Note that FAB dest is assumed to be identical on each process. More...
 
void shift (const IntVect &v)
 Shift the boxarray by vector v. More...
 
bool defined (int K) const noexcept
 
bool defined (const MFIter &mfi) const noexcept
 
void FillBoundary (bool cross=false)
 Copy on intersection within a FabArray. Data is copied from valid regions to intersecting regions of definition. The purpose is to fill in the boundary regions of each FAB in the FabArray. If cross=true, corner cells are not filled. If the length of periodic is provided, periodic boundaries are also filled. Note that FabArray itself does not contains any periodicity information. FillBoundary expects that its cell-centered version of its BoxArray is non-overlapping. More...
 
void FillBoundary (const Periodicity &period, bool cross=false)
 
void FillBoundary (const IntVect &nghost, const Periodicity &period, bool cross=false)
 
void FillBoundary (int scomp, int ncomp, bool cross=false)
 Same as FillBoundary(), but only copies ncomp components starting at scomp. More...
 
void FillBoundary (int scomp, int ncomp, const Periodicity &period, bool cross=false)
 
void FillBoundary (int scomp, int ncomp, const IntVect &nghost, const Periodicity &period, bool cross=false)
 
void FillBoundary_nowait (bool cross=false)
 
void FillBoundary_nowait (const Periodicity &period, bool cross=false)
 
void FillBoundary_nowait (const IntVect &nghost, const Periodicity &period, bool cross=false)
 
void FillBoundary_nowait (int scomp, int ncomp, bool cross=false)
 
void FillBoundary_nowait (int scomp, int ncomp, const Periodicity &period, bool cross=false)
 
void FillBoundary_nowait (int scomp, int ncomp, const IntVect &nghost, const Periodicity &period, bool cross=false)
 
void FillBoundary_finish ()
 
void FillBoundary_finish ()
 
void FillBoundary_test ()
 
void FillBoundaryAndSync (const Periodicity &period=Periodicity::NonPeriodic())
 Fill ghost cells and synchronize nodal data. Ghost regions are filled with data from the intersecting valid regions. The synchronization will override valid regions by the intersecting valid regions with a higher precedence. The smaller the global box index is, the higher precedence the box has. With periodic boundaries, for cells in the same box, those near the lower corner have higher precedence than those near the upper corner. More...
 
void FillBoundaryAndSync (int scomp, int ncomp, const IntVect &nghost, const Periodicity &period)
 Fill ghost cells and synchronize nodal data. Ghost regions are filled with data from the intersecting valid regions. The synchronization will override valid regions by the intersecting valid regions with a higher precedence. The smaller the global box index is, the higher precedence the box has. With periodic boundaries, for cells in the same box, those near the lower corner have higher precedence than those near the upper corner. More...
 
void FillBoundaryAndSync_nowait (const Periodicity &period=Periodicity::NonPeriodic())
 
void FillBoundaryAndSync_nowait (int scomp, int ncomp, const IntVect &nghost, const Periodicity &period)
 
void FillBoundaryAndSync_finish ()
 
void OverrideSync (const Periodicity &period=Periodicity::NonPeriodic())
 Synchronize nodal data. The synchronization will override valid regions by the intersecting valid regions with a higher precedence. The smaller the global box index is, the higher precedence the box has. With periodic boundaries, for cells in the same box, those near the lower corner have higher precedence than those near the upper corner. More...
 
void OverrideSync (int scomp, int ncomp, const Periodicity &period)
 Synchronize nodal data. The synchronization will override valid regions by the intersecting valid regions with a higher precedence. The smaller the global box index is, the higher precedence the box has. With periodic boundaries, for cells in the same box, those near the lower corner have higher precedence than those near the upper corner. More...
 
void OverrideSync_nowait (const Periodicity &period=Periodicity::NonPeriodic())
 
void OverrideSync_nowait (int scomp, int ncomp, const Periodicity &period)
 
void OverrideSync_finish ()
 
void SumBoundary (const Periodicity &period=Periodicity::NonPeriodic())
 Sum values in overlapped cells. The destination is limited to valid cells. More...
 
void SumBoundary (int scomp, int ncomp, const Periodicity &period=Periodicity::NonPeriodic())
 
void SumBoundary (int scomp, int ncomp, IntVect const &nghost, const Periodicity &period=Periodicity::NonPeriodic())
 Sum values in overlapped cells. The destination is limited to valid + ngrow cells. More...
 
void SumBoundary (int scomp, int ncomp, IntVect const &src_nghost, IntVect const &dst_nghost, const Periodicity &period=Periodicity::NonPeriodic())
 Sum values in overlapped cells. For computing the overlap, the dst is grown by dst_ngrow, while the src uses src_ngrow. More...
 
void SumBoundary_nowait (const Periodicity &period=Periodicity::NonPeriodic())
 
void SumBoundary_nowait (int scomp, int ncomp, const Periodicity &period=Periodicity::NonPeriodic())
 
void SumBoundary_nowait (int scomp, int ncomp, IntVect const &nghost, const Periodicity &period=Periodicity::NonPeriodic())
 
void SumBoundary_nowait (int scomp, int ncomp, IntVect const &src_nghost, IntVect const &dst_nghost, const Periodicity &period=Periodicity::NonPeriodic())
 
void SumBoundary_finish ()
 
void EnforcePeriodicity (const Periodicity &period)
 Fill ghost cells with values from their corresponding cells across periodic boundaries, regardless of whether the corresponding cells are valid. This differs from FillBoundary, which only fills from valid cells, and does not fill from ghost cells. The BoxArray is allowed to be overlapping. More...
 
void EnforcePeriodicity (int scomp, int ncomp, const Periodicity &period)
 
void EnforcePeriodicity (int scomp, int ncomp, const IntVect &nghost, const Periodicity &period)
 
void BuildMask (const Box &phys_domain, const Periodicity &period, value_type covered, value_type notcovered, value_type physbnd, value_type interior)
 
void FBEP_nowait (int scomp, int ncomp, const IntVect &nghost, const Periodicity &period, bool cross, bool enforce_periodicity_only=false, bool override_sync=false)
 
void FBEP_nowait (int scomp, int ncomp, const IntVect &nghost, const Periodicity &period, bool cross, bool enforce_periodicity_only, bool override_sync)
 
void FB_local_copy_cpu (const FB &TheFB, int scomp, int ncomp)
 
void PC_local_cpu (const CPC &thecpc, FabArray< FArrayBox > const &src, int scomp, int dcomp, int ncomp, CpOp op)
 
LayoutData< int > RecvLayoutMask (const CommMetaData &thecmd)
 
void FB_local_copy_gpu (const FB &TheFB, int scomp, int ncomp)
 
void PC_local_gpu (const CPC &thecpc, FabArray< FArrayBox > const &src, int scomp, int dcomp, int ncomp, CpOp op)
 
void CMD_local_setVal_gpu (value_type x, const CommMetaData &thecmd, int scomp, int ncomp)
 
void CMD_remote_setVal_gpu (value_type x, const CommMetaData &thecmd, int scomp, int ncomp)
 
void pack_send_buffer_gpu (FabArray< FArrayBox > const &src, int scomp, int ncomp, Vector< char * > const &send_data, Vector< std::size_t > const &send_size, Vector< CopyComTagsContainer const * > const &send_cctc)
 
void unpack_recv_buffer_gpu (FabArray< FArrayBox > &dst, int dcomp, int ncomp, Vector< char * > const &recv_data, Vector< std::size_t > const &recv_size, Vector< CopyComTagsContainer const * > const &recv_cctc, CpOp op, bool is_thread_safe)
 
void pack_send_buffer_cpu (FabArray< FArrayBox > const &src, int scomp, int ncomp, Vector< char * > const &send_data, Vector< std::size_t > const &send_size, Vector< CopyComTagsContainer const * > const &send_cctc)
 
void unpack_recv_buffer_cpu (FabArray< FArrayBox > &dst, int dcomp, int ncomp, Vector< char * > const &recv_data, Vector< std::size_t > const &recv_size, Vector< CopyComTagsContainer const * > const &recv_cctc, CpOp op, bool is_thread_safe)
 
F::value_type norminf (int comp, int ncomp, IntVect const &nghost, bool local=false,[[maybe_unused]] bool ignore_covered=false) const
 Return infinity norm. More...
 
F::value_type norminf (FabArray< IFAB > const &mask, int comp, int ncomp, IntVect const &nghost, bool local=false) const
 Return infinity norm in masked region. More...
 
- Public Member Functions inherited from amrex::FabArrayBase
 FabArrayBase ()=default
 
 ~FabArrayBase ()=default
 
 FabArrayBase (const BoxArray &bxs, const DistributionMapping &dm, int nvar, int ngrow)
 
 FabArrayBase (const BoxArray &bxs, const DistributionMapping &dm, int nvar, const IntVect &ngrow)
 
 FabArrayBase (FabArrayBase &&rhs) noexcept=default
 
 FabArrayBase (const FabArrayBase &rhs)=default
 
FabArrayBaseoperator= (const FabArrayBase &rhs)=default
 
FabArrayBaseoperator= (FabArrayBase &&rhs)=default
 
void define (const BoxArray &bxs, const DistributionMapping &dm, int nvar, int ngrow)
 
void define (const BoxArray &bxs, const DistributionMapping &dm, int nvar, const IntVect &ngrow)
 
int nGrow (int direction=0) const noexcept
 Return the grow factor that defines the region of definition. More...
 
IntVect nGrowVect () const noexcept
 
int nComp () const noexcept
 Return number of variables (aka components) associated with each point. More...
 
IndexType ixType () const noexcept
 Return index type. More...
 
bool empty () const noexcept
 
const BoxArrayboxArray () const noexcept
 Return a constant reference to the BoxArray that defines the valid region associated with this FabArray. More...
 
Box box (int K) const noexcept
 Return the Kth Box in the BoxArray. That is, the valid region of the Kth grid. More...
 
Box fabbox (int K) const noexcept
 Return the Kth FABs Box in the FabArray. That is, the region the Kth fab is actually defined on. More...
 
int size () const noexcept
 Return the number of FABs in the FabArray. More...
 
int local_size () const noexcept
 Return the number of local FABs in the FabArray. More...
 
const Vector< int > & IndexArray () const noexcept
 Return constant reference to indices in the FabArray that we have access. More...
 
int localindex (int K) const noexcept
 Return local index in the vector of FABs. More...
 
const DistributionMappingDistributionMap () const noexcept
 Return constant reference to associated DistributionMapping. More...
 
bool is_nodal () const noexcept
 This tests on whether the FabArray is fully nodal. More...
 
bool is_nodal (int dir) const noexcept
 This tests on whether the FabArray is nodal in direction dir. More...
 
bool is_cell_centered () const noexcept
 This tests on whether the FabArray is cell-centered. More...
 
void setMultiGhost (bool a_multi_ghost)
 
IntVect nGrowFilled () const noexcept
 
void setNGrowFilled (IntVect const &ng) noexcept
 
bool isFusingCandidate () const noexcept
 Is this a good candidate for kernel fusing? More...
 
BDKey getBDKey () const noexcept
 
void updateBDKey ()
 
void flushFPinfo (bool no_assertion=false) const
 
void flushCFinfo (bool no_assertion=false) const
 
const TileArraygetTileArray (const IntVect &tilesize) const
 
void clear ()
 
const std::vector< bool > & OwnerShip () const noexcept
 Return owenership of fabs. The concept of ownership only applies when UPC++ team is used. In that case, each fab is shared by team workers, with one taking the ownership. More...
 
bool isOwner (int li) const noexcept
 
void buildTileArray (const IntVect &tilesize, TileArray &ta) const
 
void flushTileArray (const IntVect &tilesize=IntVect::TheZeroVector(), bool no_assertion=false) const
 
void define_fb_metadata (CommMetaData &cmd, const IntVect &nghost, bool cross, const Periodicity &period, bool multi_ghost) const
 
const FBgetFB (const IntVect &nghost, const Periodicity &period, bool cross=false, bool enforce_periodicity_only=false, bool override_sync=false) const
 
void flushFB (bool no_assertion=false) const
 This flushes its own FB. More...
 
const CPCgetCPC (const IntVect &dstng, const FabArrayBase &src, const IntVect &srcng, const Periodicity &period, bool to_ghost_cells_only=false) const
 
void flushCPC (bool no_assertion=false) const
 This flushes its own CPC. More...
 
const RB90getRB90 (const IntVect &nghost, const Box &domain) const
 
void flushRB90 (bool no_assertion=false) const
 This flushes its own RB90. More...
 
const RB180getRB180 (const IntVect &nghost, const Box &domain) const
 
void flushRB180 (bool no_assertion=false) const
 This flushes its own RB180. More...
 
const PolarBgetPolarB (const IntVect &nghost, const Box &domain) const
 
void flushPolarB (bool no_assertion=false) const
 This flushes its own PolarB. More...
 
ParForInfo const & getParForInfo (const IntVect &nghost, int nthreads) const
 
void flushParForInfo (bool no_assertion=false) const
 
void clearThisBD (bool no_assertion=false) const
 clear BD count and caches associated with this BD, if no other is using this BD. More...
 
void addThisBD ()
 add the current BD into BD count database More...
 

Static Public Member Functions

static Real Dot (const MultiFab &x, int xcomp, const MultiFab &y, int ycomp, int numcomp, int nghost, bool local=false)
 Returns the dot product of two MultiFabs. More...
 
static Real Dot (const MultiFab &x, int xcomp, int numcomp, int nghost, bool local=false)
 Returns the dot product of a MultiFab with itself. More...
 
static Real Dot (const iMultiFab &mask, const MultiFab &x, int xcomp, const MultiFab &y, int ycomp, int numcomp, int nghost, bool local=false)
 
static void Add (MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
 Add src to dst including nghost ghost cells. The two MultiFabs MUST have the same underlying BoxArray. More...
 
static void Add (MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, const IntVect &nghost)
 
static void Copy (MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
 Copy from src to dst including nghost ghost cells. The two MultiFabs MUST have the same underlying BoxArray. The copy is local. The parallel copy function is in FabArray.H. More...
 
static void Copy (MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, const IntVect &nghost)
 
static void Swap (MultiFab &dst, MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
 Swap from src to dst including nghost ghost cells. The two MultiFabs MUST have the same underlying BoxArray. The swap is local. More...
 
static void Swap (MultiFab &dst, MultiFab &src, int srccomp, int dstcomp, int numcomp, const IntVect &nghost)
 
static void Subtract (MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
 Subtract src from dst including nghost ghost cells. The two MultiFabs MUST have the same underlying BoxArray. More...
 
static void Subtract (MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, const IntVect &nghost)
 
static void Multiply (MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
 Multiply dst by src including nghost ghost cells. The two MultiFabs MUST have the same underlying BoxArray. More...
 
static void Multiply (MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, const IntVect &nghost)
 
static void Divide (MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
 Divide dst by src including nghost ghost cells. The two MultiFabs MUST have the same underlying BoxArray. More...
 
static void Divide (MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, const IntVect &nghost)
 
static void Saxpy (MultiFab &dst, Real a, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
 dst += a*src More...
 
static void Xpay (MultiFab &dst, Real a, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
 dst = src + a*dst More...
 
static void LinComb (MultiFab &dst, Real a, const MultiFab &x, int xcomp, Real b, const MultiFab &y, int ycomp, int dstcomp, int numcomp, int nghost)
 dst = a*x + b*y More...
 
static void AddProduct (MultiFab &dst, const MultiFab &src1, int comp1, const MultiFab &src2, int comp2, int dstcomp, int numcomp, int nghost)
 dst += src1*src2 More...
 
static void AddProduct (MultiFab &dst, const MultiFab &src1, int comp1, const MultiFab &src2, int comp2, int dstcomp, int numcomp, const IntVect &nghost)
 
static void Initialize ()
 
static void Finalize ()
 
- Static Public Member Functions inherited from amrex::FabArray< FArrayBox >
static void pack_send_buffer_gpu (FabArray< FArrayBox > const &src, int scomp, int ncomp, Vector< char * > const &send_data, Vector< std::size_t > const &send_size, Vector< const CopyComTagsContainer * > const &send_cctc)
 
static void unpack_recv_buffer_gpu (FabArray< FArrayBox > &dst, int dcomp, int ncomp, Vector< char * > const &recv_data, Vector< std::size_t > const &recv_size, Vector< const CopyComTagsContainer * > const &recv_cctc, CpOp op, bool is_thread_safe)
 
static void pack_send_buffer_cpu (FabArray< FArrayBox > const &src, int scomp, int ncomp, Vector< char * > const &send_data, Vector< std::size_t > const &send_size, Vector< const CopyComTagsContainer * > const &send_cctc)
 
static void unpack_recv_buffer_cpu (FabArray< FArrayBox > &dst, int dcomp, int ncomp, Vector< char * > const &recv_data, Vector< std::size_t > const &recv_size, Vector< const CopyComTagsContainer * > const &recv_cctc, CpOp op, bool is_thread_safe)
 
static void Saxpy (FabArray< FArrayBox > &y, value_type a, FabArray< FArrayBox > const &x, int xcomp, int ycomp, int ncomp, IntVect const &nghost)
 y += a*x More...
 
static void Xpay (FabArray< FArrayBox > &y, value_type a, FabArray< FArrayBox > const &x, int xcomp, int ycomp, int ncomp, IntVect const &nghost)
 y = x + a*y More...
 
static void LinComb (FabArray< FArrayBox > &dst, value_type a, const FabArray< FArrayBox > &x, int xcomp, value_type b, const FabArray< FArrayBox > &y, int ycomp, int dstcomp, int numcomp, const IntVect &nghost)
 dst = a*x + b*y More...
 
- Static Public Member Functions inherited from amrex::FabArrayBase
static Long bytesOfMapOfCopyComTagContainers (const MapOfCopyComTagContainers &)
 
static void Initialize ()
 Initialize from ParmParse with "fabarray" prefix. More...
 
static void Finalize ()
 
static const FPinfoTheFPinfo (const FabArrayBase &srcfa, const FabArrayBase &dstfa, const IntVect &dstng, const BoxConverter &coarsener, const Geometry &fgeom, const Geometry &cgeom, const EB2::IndexSpace *)
 
static const CFinfoTheCFinfo (const FabArrayBase &finefa, const Geometry &finegm, const IntVect &ng, bool include_periodic, bool include_physbndry)
 
static void updateMemUsage (std::string const &tag, Long nbytes, Arena const *ar)
 
static void printMemUsage ()
 
static Long queryMemUsage (const std::string &tag=std::string("All"))
 
static Long queryMemUsageHWM (const std::string &tag=std::string("All"))
 
static void pushRegionTag (const char *t)
 
static void pushRegionTag (std::string t)
 
static void popRegionTag ()
 
static void flushTileArrayCache ()
 This flushes the entire cache. More...
 
static void flushFBCache ()
 This flushes the entire cache. More...
 
static void flushCPCache ()
 This flusheds the entire cache. More...
 
static void flushRB90Cache ()
 This flushes the entire cache. More...
 
static void flushRB180Cache ()
 This flushes the entire cache. More...
 
static void flushPolarBCache ()
 This flushes the entire cache. More...
 
static void flushParForCache ()
 
static bool getAllocSingleChunk ()
 

Private Types

using CopyComTagsContainer = FabArrayBase::CopyComTagsContainer
 Some useful typedefs. More...
 
using MapOfCopyComTagContainers = FabArrayBase::MapOfCopyComTagContainers
 

Private Member Functions

void initVal ()
 

Additional Inherited Members

- Public Types inherited from amrex::FabArray< FArrayBox >
using value_type = typename std::conditional_t< IsBaseFab< FArrayBox >::value, FArrayBox, FABType >::value_type
 
using fab_type = FArrayBox
 
- Public Types inherited from amrex::FabArrayBase
enum  CpOp { COPY = 0 , ADD = 1 }
 parallel copy or add More...
 
using CopyComTagsContainer = CopyComTag::CopyComTagsContainer
 
using MapOfCopyComTagContainers = CopyComTag::MapOfCopyComTagContainers
 
using FPinfoCache = std::multimap< BDKey, FabArrayBase::FPinfo * >
 
using FPinfoCacheIter = FPinfoCache::iterator
 
using CFinfoCache = std::multimap< BDKey, FabArrayBase::CFinfo * >
 
using CFinfoCacheIter = CFinfoCache::iterator
 
using TAMap = std::map< std::pair< IntVect, IntVect >, TileArray >
 
using TACache = std::map< BDKey, TAMap >
 
using FBCache = std::multimap< BDKey, FabArrayBase::FB * >
 
using FBCacheIter = FBCache::iterator
 
using CPCache = std::multimap< BDKey, FabArrayBase::CPC * >
 
using CPCacheIter = CPCache::iterator
 
using RB90Cache = std::multimap< BDKey, FabArrayBase::RB90 * >
 
using RB90CacheIter = RB90Cache::iterator
 
using RB180Cache = std::multimap< BDKey, FabArrayBase::RB180 * >
 
using RB180CacheIter = RB180Cache::iterator
 
using PolarBCache = std::multimap< BDKey, FabArrayBase::PolarB * >
 
using PolarBCacheIter = PolarBCache::iterator
 
- Public Attributes inherited from amrex::FabArray< FArrayBox >
std::unique_ptr< FBData< FArrayBox > > fbd
 
std::unique_ptr< PCData< FArrayBox > > pcd
 
std::unique_ptr< FabArray< FArrayBox > > os_temp
 
- Public Attributes inherited from amrex::FabArrayBase
BoxArray boxarray
 
DistributionMapping distributionMap
 
Vector< intindexArray
 
std::vector< bool > ownership
 
IntVect n_grow
 
int n_comp
 
BDKey m_bdkey
 
IntVect n_filled
 
bool m_multi_ghost = false
 
- Static Public Attributes inherited from amrex::FabArrayBase
static AMREX_EXPORT IntVect mfiter_tile_size
 Default tilesize in MFIter. More...
 
static AMREX_EXPORT int MaxComp
 The maximum number of components to copy() at a time. More...
 
static AMREX_EXPORT IntVect comm_tile_size
 communication tile size More...
 
static FPinfoCache m_TheFillPatchCache
 
static CacheStats m_FPinfo_stats
 
static CFinfoCache m_TheCrseFineCache
 
static CacheStats m_CFinfo_stats
 
static std::map< std::string, meminfom_mem_usage
 
static AMREX_EXPORT std::vector< std::string > m_region_tag
 
static TACache m_TheTileArrayCache
 
static CacheStats m_TAC_stats
 
static FBCache m_TheFBCache
 
static CacheStats m_FBC_stats
 
static CPCache m_TheCPCache
 
static CacheStats m_CPC_stats
 
static RB90Cache m_TheRB90Cache
 
static RB180Cache m_TheRB180Cache
 
static PolarBCache m_ThePolarBCache
 
static std::multimap< BDKey, ParForInfo * > m_TheParForCache
 
static std::map< BDKey, intm_BD_count
 Keep track of how many FabArrays are built with the same BDKey. More...
 
static AMREX_EXPORT FabArrayStats m_FA_stats
 
static AMREX_EXPORT bool m_alloc_single_chunk
 
- Protected Member Functions inherited from amrex::FabArray< FArrayBox >
bool SharedMemory () const noexcept
 
- Protected Attributes inherited from amrex::FabArray< FArrayBox >
std::unique_ptr< FabFactory< FArrayBox > > m_factory
 
DataAllocator m_dallocator
 
std::unique_ptr< detail::SingleChunkArenam_single_chunk_arena
 
Long m_single_chunk_size
 
bool define_function_called
 has define() been called? More...
 
std::vector< FArrayBox * > m_fabs_v
 The data. More...
 
void * m_dp_arrays
 
void * m_hp_arrays
 
MultiArray4< value_typem_arrays
 
MultiArray4< value_type const > m_const_arrays
 
Vector< std::string > m_tags
 
ShMem shmem
 

Detailed Description

A collection (stored as an array) of FArrayBox objects.

This class is useful for storing floating point data on a domain defined by a union of rectangular regions embedded in a uniform index space. MultiFab class extends the function of the underlying FabArray class just as the FArrayBox class extends the function of BaseFab<Real>. Additional member functions are defined for I/O and simple arithmetic operations on these aggregate objects.

This class does NOT provide a copy constructor or assignment operator.

Member Typedef Documentation

◆ CopyComTagsContainer

◆ MapOfCopyComTagContainers

Constructor & Destructor Documentation

◆ MultiFab() [1/7]

amrex::MultiFab::MultiFab ( )
noexcept

Constructs an empty MultiFab.

Data can be defined at a later time using the define member functions inherited from FabArray.

◆ MultiFab() [2/7]

amrex::MultiFab::MultiFab ( Arena a)
explicitnoexcept

Constructs an empty MultiFab.

Data can be defined at a later time using the define member functions. If define is called later with a nullptr as MFInfo's arena, the default Arena a will be used. If the arena in MFInfo is not a nullptr, the MFInfo's arena will be used.

◆ MultiFab() [3/7]

amrex::MultiFab::MultiFab ( const BoxArray bxs,
const DistributionMapping dm,
int  ncomp,
int  ngrow,
const MFInfo info = MFInfo(),
const FabFactory< FArrayBox > &  factory = FArrayBoxFactory() 
)

Constructs a MultiFab.

The size of the FArrayBox is given by the Box grown by ngrow, and the number of components is given by ncomp. If info is set to not allocating memory, then no FArrayBoxes are allocated at this time but can be defined later.

Parameters
bxsa valid region
dma DistribuionMapping
ncompnumber of components
ngrownumber of cells the region grows
infoMFInfo

◆ MultiFab() [4/7]

amrex::MultiFab::MultiFab ( const BoxArray bxs,
const DistributionMapping dm,
int  ncomp,
const IntVect ngrow,
const MFInfo info = MFInfo(),
const FabFactory< FArrayBox > &  factory = FArrayBoxFactory() 
)

◆ MultiFab() [5/7]

amrex::MultiFab::MultiFab ( const MultiFab rhs,
MakeType  maketype,
int  scomp,
int  ncomp 
)

Make an alias MultiFab.

Note that maketype must be amrex::make_alias, scomp is the starting component of the alias, and ncomp is the number of components in the new aliasing MultiFab.

◆ ~MultiFab()

amrex::MultiFab::~MultiFab ( )

◆ MultiFab() [6/7]

amrex::MultiFab::MultiFab ( MultiFab &&  rhs)
noexcept

◆ MultiFab() [7/7]

amrex::MultiFab::MultiFab ( const MultiFab rhs)
delete

Member Function Documentation

◆ Add() [1/2]

void amrex::MultiFab::Add ( MultiFab dst,
const MultiFab src,
int  srccomp,
int  dstcomp,
int  numcomp,
const IntVect nghost 
)
static

◆ Add() [2/2]

void amrex::MultiFab::Add ( MultiFab dst,
const MultiFab src,
int  srccomp,
int  dstcomp,
int  numcomp,
int  nghost 
)
static

Add src to dst including nghost ghost cells. The two MultiFabs MUST have the same underlying BoxArray.

◆ AddProduct() [1/2]

void amrex::MultiFab::AddProduct ( MultiFab dst,
const MultiFab src1,
int  comp1,
const MultiFab src2,
int  comp2,
int  dstcomp,
int  numcomp,
const IntVect nghost 
)
static

◆ AddProduct() [2/2]

void amrex::MultiFab::AddProduct ( MultiFab dst,
const MultiFab src1,
int  comp1,
const MultiFab src2,
int  comp2,
int  dstcomp,
int  numcomp,
int  nghost 
)
static

dst += src1*src2

◆ AverageSync()

void amrex::MultiFab::AverageSync ( const Periodicity period = Periodicity::NonPeriodic())

Sync up nodal data via averaging.

◆ contains_inf() [1/3]

bool amrex::MultiFab::contains_inf ( bool  local = false) const

Are there any Infs in the MF? This may return false, even if the MF contains Infs, if the machine doesn't support the appropriate Inf testing functions. This version checks all components and grow cells.

◆ contains_inf() [2/3]

bool amrex::MultiFab::contains_inf ( int  scomp,
int  ncomp,
const IntVect ngrow,
bool  local = false 
) const

◆ contains_inf() [3/3]

bool amrex::MultiFab::contains_inf ( int  scomp,
int  ncomp,
int  ngrow = 0,
bool  local = false 
) const

◆ contains_nan() [1/3]

bool amrex::MultiFab::contains_nan ( bool  local = false) const

Are there any NaNs in the MF? This may return false, even if the MF contains NaNs, if the machine doesn't support the appropriate NaN testing functions.

This version checks all components and grow cells.

◆ contains_nan() [2/3]

bool amrex::MultiFab::contains_nan ( int  scomp,
int  ncomp,
const IntVect ngrow,
bool  local = false 
) const

◆ contains_nan() [3/3]

bool amrex::MultiFab::contains_nan ( int  scomp,
int  ncomp,
int  ngrow = 0,
bool  local = false 
) const

◆ Copy() [1/2]

void amrex::MultiFab::Copy ( MultiFab dst,
const MultiFab src,
int  srccomp,
int  dstcomp,
int  numcomp,
const IntVect nghost 
)
static

◆ Copy() [2/2]

void amrex::MultiFab::Copy ( MultiFab dst,
const MultiFab src,
int  srccomp,
int  dstcomp,
int  numcomp,
int  nghost 
)
static

Copy from src to dst including nghost ghost cells. The two MultiFabs MUST have the same underlying BoxArray. The copy is local. The parallel copy function is in FabArray.H.

◆ deepCopy()

MultiFab amrex::MultiFab::deepCopy ( ) const

Create a deep copy of this MultiFab and its data.

This uses the same arena and factory in the copy and copies all data.

◆ define() [1/2]

void amrex::MultiFab::define ( const BoxArray bxs,
const DistributionMapping dm,
int  nvar,
const IntVect ngrow,
const MFInfo info = MFInfo(),
const FabFactory< FArrayBox > &  factory = FArrayBoxFactory() 
)

◆ define() [2/2]

void amrex::MultiFab::define ( const BoxArray bxs,
const DistributionMapping dm,
int  nvar,
int  ngrow,
const MFInfo info = MFInfo(),
const FabFactory< FArrayBox > &  factory = FArrayBoxFactory() 
)

◆ divide()

void amrex::MultiFab::divide ( const MultiFab mf,
int  strt_comp,
int  num_comp,
int  nghost 
)

This function divides the values of the cells in mf from the corresponding cells of this MultiFab. mf is required to have the same BoxArray or "valid region" as this MultiFab. The division is done only to num_comp components, starting with component number strt_comp. The parameter nghost specifies the number of boundary cells that will be modified. If nghost == 0, only the valid region of each FArrayBox will be modified. Note, nothing is done to protect against divide by zero.

◆ Divide() [1/2]

void amrex::MultiFab::Divide ( MultiFab dst,
const MultiFab src,
int  srccomp,
int  dstcomp,
int  numcomp,
const IntVect nghost 
)
static

◆ Divide() [2/2]

void amrex::MultiFab::Divide ( MultiFab dst,
const MultiFab src,
int  srccomp,
int  dstcomp,
int  numcomp,
int  nghost 
)
static

Divide dst by src including nghost ghost cells. The two MultiFabs MUST have the same underlying BoxArray.

◆ Dot() [1/3]

Real amrex::MultiFab::Dot ( const iMultiFab mask,
const MultiFab x,
int  xcomp,
const MultiFab y,
int  ycomp,
int  numcomp,
int  nghost,
bool  local = false 
)
static

◆ Dot() [2/3]

Real amrex::MultiFab::Dot ( const MultiFab x,
int  xcomp,
const MultiFab y,
int  ycomp,
int  numcomp,
int  nghost,
bool  local = false 
)
static

Returns the dot product of two MultiFabs.

◆ Dot() [3/3]

Real amrex::MultiFab::Dot ( const MultiFab x,
int  xcomp,
int  numcomp,
int  nghost,
bool  local = false 
)
static

Returns the dot product of a MultiFab with itself.

◆ Finalize()

void amrex::MultiFab::Finalize ( )
static

◆ Initialize()

void amrex::MultiFab::Initialize ( )
static

◆ initVal()

void amrex::MultiFab::initVal ( )
private

◆ invert() [1/4]

void amrex::MultiFab::invert ( Real  numerator,
const Box region,
int  comp,
int  num_comp,
int  nghost = 0 
)

Identical to the previous version of invert(), with the restriction that the subregion is further constrained to the intersection with Box region. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified.

◆ invert() [2/4]

void amrex::MultiFab::invert ( Real  numerator,
const Box region,
int  nghost 
)

Replaces the value of each cell in the specified subregion of the MultiFab, that also intersects the Box region, with its reciprocal multiplied by the value of numerator. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified.

◆ invert() [3/4]

void amrex::MultiFab::invert ( Real  numerator,
int  comp,
int  num_comp,
int  nghost = 0 
)

Replaces the value of each cell in the specified subregion of the MultiFab with its reciprocal multiplied by the value of numerator. The subregion consists of the num_comp components starting at component comp. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified.

◆ invert() [4/4]

void amrex::MultiFab::invert ( Real  numerator,
int  nghost 
)

Replaces the value of each cell in the specified subregion of the MultiFab with its reciprocal multiplied by the value of numerator. The value of nghost specifies the number of cells in the boundary region that should be modified.

◆ is_finite() [1/3]

bool amrex::MultiFab::is_finite ( bool  local = false) const

Are the numbers in the MF finite (i.e., neither nan nor inf)? This may return true, even if the MF contains NaNs or Infs, if the machine doesn't support the appropriate NaN and Inf testing functions.

This version checks all components and grow cells.

◆ is_finite() [2/3]

bool amrex::MultiFab::is_finite ( int  scomp,
int  ncomp,
const IntVect ngrow,
bool  local = false 
) const

◆ is_finite() [3/3]

bool amrex::MultiFab::is_finite ( int  scomp,
int  ncomp,
int  ngrow = 0,
bool  local = false 
) const

◆ LinComb()

void amrex::MultiFab::LinComb ( MultiFab dst,
Real  a,
const MultiFab x,
int  xcomp,
Real  b,
const MultiFab y,
int  ycomp,
int  dstcomp,
int  numcomp,
int  nghost 
)
static

dst = a*x + b*y

◆ max() [1/2]

Real amrex::MultiFab::max ( const Box region,
int  comp,
int  nghost = 0,
bool  local = false 
) const

Identical to the previous max() function, but confines its search to intersection of Box b and the MultiFab.

◆ max() [2/2]

Real amrex::MultiFab::max ( int  comp,
int  nghost = 0,
bool  local = false 
) const

Returns the maximum value contained in component comp of the MultiFab.

The parameter nghost determines the number of boundary cells to search for the maximum. The default is to search only the valid regions of the FArrayBoxes.

◆ maxIndex()

IntVect amrex::MultiFab::maxIndex ( int  comp,
int  nghost = 0 
) const

◆ min() [1/2]

Real amrex::MultiFab::min ( const Box region,
int  comp,
int  nghost = 0,
bool  local = false 
) const

Identical to min() function, but confines its search to intersection of Box b and the MultiFab.

◆ min() [2/2]

Real amrex::MultiFab::min ( int  comp,
int  nghost = 0,
bool  local = false 
) const

Returns the minimum value contained in component comp of the MultiFab.

The parameter nghost determines the number of boundary cells to search for the minimum. The default is to search only the valid regions of the FArrayBoxes.

◆ minIndex()

IntVect amrex::MultiFab::minIndex ( int  comp,
int  nghost = 0 
) const

◆ minus()

void amrex::MultiFab::minus ( const MultiFab mf,
int  strt_comp,
int  num_comp,
int  nghost 
)

This function subtracts the values of the cells in mf from the corresponding cells of this MultiFab. mf is required to have the same BoxArray or "valid region" as this MultiFab. The subtraction is done only to num_comp components, starting with component number strt_comp. The parameter nghost specifies the number of boundary cells that will be modified. If nghost == 0, only the valid region of each FArrayBox will be modified.

◆ mult() [1/4]

void amrex::MultiFab::mult ( Real  val,
const Box region,
int  comp,
int  num_comp,
int  nghost = 0 
)

Identical to the previous version of mult(), with the restriction that the subregion is further constrained to the intersection with Box region. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified.

◆ mult() [2/4]

void amrex::MultiFab::mult ( Real  val,
const Box region,
int  nghost = 0 
)

Scales the value of each cell in the valid region of each component of the MultiFab by the scalar val (a[i] <- a[i]*val), that also intersects the Box region. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified.

◆ mult() [3/4]

void amrex::MultiFab::mult ( Real  val,
int  comp,
int  num_comp,
int  nghost = 0 
)

Scales the value of each cell in the specified subregion of the MultiFab by the scalar val (a[i] <- a[i]*val). The subregion consists of the num_comp components starting at component comp. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified.

◆ mult() [4/4]

void amrex::MultiFab::mult ( Real  val,
int  nghost = 0 
)

Scales the value of each cell in the valid region of each component of the MultiFab by the scalar val (a[i] <- a[i]*val). The value of nghost specifies the number of cells in the boundary region that should be modified.

◆ Multiply() [1/2]

void amrex::MultiFab::Multiply ( MultiFab dst,
const MultiFab src,
int  srccomp,
int  dstcomp,
int  numcomp,
const IntVect nghost 
)
static

◆ Multiply() [2/2]

void amrex::MultiFab::Multiply ( MultiFab dst,
const MultiFab src,
int  srccomp,
int  dstcomp,
int  numcomp,
int  nghost 
)
static

Multiply dst by src including nghost ghost cells. The two MultiFabs MUST have the same underlying BoxArray.

◆ negate() [1/4]

void amrex::MultiFab::negate ( const Box region,
int  comp,
int  num_comp,
int  nghost = 0 
)

Identical to the previous version of negate(), with the restriction that the subregion is further constrained to the intersection with Box region.

◆ negate() [2/4]

void amrex::MultiFab::negate ( const Box region,
int  nghost = 0 
)

Negates the value of each cell in the valid region of the MultiFab that also intersects the Box region. The value of nghost specifies the number of cells in the boundary region that should be modified.

◆ negate() [3/4]

void amrex::MultiFab::negate ( int  comp,
int  num_comp,
int  nghost = 0 
)

Negates the value of each cell in the specified subregion of the MultiFab. The subregion consists of the num_comp components starting at component comp. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified.

◆ negate() [4/4]

void amrex::MultiFab::negate ( int  nghost = 0)

Negates the value of each cell in the valid region of the MultiFab. The value of nghost specifies the number of cells in the boundary region that should be modified.

◆ norm0() [1/4]

Real amrex::MultiFab::norm0 ( const iMultiFab mask,
int  comp = 0,
int  nghost = 0,
bool  local = false 
) const

◆ norm0() [2/4]

Vector< Real > amrex::MultiFab::norm0 ( const Vector< int > &  comps,
int  nghost = 0,
bool  local = false,
bool  ignore_covered = false 
) const

Returns the maximum absolute values contained in each component of comps of the MultiFab. nghost ghost cells are used.

◆ norm0() [3/4]

Real amrex::MultiFab::norm0 ( int  comp,
int  ncomp,
IntVect const &  nghost,
bool  local = false,
bool  ignore_covered = false 
) const

◆ norm0() [4/4]

Real amrex::MultiFab::norm0 ( int  comp = 0,
int  nghost = 0,
bool  local = false,
bool  ignore_covered = false 
) const

Returns the maximum absolute value contained in component comp of the MultiFab.

◆ norm1() [1/3]

Vector< Real > amrex::MultiFab::norm1 ( const Vector< int > &  comps,
int  ngrow = 0,
bool  local = false 
) const

Returns the L1 norm of each component of "comps" over the MultiFab. ngrow ghost cells are used.

◆ norm1() [2/3]

Real amrex::MultiFab::norm1 ( int  comp,
const Periodicity period,
bool  ignore_covered = false 
) const

Returns the L1 norm of component comp over the MultiFab.

No ghost cells are used. This version has no double counting for nodal data.

◆ norm1() [3/3]

Real amrex::MultiFab::norm1 ( int  comp = 0,
int  ngrow = 0,
bool  local = false 
) const

Returns the L1 norm of component comp over the MultiFab. ngrow ghost cells are used.

◆ norm2() [1/3]

Vector< Real > amrex::MultiFab::norm2 ( const Vector< int > &  comps) const

Returns the L2 norm of each component of "comps" over the MultiFab. No ghost cells are used.

◆ norm2() [2/3]

Real amrex::MultiFab::norm2 ( int  comp,
const Periodicity period 
) const

Returns the L2 norm of component comp over the MultiFab. No ghost cells are used. This version has no double counting for nodal data.

◆ norm2() [3/3]

Real amrex::MultiFab::norm2 ( int  comp = 0) const

Returns the L2 norm of component comp over the MultiFab. No ghost cells are used.

◆ norminf() [1/3]

Real amrex::MultiFab::norminf ( const iMultiFab mask,
int  comp = 0,
int  nghost = 0,
bool  local = false 
) const
inline

◆ norminf() [2/3]

Vector<Real> amrex::MultiFab::norminf ( const Vector< int > &  comps,
int  nghost = 0,
bool  local = false,
bool  ignore_covered = false 
) const
inline

◆ norminf() [3/3]

Real amrex::MultiFab::norminf ( int  comp = 0,
int  nghost = 0,
bool  local = false,
bool  ignore_covered = false 
) const
inline

◆ operator=() [1/3]

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

◆ operator=() [2/3]

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

◆ operator=() [3/3]

MultiFab & amrex::MultiFab::operator= ( Real  r)

◆ OverlapMask()

std::unique_ptr< MultiFab > amrex::MultiFab::OverlapMask ( const Periodicity period = Periodicity::NonPeriodic()) const

Return a mask indicating how many duplicates are in each point.

◆ OverrideSync()

void amrex::MultiFab::OverrideSync ( const iMultiFab msk,
const Periodicity period = Periodicity::NonPeriodic() 
)

Sync up nodal data with owners overriding non-owners.

◆ OwnerMask()

std::unique_ptr< iMultiFab > amrex::MultiFab::OwnerMask ( const Periodicity period = Periodicity::NonPeriodic()) const

Owner is the grid with the lowest grid number containing the data.

◆ plus() [1/5]

void amrex::MultiFab::plus ( const MultiFab mf,
int  strt_comp,
int  num_comp,
int  nghost 
)

This function adds the values of the cells in mf to the corresponding cells of this MultiFab. mf is required to have the same BoxArray or "valid region" as this MultiFab. The addition is done only to num_comp components, starting with component number strt_comp. The parameter nghost specifies the number of boundary cells that will be modified. If nghost == 0, only the valid region of each FArrayBox will be modified.

◆ plus() [2/5]

void amrex::MultiFab::plus ( Real  val,
const Box region,
int  comp,
int  num_comp,
int  nghost = 0 
)

Identical to the previous version of plus(), with the restriction that the subregion is further constrained to the intersection with Box region.

◆ plus() [3/5]

void amrex::MultiFab::plus ( Real  val,
const Box region,
int  nghost 
)

Adds the scalar value val to the value of each cell in the valid region of each component of the MultiFab, that also intersects the Box region. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified.

◆ plus() [4/5]

void amrex::MultiFab::plus ( Real  val,
int  comp,
int  num_comp,
int  nghost = 0 
)

Adds the scalar value val to the value of each cell in the specified subregion of the MultiFab.

The subregion consists of the num_comp components starting at component comp. The value of nghost specifies the number of cells in the boundary region of each FArrayBox in the subregion that should be modified.

◆ plus() [5/5]

void amrex::MultiFab::plus ( Real  val,
int  nghost 
)

Adds the scalar value val to the value of each cell in the valid region of each component of the MultiFab. The value of nghost specifies the number of cells in the boundary region that should be modified.

◆ Saxpy()

void amrex::MultiFab::Saxpy ( MultiFab dst,
Real  a,
const MultiFab src,
int  srccomp,
int  dstcomp,
int  numcomp,
int  nghost 
)
static

dst += a*src

◆ Subtract() [1/2]

void amrex::MultiFab::Subtract ( MultiFab dst,
const MultiFab src,
int  srccomp,
int  dstcomp,
int  numcomp,
const IntVect nghost 
)
static

◆ Subtract() [2/2]

void amrex::MultiFab::Subtract ( MultiFab dst,
const MultiFab src,
int  srccomp,
int  dstcomp,
int  numcomp,
int  nghost 
)
static

Subtract src from dst including nghost ghost cells. The two MultiFabs MUST have the same underlying BoxArray.

◆ sum() [1/2]

Real amrex::MultiFab::sum ( Box const &  region,
int  comp = 0,
bool  local = false 
) const

Returns the sum of component "comp" in the given "region". – no ghost cells are included.

◆ sum() [2/2]

Real amrex::MultiFab::sum ( int  comp = 0,
bool  local = false 
) const

Returns the sum of component "comp" over the MultiFab – no ghost cells are included.

◆ sum_unique() [1/2]

Real amrex::MultiFab::sum_unique ( Box const &  region,
int  comp = 0,
bool  local = false 
) const

Returns the unique sum of component "comp" in the given region. Non-unique points owned by multiple boxes in the MultiFab are only added once. No ghost cells are included. This function does not take periodicity into account in the determination of uniqueness of points.

◆ sum_unique() [2/2]

Real amrex::MultiFab::sum_unique ( int  comp = 0,
bool  local = false,
const Periodicity period = Periodicity::NonPeriodic() 
) const

Same as sum with local =false, but for non-cell-centered data, this only adds non-unique points that are owned by multiple boxes once.

◆ Swap() [1/2]

void amrex::MultiFab::Swap ( MultiFab dst,
MultiFab src,
int  srccomp,
int  dstcomp,
int  numcomp,
const IntVect nghost 
)
static

◆ Swap() [2/2]

void amrex::MultiFab::Swap ( MultiFab dst,
MultiFab src,
int  srccomp,
int  dstcomp,
int  numcomp,
int  nghost 
)
static

Swap from src to dst including nghost ghost cells. The two MultiFabs MUST have the same underlying BoxArray. The swap is local.

◆ WeightedSync()

void amrex::MultiFab::WeightedSync ( const MultiFab wgt,
const Periodicity period = Periodicity::NonPeriodic() 
)

Sync up nodal data with weights.

◆ Xpay()

void amrex::MultiFab::Xpay ( MultiFab dst,
Real  a,
const MultiFab src,
int  srccomp,
int  dstcomp,
int  numcomp,
int  nghost 
)
static

dst = src + a*dst


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