|
| | PatchCellDataMiscellaneousOpsReal () |
| |
| virtual | ~PatchCellDataMiscellaneousOpsReal () |
| |
| int | computeConstrProdPos (const tbox::Pointer< pdat::CellData< DIM, TYPE > > &data1, const tbox::Pointer< pdat::CellData< DIM, TYPE > > &data2, const hier::Box< DIM > &box, const tbox::Pointer< pdat::CellData< DIM, double > > cvol=((pdat::CellData< DIM, double > *) NULL)) const |
| |
| void | compareToScalar (tbox::Pointer< pdat::CellData< DIM, TYPE > > &dst, const tbox::Pointer< pdat::CellData< DIM, TYPE > > &src, const TYPE &alpha, const hier::Box< DIM > &box, const tbox::Pointer< pdat::CellData< DIM, double > > cvol=((pdat::CellData< DIM, double > *) NULL)) const |
| |
| int | testReciprocal (tbox::Pointer< pdat::CellData< DIM, TYPE > > &dst, const tbox::Pointer< pdat::CellData< DIM, TYPE > > &src, const hier::Box< DIM > &box, const tbox::Pointer< pdat::CellData< DIM, double > > cvol=((pdat::CellData< DIM, double > *) NULL)) const |
| |
| TYPE | maxPointwiseDivide (const tbox::Pointer< pdat::CellData< DIM, TYPE > > &numer, const tbox::Pointer< pdat::CellData< DIM, TYPE > > &denom, const hier::Box< DIM > &box) const |
| | Compute max of "conditional" quotients of two arrays. More...
|
| |
| TYPE | minPointwiseDivide (const tbox::Pointer< pdat::CellData< DIM, TYPE > > &numer, const tbox::Pointer< pdat::CellData< DIM, TYPE > > &denom, const hier::Box< DIM > &box) const |
| | Compute min of quotients of two arrays. More...
|
| |
template<int DIM, class TYPE>
class SAMRAI::math::PatchCellDataMiscellaneousOpsReal< DIM, TYPE >
Class PatchCellDataMiscellaneousOpsReal<DIM> provides access to a collection of operations that may be applied to numerical cell-centered patch data of type double and float. The primary intent of this class is to provide the interface to these operations for the class PatchCellDataOpsReal<DIM> which provides access to a more complete set of operations that may be used to manipulate cell-centered patch data. Each member function accepts a box argument indicating the region of index space on which the operation should be performed. The operation will be performed on the intersection of this box and those boxes corresponding to the patch data objects. Also, each operation allows an additional cell-centered patch data object to be used to represent a control volume that weights the contribution of each data entry in the given norm calculation. Note that the control volume patch data must be of type double and have cell-centered geometry (i.e., the same as the data itself). The use of control volumes is important when these operations are used in vector kernels where the data resides over multiple levels of spatial resolution in an AMR hierarchy. If the control volume is not given in the function call, it will be ignored in the calculation. Also, note that the depth of the control volume patch data object must be either 1 or be equal to the depth of the other data objects.
Since these operations are used only by the vector kernels for the KINSOL and CVODE solver packages at this time, they are intended to be instantiated for the standard built-in types double and float (since those solvers only treat double and float data). To extend this class to other data types or to include other operations, the member functions must be specialized or the new operations must be added.
- See also
- math::ArrayDataMiscellaneousOpsReal
template<int DIM, class TYPE >
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!