SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE > Class Template Reference

#include <source/mathops/hierarchy/HierarchyDataOpsReal.h>

Inheritance diagram for SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >:

Inheritance graph
[legend]
List of all members.

Public Member Functions

 HierarchyDataOpsReal ()
virtual ~HierarchyDataOpsReal ()
virtual void setPatchHierarchy (tbox::Pointer< hier::PatchHierarchy< DIM > > hierarchy)=0
virtual void resetLevels (const int coarsest_level, const int finest_level)=0
virtual const tbox::Pointer<
hier::PatchHierarchy< DIM > > 
getPatchHierarchy () const=0
virtual void copyData (const int dst_id, const int src_id, const bool interior_only=true) const =0
virtual void swapData (const int data1_id, const int data2_id) const=0
virtual void printData (const int data_id, std::ostream &s, const bool interior_only=true) const =0
virtual void setToScalar (const int data_id, const TYPE &alpha, const bool interior_only=true) const =0
virtual void scale (const int dst_id, const TYPE &alpha, const int src_id, const bool interior_only=true) const =0
virtual void addScalar (const int dst_id, const int src_id, const TYPE &alpha, const bool interior_only=true) const =0
virtual void add (const int dst_id, const int src1_id, const int src2_id, const bool interior_only=true) const =0
virtual void subtract (const int dst_id, const int src1_id, const int src2_id, const bool interior_only=true) const =0
virtual void multiply (const int dst_id, const int src1_id, const int src2_id, const bool interior_only=true) const =0
virtual void divide (const int dst_id, const int src1_id, const int src2_id, const bool interior_only=true) const =0
virtual void reciprocal (const int dst_id, const int src_id, const bool interior_only=true) const =0
virtual void linearSum (const int dst_id, const TYPE &alpha, const int src1_id, const TYPE &beta, const int src2_id, const bool interior_only=true) const =0
virtual void axpy (const int dst_id, const TYPE &alpha, const int src1_id, const int src2_id, const bool interior_only=true) const =0
virtual void axmy (const int dst_id, const TYPE &alpha, const int src1_id, const int src2_id, const bool interior_only=true) const =0
virtual void abs (const int dst_id, const int src_id, const bool interior_only=true) const =0
virtual TYPE min (const int data_id, const bool interior_only=true) const =0
virtual TYPE max (const int data_id, const bool interior_only=true) const =0
virtual void setRandomValues (const int data_id, const TYPE &width, const TYPE &low, const bool interior_only=true) const =0
virtual int numberOfEntries (const int data_id, const bool interior_only=true) const =0
virtual double sumControlVolumes (const int data_id, const int vol_id) const =0
virtual double L1Norm (const int data_id, const int vol_id=-1, bool local_only=false) const =0
virtual double L2Norm (const int data_id, const int vol_id=-1, bool local_only=false) const =0
virtual double weightedL2Norm (const int data_id, const int weight_id, const int vol_id=-1) const =0
virtual double RMSNorm (const int data_id, const int vol_id=-1) const=0
virtual double weightedRMSNorm (const int data_id, const int weight_id, const int vol_id=-1) const =0
virtual double maxNorm (const int data_id, const int vol_id=-1, bool local_only=false) const =0
virtual TYPE dot (const int data1_id, const int data2_id, const int vol_id=-1, bool local_only=false) const=0
virtual int computeConstrProdPos (const int data1_id, const int data2_id, const int vol_id=-1) const =0
virtual void compareToScalar (const int dst_id, const int src_id, const TYPE &alpha, const int vol_id=-1) const=0
virtual int testReciprocal (const int dst_id, const int src_id, const int vol_id=-1) const =0
virtual TYPE maxPointwiseDivide (const int numer_id, const int denom_id, bool local_only=false) const =0
 Compute max of "conditional" quotients of two arrays.
virtual TYPE minPointwiseDivide (const int numer_id, const int denom_id, bool local_only=false) const =0
 Compute min of quotients of two arrays.

Detailed Description

template<int DIM, class TYPE>
class SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >

