FabArray< FAB > Class Template Reference

#include <FabArray.H>

Inheritance diagram for FabArray< FAB >:

Inheritance graph
[legend]
Collaboration diagram for FabArray< FAB >:

Collaboration graph
[legend]

List of all members.

Public Types

typedef FAB::value_type value_type

Public Member Functions

 FabArray ()
 Constructs an empty FabArray<FAB>.
 FabArray (const BoxArray &bxs, int nvar, int ngrow, FabAlloc mem_mode=Fab_allocate)
 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.
 FabArray (const BoxArray &bxs, int nvar, int ngrow, const DistributionMapping &dm, FabAlloc mem_mode)
virtual ~FabArray ()
 The (virtual) destructor -- deletes all FABs in the array.
void define (const BoxArray &bxs, int nvar, int ngrow, FabAlloc mem_mode)
 Define this FabArray identically to that performed by the constructor having an analogous function signature.
void define (const BoxArray &bxs, int nvar, int ngrow, const DistributionMapping &dm, FabAlloc mem_mode)
bool ok () const
 Returns true if the FabArray is well-defined.
const FAB & operator[] (const MFIter &mfi) const
 Returns a constant reference to the FAB associated with the Kth element.
const FAB & get (const MFIter &mfi) const
FAB & operator[] (const MFIter &mfi)
 Returns a reference to the FAB associated with the Kth element.
FAB & get (const MFIter &mfi)
const FAB & operator[] (int K) const
 Returns a constant reference to the FAB associated with the Kth element.
FAB & operator[] (int K)
 Returns a reference to the FAB associated with the Kth element.
void setFab (int K, FAB *elem)
 Explicitly set the Kth FAB in the FabArray to point to elem.
void clear ()
 Releases FAB memory in the FabArray.
void setVal (value_type val)
 Set all components in the valid region of each FAB to val.
void operator= (const value_type &val)
void setVal (value_type val, int comp, int num_comp, int nghost=0)
 Set the value of num components in the valid region of each FAB in the FabArray, starting at component comp to val.
void setVal (value_type val, const Box &region, int comp, int num_comp, int nghost=0)
 Set the value of num 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.
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.
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.
void setBndry (value_type val)
 Set all values in the boundary region to val.
void setBndry (value_type val, int strt_comp, int ncomp)
 Set ncomp values in the boundary region, starting at start to val.
void copy (const FabArray< FAB > &fa)
 This function copies data from fa to this FabArray.
void copy (const FabArray< FAB > &src, int src_comp, int dest_comp, int num_comp)
 This function copies data from src to this FabArray.
void copy (FAB &dest) const
 Copies the values contained in the intersection of the valid region of this FabArray with the FAB dest into dest.
void copy (FAB &dest, const Box &subbox) const
 Copies the values contained in the intersection of the valid region of this FabArray with the FAB dest and the Box subbox into that subregion of dest.
void copy (FAB &dest, int src_comp, int dest_comp, int num_comp) const
 Copies the values contained in the intersection of the num component valid region of this FabArray, starting at component src, with the FAB dest into dest, starting at component dest in dest.
void copy (FAB &dest, const Box &subbox, int src_comp, int dest_comp, int num_comp) const
 Copies the values contained in the intersection of the num component valid region of this FabArray, starting at component src, with the FAB dest and the Box subbox, into dest, starting at component dest in dest.
void shift (const IntVect &v)
 Perform shifts on the FabArray.
int nGrow () const
 Returns the grow factor that defines the region of definition.
int nComp () const
 Returns number of variables associated with each point (nvar).
const BoxArrayboxArray () const
 Returns a constant reference to the BoxArray that defines the valid region associated with this FabArray.
virtual const Boxbox (int K) const
 Returns a constant reference to the Kth Box in the BoxArray.
virtual Box fabbox (int K) const
 Returns the Kth FABs Box in the FabArray.
int size () const
 Returns the number of FABs in the FabArray..
const DistributionMappingDistributionMap () const
 Returns constant reference to associated DistributionMapping.

Protected Attributes

PArray< FAB > fabparray
BoxArray boxarray
DistributionMapping distributionMap
int n_grow
int n_comp

template<class FAB>
class FabArray< FAB >


Member Typedef Documentation

template<class FAB>
typedef FAB::value_type FabArray< FAB >::value_type


Constructor & Destructor Documentation

template<class FAB>
FabArray< FAB >::FabArray (  )  [inline]

Constructs an empty FabArray<FAB>.

template<class FAB>
FabArray< FAB >::FabArray ( const BoxArray bxs,
int  nvar,
int  ngrow,
FabAlloc  mem_mode = Fab_allocate 
) [inline]

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.

If mem is defined to be Fab then FABs are allocated for each Box in the BoxArray. The size of the Kth FAB is given by bxs[K] grown by ngrow. If mem is defined to be Fab, then no FABs are allocated at this time, but can be defined later. The number of components in each FAB is not specified and is expected to be implicit in the definition of the FAB class. That is, the FAB constructor will take only a Box argument. Call this constructor number two.

References FabArray< FAB >::define().