Class HierarchyDataOpsReal<DIM> defines the interface to a collection of operations that may be used to manipulate real numerical (float or double) patch data components over multiple levels in an AMR hierarchy. It serves as a base class for subclasses which implement the operations for cell-centered, face-centered, or node-centered data types. The patch hierarchy and set of levels within that hierarcy over which the operations will be performed are set in the constructor of the subclass. However, these data members may be changed at any time via the virtual access functions setPatchHierarchy() and resetLevels() below. The operations include basic arithmetic, norms and ordering, and assorted miscellaneous operations.

Note that, when it makes sense, an operation accept a boolean argument which indicates whether the operation should be performed on all of the data or just those data elements corresponding to the patch interiors. If no boolean argument is provided, the default behavior is to treat only the patch interiors. Also, interfaces for similar sets of operations for complex and integer hierarchy data are defined in the classes HierarchyDataOpsComplex<DIM> and HierarchyDataOpsInteger<DIM>, respectively.


Constructor & Destructor Documentation

template<int DIM, class TYPE>
SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::HierarchyDataOpsReal (  ) 

The constructor for the HierarchyDataOpsReal<DIM> class.

template<int DIM, class TYPE>
SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::~HierarchyDataOpsReal (  )  [virtual]

Virtual destructor for the HierarchyDataOpsReal<DIM> class.


Member Function Documentation

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::setPatchHierarchy ( tbox::Pointer< hier::PatchHierarchy< DIM > >  hierarchy  )  [pure virtual]

Reset patch hierarchy over which operations occur.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::resetLevels ( const int  coarsest_level,
const int  finest_level 
) [pure virtual]

Reset range of patch levels over which operations occur. Typically, levels must exist in hierarchy or an assertion will result.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual const tbox::Pointer< hier::PatchHierarchy<DIM> > SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::getPatchHierarchy (  )  const [pure virtual]

Return const pointer to patch hierarchy associated with operations.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::copyData ( const int  dst_id,
const int  src_id,
const bool  interior_only = true 
) const [pure virtual]

Copy source data to destination data.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::swapData ( const int  data1_id,
const int  data2_id 
) const [pure virtual]

Swap data pointers (i.e., storage) between two data components.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::printData ( const int  data_id,
std::ostream &  s,
const bool  interior_only = true 
) const [pure virtual]

Print data over multiple levels to specified output stream.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::setToScalar ( const int  data_id,
const TYPE &  alpha,
const bool  interior_only = true 
) const [pure virtual]

Set data component to given scalar.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::scale ( const int  dst_id,
const TYPE &  alpha,
const int  src_id,
const bool  interior_only = true 
) const [pure virtual]

Set destination to source multiplied by given scalar, pointwise.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::addScalar ( const int  dst_id,
const int  src_id,
const TYPE &  alpha,
const bool  interior_only = true 
) const [pure virtual]

Add scalar to each entry in source data and set destination to result.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::add ( const int  dst_id,
const int  src1_id,
const int  src2_id,
const bool  interior_only = true 
) const [pure virtual]

Set destination to sum of two source components, pointwise.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::subtract ( const int  dst_id,
const int  src1_id,
const int  src2_id,
const bool  interior_only = true 
) const [pure virtual]

Subtract second source component from first source component pointwise and set destination data component to result.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::multiply ( const int  dst_id,
const int  src1_id,
const int  src2_id,
const bool  interior_only = true 
) const [pure virtual]

Set destination component to product of two source components, pointwise.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::divide ( const int  dst_id,
const int  src1_id,
const int  src2_id,
const bool  interior_only = true 
) const [pure virtual]

Divide first data component by second source component pointwise and set destination data component to result.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::reciprocal ( const int  dst_id,
const int  src_id,
const bool  interior_only = true 
) const [pure virtual]

Set each entry of destination component to reciprocal of corresponding source data component entry.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::linearSum ( const int  dst_id,
const TYPE &  alpha,
const int  src1_id,
const TYPE &  beta,
const int  src2_id,
const bool  interior_only = true 
) const [pure virtual]

Set $d = \alpha s_1 + \beta s_2$, where $d$ is the destination patch data component and $s_1, s_2$ are the first and second source components, respectively. Here $\alpha, \beta$ are scalar values.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::axpy ( const int  dst_id,
const TYPE &  alpha,
const int  src1_id,
const int  src2_id,
const bool  interior_only = true 
) const [pure virtual]

Set $d = \alpha s_1 + s_2$, where $d$ is the destination patch data component and $s_1, s_2$ are the first and second source components, respectively. Here $\alpha$ is a scalar.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::axmy ( const int  dst_id,
const TYPE &  alpha,
const int  src1_id,
const int  src2_id,
const bool  interior_only = true 
) const [pure virtual]

Set $d = \alpha s_1 - s_2$, where $d$ is the destination patch data component and $s_1, s_2$ are the first and second source components, respectively. Here $\alpha$ is a scalar.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::abs ( const int  dst_id,
const int  src_id,
const bool  interior_only = true 
) const [pure virtual]

Set destination data to absolute value of source data, pointwise.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual TYPE SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::min ( const int  data_id,
const bool  interior_only = true 
) const [pure virtual]

Return minimum data value over all patches in the collection of levels.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual TYPE SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::max ( const int  data_id,
const bool  interior_only = true 
) const [pure virtual]

Return maximum data value over all patches in the collection of levels.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::setRandomValues ( const int  data_id,
const TYPE &  width,
const TYPE &  low,
const bool  interior_only = true 
) const [pure virtual]

Set data entries to random values. See the operations in the array data operation classes for each data type for details on the generation of the random values.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual int SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::numberOfEntries ( const int  data_id,
const bool  interior_only = true 
) const [pure virtual]

Return the total number of data values for the component on the set of hierarchy levels. If the boolean argument is true, the number of elements will be summed over patch interiors. If the boolean argument is false, all elements will be counted (including ghost values) over all patches.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual double SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::sumControlVolumes ( const int  data_id,
const int  vol_id 
) const [pure virtual]

Return sum of the control volumes associated with the data component. Note that if the ontrol volumes are set propery, this is equivalent to integrating a data component containing all ones over the collection of hierarchy levels.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual double SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::L1Norm ( const int  data_id,
const int  vol_id = -1,
bool  local_only = false 
) const [pure virtual]

Return discrete $L_1$-norm of the data using the control volume to weight the contribution of each data entry to the sum. That is, the return value is the sum $\sum_i ( \| data_i \| cvol_i )$. If the control volume is not defined (vol_id < 0), the return value is $\sum_i ( \| data_i \| )$. If local_only is true, the global reduction is not performed (thus each process will get only local results).

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual double SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::L2Norm ( const int  data_id,
const int  vol_id = -1,
bool  local_only = false 
) const [pure virtual]

Return discrete $L_2$-norm of the data using the control volume to weight the contribution of each data entry to the sum. That is, the return value is the sum $\sqrt{ \sum_i ( (data_i)^2 cvol_i ) }$. If the control volume is not defined (vol_id < 0), the return value is $\sqrt{ \sum_i ( (data_i)^2 cvol_i ) }$. If local_only is true, the global reduction is not performed (thus each process will get only local results).

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual double SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::weightedL2Norm ( const int  data_id,
const int  weight_id,
const int  vol_id = -1 
) const [pure virtual]

Return discrete weighted $L_2$-norm of the data using the control volume to weight the contribution of the data and weight entries to the sum. That is, the return value is the sum $\sqrt{ \sum_i ( (data_i * weight_i)^2 cvol_i ) }$. If the control volume is not defined (vol_id < 0), the return value is $\sqrt{ \sum_i ( (data_i * weight_i)^2 ) }$.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual double SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::RMSNorm ( const int  data_id,
const int  vol_id = -1 
) const [pure virtual]

Return discrete root mean squared norm of the data. If the control volume is defined, the return value is the $L_2$-norm divided by the square root of the sum of the control volumes. Otherwise, the return value is the $L_2$-norm divided by the square root of the number of data entries.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual double SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::weightedRMSNorm ( const int  data_id,
const int  weight_id,
const int  vol_id = -1 
) const [pure virtual]