template<class FAB>
FabArray< FAB >::FabArray ( const BoxArray bxs,
int  nvar,
int  ngrow,
const DistributionMapping dm,
FabAlloc  mem_mode 
) [inline]

template<class FAB>
FabArray< FAB >::~FabArray (  )  [inline, virtual]

The (virtual) destructor -- deletes all FABs in the array.


Member Function Documentation

template<class FAB>
void FabArray< FAB >::define ( const BoxArray bxs,
int  nvar,
int  ngrow,
FabAlloc  mem_mode 
) [inline]

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.

References BL_ASSERT, FabArrayBase::boxarray, DistributionMapping::define(), BoxArray::define(), FabArrayBase::distributionMap, Fab_allocate, FabArray< FAB >::fabparray, FabArrayBase::n_comp, FabArrayBase::n_grow, ParallelDescriptor::NProcsCFD(), PArray< T >::resize(), and BoxArray::size().

Referenced by FabArray< FAB >::FabArray(), and VisMF::Read().

template<class FAB>
void FabArray< FAB >::define ( const BoxArray bxs,
int  nvar,
int  ngrow,
const DistributionMapping dm,
FabAlloc  mem_mode 
) [inline]

template<class FAB>
bool FabArray< FAB >::ok (  )  const [inline]

Returns true if the FabArray is well-defined.

That is, if 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.

References FabArrayBase::box(), PArray< T >::defined(), FabArray< FAB >::fabparray, BoxLib::grow(), FabArrayBase::n_grow, and ParallelDescriptor::ReduceLongAnd().

Referenced by VisMF::Read().

template<class FAB>
const FAB & FabArray< FAB >::operator[] ( const MFIter mfi  )  const [inline]

Returns a constant reference to the FAB associated with the Kth element.

References FabArray< FAB >::fabparray, and MFIter::index().

template<class FAB>
const FAB & FabArray< FAB >::get ( const MFIter mfi  )  const [inline]

template<class FAB>
FAB & FabArray< FAB >::operator[] ( const MFIter mfi  )  [inline]

Returns a reference to the FAB associated with the Kth element.

References FabArray< FAB >::fabparray, and MFIter::index().

template<class FAB>
FAB & FabArray< FAB >::get ( const MFIter mfi  )  [inline]

template<class FAB>
const FAB & FabArray< FAB >::operator[] ( int  K  )  const [inline]

Returns a constant reference to the FAB associated with the Kth element.

References FabArray< FAB >::fabparray.

template<class FAB>
FAB & FabArray< FAB >::operator[] ( int  K  )  [inline]

Returns a reference to the FAB associated with the Kth element.

References FabArray< FAB >::fabparray.

template<class FAB>
void FabArray< FAB >::setFab ( int  K,
FAB *  elem 
) [inline]

template<class FAB>
void FabArray< FAB >::clear (  )  [inline]

Releases FAB memory in the FabArray.

References PArray< T >::clear(), and FabArray< FAB >::fabparray.

template<class FAB>
void FabArray< FAB >::setVal ( value_type  val  )  [inline]

Set all components in the valid region of each FAB to val.

Referenced by FabArray< FAB >::operator=(), and FabArray< FAB >::setVal().

template<class FAB>
void FabArray< FAB >::operator= ( const value_type val  )  [inline]

template<class FAB>
void FabArray< FAB >::setVal ( value_type  val,
int  comp,
int  num_comp,
int  nghost = 0 
) [inline]

Set the value of num 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.

References BL_ASSERT, BoxLib::grow(), FabArrayBase::n_comp, FabArrayBase::n_grow, and FabArray< FAB >::setVal().

template<class FAB>
void FabArray< FAB >::setVal ( value_type  val,
const Box region,
int  comp,
int  num_comp,
int  nghost = 0 
) [inline]

Set the value of num 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.

References BL_ASSERT, BoxLib::grow(), FabArrayBase::n_comp, FabArrayBase::n_grow, and Box::ok().

template<class FAB>
void FabArray< FAB >::setVal ( value_type  val,
int  nghost 
) [inline]

Set all components in the valid region of each FAB in the FabArray to val, including nghost boundary cells.

References FabArrayBase::n_comp, and FabArray< FAB >::setVal().

template<class FAB>
void FabArray< FAB >::setVal ( value_type  val,
const Box region,
int  nghost 
) [inline]

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.

References FabArrayBase::n_comp, and FabArray< FAB >::setVal().

template<class FAB>
void FabArray< FAB >::setBndry ( value_type  val  )  [inline]

Set all values in the boundary region to val.

References FabArrayBase::n_comp.

template<class FAB>
void FabArray< FAB >::setBndry ( value_type  val,
int  strt_comp,
int  ncomp 
) [inline]

Set ncomp values in the boundary region, starting at start to val.

References FabArrayBase::n_grow.

template<class FAB>
void FabArray< FAB >::copy ( const FabArray< FAB > &  fa  )  [inline]

This function copies data from fa 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 region of each FAB.

References FabArrayBase::nComp().

Referenced by FabArray< FAB >::copy().

template<class FAB>
void FabArray< FAB >::copy ( const FabArray< FAB > &  src,
int  src_comp,
int  dest_comp,
int  num_comp 
) [inline]

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 components starting at src in the FabArray src, with the destination components in this FabArray starting at dest. This assumes that the source and destination FabArray have identical valid regions.

References Arena::alloc(), ParallelDescriptor::Arecv(), BL_ASSERT, FabArrayBase::box(), FabComTag::box, FabArrayBase::boxarray, FabArray< FAB >::copy(), FabArrayBase::distributionMap, FabComTag::fabIndex, FabArray< FAB >::fabparray, Arena::free(), BoxLib::max(), MPI_REQUEST_NULL, ParallelDescriptor::MyProc(), ParallelDescriptor::NProcs(), Box::numPts(), ParallelDescriptor::Send(), ParallelDescriptor::SeqNum(), Array< T >::size(), BoxArray::size(), FabArrayBase::size(), BoxLib::The_Arena(), and ParallelDescriptor::Waitsome().

template<class FAB>
void FabArray< FAB >::copy ( FAB &  dest  )  const [inline]

Copies the values contained in the intersection of the valid region of this FabArray with the FAB dest into dest.

References FabArray< FAB >::copy().

template<class FAB>
void FabArray< FAB >::copy ( FAB &  dest,
const Box subbox 
) const [inline]

Copies the values contained in the intersection of the valid region of this FabArray with the FAB dest and the Box subbox into that subregion of dest.

References FabArray< FAB >::copy().

template<class FAB>
void FabArray< FAB >::copy ( FAB &  dest,
int  src_comp,
int  dest_comp,
int  num_comp 
) const [inline]

Copies the values contained in the intersection of the num component valid region of this FabArray, starting at component src, with the FAB dest into dest, starting at component dest in dest.

References FabArray< FAB >::copy().

template<class FAB>
void FabArray< FAB >::copy ( FAB &  dest,
const Box subbox,
int  src_comp,
int  dest_comp,
int  num_comp 
) const [inline]

Copies the values contained in the intersection of the num component valid region of this FabArray, starting at component src, with the FAB dest and the Box subbox, into dest, starting at component dest in dest.

References ParallelDescriptor::Bcast(), BL_ASSERT, FabArrayBase::boxarray, FabArrayBase::distributionMap, FabArray< FAB >::fabparray, Box::intersects(), ParallelDescriptor::MyProc(), ParallelDescriptor::NProcs(), Box::numPts(), BoxArray::resize(), and FabArrayBase::size().

template<class FAB>
void FabArray< FAB >::shift ( const IntVect v  )  [inline]

Perform shifts on the FabArray.

References FabArrayBase::boxarray, and BoxArray::shift().

int FabArrayBase::nGrow (  )  const [inline, inherited]

Returns the grow factor that defines the region of definition.

References FabArrayBase::n_grow.

Referenced by BuildFBsirec(), MultiFab::Copy(), and TheFBsirec().

int FabArrayBase::nComp (  )  const [inline, inherited]

Returns number of variables associated with each point (nvar).

References FabArrayBase::n_comp.

Referenced by FabArray< FAB >::copy(), and VisMF::Write().

const BoxArray & FabArrayBase::boxArray (  )  const [inline, inherited]

Returns a constant reference to the BoxArray that defines the valid region associated with this FabArray.

References FabArrayBase::boxarray.

Referenced by BuildFBsirec(), MultiFab::Copy(), TheFBsirec(), and VisMF::Write().

const Box & FabArrayBase::box ( int  K  )  const [inline, virtual, inherited]

Returns a constant reference to the Kth Box in the BoxArray.

That is, the valid region of the Kth grid.

References FabArrayBase::boxarray.

Referenced by FabArrayCopyDescriptor< FAB >::AddBoxDoIt(), BuildFBsirec(), FabArray< FAB >::copy(), FabArray< FAB >::ok(), and MFIter::validbox().

Box FabArrayBase::fabbox ( int  K  )  const [virtual, inherited]

Returns the Kth FABs Box in the FabArray.

That is, the region the Kth fab is actually defined on.

References FabArrayBase::boxarray, BoxLib::grow(), and FabArrayBase::n_grow.

Referenced by FabArrayCopyDescriptor< FAB >::AddBoxDoIt(), and MFIter::fabbox().

int FabArrayBase::size (  )  const [inline, inherited]

Returns the number of FABs in the FabArray..

References FabArrayBase::boxarray, and BoxArray::size().

Referenced by BuildFBsirec(), FabArray< FAB >::copy(), and MFIter::isValid().

const DistributionMapping & FabArrayBase::DistributionMap (  )  const [inline, inherited]


Member Data Documentation

template<class FAB>
PArray<FAB> FabArray< FAB >::fabparray [protected]

BoxArray FabArrayBase::boxarray [mutable, protected, inherited]

int FabArrayBase::n_grow [protected, inherited]

int FabArrayBase::n_comp [protected, inherited]


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

Generated on Fri Nov 21 10:11:01 2008 for AMRParticlePaths by  doxygen 1.5.5