Return discrete weighted root mean squared norm of the data. If the control volume is define, the return value is the weighted $L_2$-norm divided by the square root of the sum of the control volumes. Otherwise, the return value is the weighted $L_2$-norm divided by the square root of the number of data entries.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual double SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::maxNorm ( const int  data_id,
const int  vol_id = -1,
bool  local_only = false 
) const [pure virtual]

Return the $\max$-norm of the data using the control volume to weight the contribution of each data entry to the maximum. That is, the return value is $\max_i ( \| data_i \| )$, where the max is over the data elements where $cvol_i > 0$. If the control volume is undefined (vol_id < 0), it is ignored during the computation of the maximum. If local_only is true, the global reduction is not performed (thus each process will get only local results).

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual TYPE SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::dot ( const int  data1_id,
const int  data2_id,
const int  vol_id = -1,
bool  local_only = false 
) const [pure virtual]

Return the dot product of the two data arrays using the control volume to weight the contribution of each product to the sum. That is, the return value is the sum $\sum_i ( data1_i * data2_i * cvol_i )$. If the control volume is undefined (vol_id < 0), it is ignored during the summation. If local_only is true, the global reduction is not performed (thus each process will get only local results).

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual int SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::computeConstrProdPos ( const int  data1_id,
const int  data2_id,
const int  vol_id = -1 
) const [pure virtual]

Return 1 if $\|data2_i\| > 0$ and $data1_i * data2_i \leq 0$, for any $i$ in the set of patch data indices, where $cvol_i > 0$. Otherwise, return 0. If the control volume is undefined (vol_id < 0), all values on the patch interiors are considered.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual void SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::compareToScalar ( const int  dst_id,
const int  src_id,
const TYPE &  alpha,
const int  vol_id = -1 
) const [pure virtual]

Wherever $cvol_i > 0$ in the set of patch data indices, set $dst_i = 1$ if $\|src_i\| > \alpha$, and $dst_i = 0$ otherwise. If the control volume is undefined (vol_id < 0), all values on the patch interiors are considered.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual int SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::testReciprocal ( const int  dst_id,
const int  src_id,
const int  vol_id = -1 
) const [pure virtual]

Wherever $cvol_i > 0$ in the set of patch data indices, set $dst_i = 1/src_i$ if $src_i \neq 0$, and $dst_i = 0$ otherwise. If $dst_i = 0$ anywhere, 0 is the return value. Otherwise 1 is returned. If the control volume is undefined (vol_id < 0), all values on the patch interiors are considered.

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual TYPE SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::maxPointwiseDivide ( const int  numer_id,
const int  denom_id,
bool  local_only = false 
) const [pure virtual]

Compute max of "conditional" quotients of two arrays.

Return the maximum of pointwise "conditional" quotients of the numerator and denominator.

The "conditional" quotient is defined as |numerator/denominator| where the denominator is nonzero. Otherwise, it is defined as |numerator|.

Note: This method is currently intended to support the PETSc-2.1.6 vector wrapper only. Please do not use it!

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.

template<int DIM, class TYPE>
virtual TYPE SAMRAI::math::HierarchyDataOpsReal< DIM, TYPE >::minPointwiseDivide ( const int  numer_id,
const int  denom_id,
bool  local_only = false 
) const [pure virtual]

Compute min of quotients of two arrays.

Return the minimum of pointwise quotients of the numerator and denominator.

The quotient is defined as (numerator/denominator) where the denominator is nonzero. When the denominator is zero, the entry is skipped. If the denominator is always zero, the value of tbox::IEEE::getFLT_MAX() is returned (see SAMRAI::tbox::IEEE).

Note: This method is currently intended to support the SUNDIALS vector wrapper only. Please do not use it!

Implemented in SAMRAI::math::HierarchyCellDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyEdgeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyFaceDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyNodeDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchySideDataOpsReal< DIM, TYPE >, SAMRAI::math::HierarchyCellDataOpsReal< DIM, double >, and SAMRAI::math::HierarchySideDataOpsReal< DIM, double >.


The documentation for this class was generated from the following files:
Generated on Thu Jun 18 11:28:36 2009 for SAMRAI by  doxygen 1.5.1