RTOp Package Browser (Single Doxygen Collection)
Version of the Day
|
Classes | |
class | ReductTarget |
Abstract base class for all reduction objects. More... | |
class | RTOpT |
Templated interface to vector reduction/transformation operators {abstract}. More... | |
class | UnknownError |
class | InvalidUsage |
class | InvalidNumVecs |
class | InvalidNumTargVecs |
class | IncompatibleVecs |
class | IncompatibleReductObj |
class | ConstSubVectorView |
Class for a non-changeable sub-vector. More... | |
class | SubVectorView |
Class for a changeable sub-vector. More... | |
class | ConstSubMultiVectorView |
Class for a non-changeable sub-multi-vector (submatrix). More... | |
class | SubMultiVectorView |
Class for a changeable sub-vector. More... | |
class | PrimitiveTypeTraits |
A templated traits class for decomposing object into an array of primitive objects. More... | |
class | PrimitiveTypeTraits< Scalar, Scalar > |
Specialization where the scalar type is the same as the concrete object type. More... | |
class | PrimitiveTypeTraits< Scalar, index_type > |
Specialization for index_type concrete object. More... | |
class | ROpGetElementEleWiseReductionOp |
class | ROpGetElement |
Returns the value of the element at index global_i . More... | |
class | ROpGetSubVector |
Reduction operator that extracts a sub-vector in the range of global zero-based indexes [l,u]. More... | |
class | ROpMaxIndexEleWiseReductionOp |
class | ROpMaxIndexReductObjReductionOp |
class | ROpMaxIndex |
Returns the maximum element and its index: result.scalar = x(k) and result.index = k such that x(k) >= x(i) for i=0...n-1 and k is the minimum index to break ties. More... | |
class | ROpMaxIndexLessThanBoundEleWiseReductionOp |
class | ROpMaxIndexLessThanBound |
Returns the maximum element less than some bound along with its index: result.scalar = x(k) and result.index = k such that x(k) >= x(i) for all i where x(i) < bound and k is the minimum index to break ties. More... | |
class | ROpMinIndexEleWiseReductionOp |
class | ROpMinIndexReductObjReductionOp |
class | ROpMinIndex |
Returns the minimum element and its index: result.scalar = x(k) and result.index = k such that x(k) <= x(i) for i=0...n-1 and k is the minimum index to break ties. More... | |
class | ROpMinIndexGreaterThanBoundEleWiseReductionOp |
class | ROpMinIndexGreaterThanBound |
Returns the minimum element greater than some bound along with its index: result.scalar = x(k) and result.index = k such that x(k) <= x(i) for all i where x(i) > bound and k is the minimum index to break ties. More... | |
class | ROpNorm2EleWiseReduction |
class | ROpNorm2 |
Two (Euclidean) norm reduction operator: result = sqrt( sum( conj(v0[i])*v0[i], i=0...n-1 ) ) . More... | |
class | ROpWeightedNorm2EleWiseReduction |
class | ROpWeightedNorm2 |
Weighted Two (Euclidean) norm reduction operator: result = sqrt( sum( v0[i]*conj(v1[i])*v1[i], i=0...n-1 ) ) . More... | |
class | TOpAddScalarEleWiseTransformation |
Element-wise transformation operator for TOpAddScalar. More... | |
class | TOpAddScalar |
Add a scalar to a vector transformation operator: z0[i] += alpha, i=0...n-1 . More... | |
class | TOpAssignScalarEleWiseTransformation |
Element-wise transformation operator for TOpAssignScalar. More... | |
class | TOpAssignScalar |
Assign a scalar to a vector transformation operator: z0[i] = alpha, i=0...n-1 . More... | |
class | TOpAXPYEleWiseTransformation |
Element-wise transformation operator for TOpAXPY. More... | |
class | TOpAXPY |
AXPY transformation operator: z0[i] += alpha*v0[i], i=0...n-1 . More... | |
class | TOpEleWiseConjProdEleWiseTransformation |
Element-wise transformation operator for TOpEleWiseConjProd. More... | |
class | TOpEleWiseConjProd |
Element-wise product transformation operator: z0[i] += alpha*conj(v0[i])*v1[i], i=0...n-1 . More... | |
class | TOpEleWiseDivideEleWiseTransformation |
Element-wise transformation operator for TOpEleWiseDivide. More... | |
class | TOpEleWiseDivide |
Element-wise division transformation operator: z0[i] += alpha*v0[i]/v1[i], i=0...n-1 . More... | |
class | TOpEleWiseProdEleWiseTransformation |
Element-wise transformation operator for TOpEleWiseProd. More... | |
class | TOpEleWiseProd |
Element-wise product transformation operator: z0[i] += alpha*v0[i]*v1[i], i=0...n-1 . More... | |
class | TOpEleWiseProdUpdateEleWiseTransformation |
Element-wise transformation operator for TOpEleWiseProdUpdate. More... | |
class | TOpEleWiseProdUpdate |
Element-wise product update transformation operator: z0[i] *= alpha*v0[i], i=0...n-1 . More... | |
class | TOpEleWiseScaleEleWiseTransformation |
Element-wise vector scaling op for TOpEleWiseScaling. More... | |
class | TOpEleWiseScale |
Element-wise vector scaling: z0[i] *= v0[i], i=0...n-1 . More... | |
class | TOpLinearCombination |
Linear combination transformation operator: z0[i] = beta*z0[i]. More... | |
class | TOpPairWiseMaxPairWiseTransformation |
Pair-wise transformation operator for TOpPairWiseMax. More... | |
class | TOpPairWiseMax |
Pair-wise Maximum transformation operator: z0[i] = alpha*max(v0[i],v1[i]), i=0...n-1 . More... | |
class | TOpPairWiseMaxUpdatePairWiseTransformation |
Pair-wise transformation operator for TOpPairWiseMaxUpdate. More... | |
class | TOpPairWiseMaxUpdate |
Pair-wise Maximum update transformation operator: z0[i] = alpha*max(z0[i],v0[i]), i=0...n-1 . More... | |
class | TOpRandomize |
Generate a random vector in the range [l,u]: z0[i] = 0.5*((u-l)*Teuchos::ScalarTraits<Scalar>::random()+(u+l)), i=0...n-1 . More... | |
class | TOpScaleVectorEleWiseTransformation |
Element-wise transformation operator for TOpScaleVector. More... | |
class | TOpScaleVector |
Simple transformation operator that scales every vector element by a scalar alpha . More... | |
class | TOpSetAssendingValuesEleWiseTransformation |
Element-wise transformation for TOpSetAssendingValues. More... | |
class | TOpSetAssendingValues |
Set the elements of a vector to: z0[i] = i+offset+1, i=0...n-1 . More... | |
class | TOpSetElementEleWiseTransformation |
Element-wise transformation for TOpSetElement. More... | |
class | TOpSetElement |
Set the elements of a vector to: z0[i] = i+global_i+1, i=0...n-1 . More... | |
class | TOpSetSubVector |
Advanced transformation operator that assigns elements from a sparse explicit vector. More... | |
class | RTOpServer |
Server for creating RTOpT objects given just the operators name. More... | |
class | RTOpSubRangeDecorator |
Decorator subclass that restricts the range of elements to apply the underlying RTOpT object to. More... | |
struct | ScalarIndex |
Simple struct for a Scalar and an Ordinal object. More... | |
class | PrimitiveTypeTraits< Scalar, ScalarIndex< Scalar > > |
Partial specialization of PrimitiveTypeTraits for ScalarIndex . More... | |
class | DefaultReductTarget |
Simple ReductTarget subclass for simple scalar objects. More... | |
class | BasicReductObjReductionOp |
class | BasicReductObjReductionOp< ConcreteReductObj, REDUCT_TYPE_SUM > |
class | BasicReductObjReductionOp< ConcreteReductObj, REDUCT_TYPE_MAX > |
class | BasicReductObjReductionOp< ConcreteReductObj, REDUCT_TYPE_MIN > |
class | SumScalarReductObjReduction |
Null reduction object reduction operator. More... | |
class | ROpScalarReductionWithOpBase |
class | ROp_1_ScalarReduction |
Base class for scalar reduction RTOps with one input vector. More... | |
class | ROp_1_CoordVariantScalarReduction |
Base class for coordinate-variant scalar reduction RTOps with one input vector. More... | |
class | ROp_2_ScalarReduction |
Base class for scalar reduction RTOps with two input vectors. More... | |
class | TOp_0_1_Base |
Base class for transformations for 0 input and 1 output vector. More... | |
class | TOp_0_1_CoordVariantBase |
Base class for coordinate variant transformations for 0 input and 1 output vector. More... | |
class | TOp_1_1_Base |
Base class for transformations for 1 input and 1 output vector. More... | |
class | TOp_2_1_Base |
Base class for transformations for 2 input and 1 output vector. More... | |
class | TOp_3_1_Base |
Base class for transformations for 3 input and 1 output vector. More... | |
class | SparseSubVectorT |
Class for a (sparse or dense) sub-vector. More... | |
class | ReductTargetSerializer |
Serializer subclass for ReductTarget objects. More... | |
class | ReductTargetReductionOp |
ReductionOp subclass for ReductTarget objects. More... | |
class | TOpUnaryFuncPtr |
RTOpT subclass for unary transformation functions using a function pointer. More... | |
Typedefs | |
typedef Teuchos_Ordinal | Ordinal |
typedef Teuchos_Ordinal | index_type |
typedef char | char_type |
Enumerations | |
enum | ETransp { NOTRANS, TRANS, CONJTRANS } |
enum | EBasicReductTypes { REDUCT_TYPE_SUM, REDUCT_TYPE_MAX, REDUCT_TYPE_MIN } |
Functions | |
template<class Scalar > | |
void | assign_entries (const Ptr< const SubVectorView< Scalar > > &msv, const ConstSubVectorView< Scalar > &sv) |
template<class Scalar > | |
void | assign_entries (const Ptr< const SubMultiVectorView< Scalar > > &msmv, const ConstSubMultiVectorView< Scalar > &smv) |
template<class Scalar > | |
void | getrf (const SubMultiVectorView< Scalar > &A, const ArrayView< int > &ipiv, const Ptr< int > &rank) |
Peform an in-place factorization of a square or rectangular matrix. More... | |
template<class Scalar > | |
void | getrs (const ConstSubMultiVectorView< Scalar > &A, const ArrayView< const int > &ipiv, const ETransp transp, const Ptr< const SubMultiVectorView< Scalar > > &BX) |
RTOP_ROP_1_REDUCT_SCALAR (ROpCountNanInf, index_type, REDUCT_TYPE_SUM) | |
Reduction operator that counts the number of entries that are NaN or Inf. More... | |
RTOP_ROP_2_REDUCT_SCALAR (ROpDotProd, Scalar, REDUCT_TYPE_SUM) | |
RTOP_ROP_1_REDUCT_SCALAR_CUSTOM_DEFAULT (ROpMax, Scalar, REDUCT_TYPE_MAX, Teuchos::as< Scalar >(-std::numeric_limits< Scalar >::max())) | |
Maximum element: result = max{ v0[i], i=0...n-1 } . More... | |
RTOP_ROP_1_REDUCT_SCALAR_CUSTOM_DEFAULT (ROpMin, Scalar, REDUCT_TYPE_MIN, std::numeric_limits< Scalar >::max()) | |
Minimum element: result = min{ v0[i], i=0...n-1 } . More... | |
RTOP_ROP_1_REDUCT_SCALAR (ROpNorm1, typename ScalarTraits< Scalar >::magnitudeType, REDUCT_TYPE_SUM) | |
Class ROpNorm1. More... | |
RTOP_ROP_1_REDUCT_SCALAR (ROpNormInf, typename ScalarTraits< Scalar >::magnitudeType, REDUCT_TYPE_MAX) | |
Norm Inf: result = max(|x[i]|, i=0...n-1). More... | |
RTOP_ROP_1_REDUCT_SCALAR (ROpSum, Scalar, REDUCT_TYPE_SUM) | |
Class ROpSum: result = sum( v0[i], i=0...n-1 ) More... | |
RTOP_TOP_1_1 (TOpAbs) | |
Transformation operator that takes absolute values of elements: z0[i] = abs(v0[i]), i=0...n-1 . More... | |
RTOP_TOP_1_1 (TOpAssignVectors) | |
VectorBase assignment transformation operator: z0[i] = v0[i], i=0...n-1 . More... | |
RTOP_TOP_1_1 (TOpReciprocal) | |
VectorBase assignment transformation operator: z0[i] = v0[i], i=0...n-1 . More... | |
std::string | version () |
Print the version of RTOp. More... | |
template<class Scalar > | |
void | validate_apply_op (const RTOpT< Scalar > &op, const int allowed_num_sub_vecs, const int allowed_num_targ_sub_vecs, const bool expect_reduct_obj, const ArrayView< const ConstSubVectorView< Scalar > > &sub_vecs, const ArrayView< const SubVectorView< Scalar > > &targ_sub_vecs, const Ptr< const ReductTarget > &reduct_obj) |
Validate the input to an apply_op(...) function. More... | |
void | set_SPMD_apply_op_dump_out (const RCP< FancyOStream > &dumpOut) |
Set up to show a dump of RTOps applied through SPMD_apply_op(). More... | |
template<class PrimitiveScalar > | |
int | serializedSize (int num_values, int num_indexes, int num_chars) |
Return the size in bytes of an external representation of a ReductTarget object. More... | |
template<class Scalar > | |
void | serialize (const RTOpT< Scalar > &op, Ordinal num_values, Ordinal num_indexes, Ordinal num_chars, const ReductTarget &reduct_obj, char reduct_obj_ext[]) |
Serialize a ReductTarget object. More... | |
template<class Scalar > | |
void | deserialize (const RTOpT< Scalar > &op, int num_values, int num_indexes, int num_chars, const char reduct_obj_ext[], ReductTarget *reduct_obj) |
Deserialize a ReductTarget object. More... | |
template<class Scalar > | |
void | SPMD_all_reduce (const Teuchos::Comm< index_type > *comm, const RTOpT< Scalar > &op, const int num_cols, const ReductTarget *const i_reduct_objs[], ReductTarget *const reduct_objs[]) |
Reduce a set of reduction objects. More... | |
template<class Scalar > | |
void | SPMD_apply_op (const Teuchos::Comm< index_type > *comm, const RTOpT< Scalar > &op, const int num_vecs, const ConstSubVectorView< Scalar > sub_vecs[], const int num_targ_vecs, const SubVectorView< Scalar > targ_sub_vecs[], ReductTarget *reduct_obj) |
Apply an RTOp in SMPD mode to a set of vectors with contiguous storage per process. More... | |
template<class Scalar > | |
void | SPMD_apply_op (const Teuchos::Comm< index_type > *comm, const RTOpT< Scalar > &op, const int num_cols, const int num_multi_vecs, const ConstSubMultiVectorView< Scalar > sub_multi_vecs[], const int num_targ_multi_vecs, const SubMultiVectorView< Scalar > targ_sub_multi_vecs[], ReductTarget *const reduct_objs[]) |
Apply an RTOp in SMPD mode to a set of columns to a set of multi-vectors with contiguous storage per process. More... | |
template<class Scalar > | |
void | SPMD_apply_op (const Teuchos::Comm< index_type > *comm, const RTOpT< Scalar > &op, const int num_cols, const int num_vecs, const ConstSubVectorView< Scalar > sub_vecs[], const int num_targ_vecs, const SubVectorView< Scalar > sub_targ_vecs[], ReductTarget *const reduct_objs[]) |
Apply an RTOp in SMPD mode to a set of columns to a set of multi-vectors with contiguous storage per process. More... | |
RCP< FancyOStream > & | spmdApplyOpDumpOut () |
template<class Scalar > | |
void | print (const ConstSubVectorView< Scalar > &v, Teuchos::FancyOStream &out_arg) |
Variables | |
const int | NUM_ETRANS_ARGS = 3 |
const Teuchos::Tuple< char, NUM_ETRANS_ARGS > | transpMap = Teuchos::tuple('N', 'T', 'C') |
typedef Teuchos_Ordinal RTOpPack::Ordinal |
Definition at line 36 of file RTOpPack_Types.hpp.
typedef Teuchos_Ordinal RTOpPack::index_type |
Definition at line 57 of file RTOpPack_Types.hpp.
typedef char RTOpPack::char_type |
Definition at line 59 of file RTOpPack_Types.hpp.
enum RTOpPack::ETransp |
Enumerator | |
---|---|
NOTRANS | |
TRANS | |
CONJTRANS |
Definition at line 24 of file RTOpPack_LapackWrappers.hpp.
Enumerator | |
---|---|
REDUCT_TYPE_SUM | |
REDUCT_TYPE_MAX | |
REDUCT_TYPE_MIN |
Definition at line 211 of file RTOpPack_RTOpTHelpers_decl.hpp.
void RTOpPack::assign_entries | ( | const Ptr< const SubVectorView< Scalar > > & | msv, |
const ConstSubVectorView< Scalar > & | sv | ||
) |
Definition at line 256 of file RTOpPack_Types.hpp.
void RTOpPack::assign_entries | ( | const Ptr< const SubMultiVectorView< Scalar > > & | msmv, |
const ConstSubMultiVectorView< Scalar > & | smv | ||
) |
Definition at line 503 of file RTOpPack_Types.hpp.
void RTOpPack::getrf | ( | const SubMultiVectorView< Scalar > & | A, |
const ArrayView< int > & | ipiv, | ||
const Ptr< int > & | rank | ||
) |
Peform an in-place factorization of a square or rectangular matrix.
A | [in/out] On input, contains the entries of the square matrix. On output, contains the L and U factors. |
ipiv | [in] On output, contains the pivots used in the factorization. Note: This will be a 1-based valued array since this is a Fortran routine! |
rank | [out] On output, gives the rank of the factorization. |
Definition at line 70 of file RTOpPack_LapackWrappers.hpp.
void RTOpPack::getrs | ( | const ConstSubMultiVectorView< Scalar > & | A, |
const ArrayView< const int > & | ipiv, | ||
const ETransp | transp, | ||
const Ptr< const SubMultiVectorView< Scalar > > & | BX | ||
) |
Definition at line 100 of file RTOpPack_LapackWrappers.hpp.
RTOpPack::RTOP_ROP_1_REDUCT_SCALAR | ( | ROpCountNanInf | , |
index_type | , | ||
REDUCT_TYPE_SUM | |||
) |
Reduction operator that counts the number of entries that are NaN or Inf.
Definition at line 23 of file RTOpPack_ROpCountNanInf.hpp.
RTOpPack::RTOP_ROP_2_REDUCT_SCALAR | ( | ROpDotProd | , |
Scalar | , | ||
REDUCT_TYPE_SUM | |||
) |
Definition at line 20 of file RTOpPack_ROpDotProd.hpp.
RTOpPack::RTOP_ROP_1_REDUCT_SCALAR_CUSTOM_DEFAULT | ( | ROpMax | , |
Scalar | , | ||
REDUCT_TYPE_MAX | , | ||
Teuchos::as< Scalar > | -std::numeric_limits< Scalar >::max() | ||
) |
Maximum element: result = max{ v0[i], i=0...n-1 }
.
Definition at line 22 of file RTOpPack_ROpMax.hpp.
RTOpPack::RTOP_ROP_1_REDUCT_SCALAR_CUSTOM_DEFAULT | ( | ROpMin | , |
Scalar | , | ||
REDUCT_TYPE_MIN | , | ||
std::numeric_limits< Scalar >:: | max() | ||
) |
Minimum element: result = min{ v0[i], i=0...n-1 }
.
Definition at line 22 of file RTOpPack_ROpMin.hpp.
RTOpPack::RTOP_ROP_1_REDUCT_SCALAR | ( | ROpNorm1 | , |
typename ScalarTraits< Scalar >::magnitudeType | , | ||
REDUCT_TYPE_SUM | |||
) |
Class ROpNorm1.
Definition at line 21 of file RTOpPack_ROpNorm1.hpp.
RTOpPack::RTOP_ROP_1_REDUCT_SCALAR | ( | ROpNormInf | , |
typename ScalarTraits< Scalar >::magnitudeType | , | ||
REDUCT_TYPE_MAX | |||
) |
Norm Inf: result = max(|x[i]|, i=0...n-1).
Definition at line 34 of file RTOpPack_ROpNormInf.hpp.
RTOpPack::RTOP_ROP_1_REDUCT_SCALAR | ( | ROpSum | , |
Scalar | , | ||
REDUCT_TYPE_SUM | |||
) |
Class ROpSum: result = sum( v0[i], i=0...n-1 )
Definition at line 20 of file RTOpPack_ROpSum.hpp.
RTOpPack::RTOP_TOP_1_1 | ( | TOpAbs | ) |
Transformation operator that takes absolute values of elements: z0[i] = abs(v0[i]), i=0...n-1
.
Definition at line 23 of file RTOpPack_TOpAbs.hpp.
RTOpPack::RTOP_TOP_1_1 | ( | TOpAssignVectors | ) |
VectorBase assignment transformation operator: z0[i] = v0[i], i=0...n-1
.
Definition at line 22 of file RTOpPack_TOpAssignVectors.hpp.
RTOpPack::RTOP_TOP_1_1 | ( | TOpReciprocal | ) |
VectorBase assignment transformation operator: z0[i] = v0[i], i=0...n-1
.
Definition at line 23 of file RTOpPack_TOpReciprocal.hpp.
std::string RTOpPack::version | ( | ) |
Print the version of RTOp.
Definition at line 13 of file RTOpPack_version.cpp.
void RTOpPack::validate_apply_op | ( | const RTOpT< Scalar > & | op, |
const int | allowed_num_sub_vecs, | ||
const int | allowed_num_targ_sub_vecs, | ||
const bool | expect_reduct_obj, | ||
const ArrayView< const ConstSubVectorView< Scalar > > & | sub_vecs, | ||
const ArrayView< const SubVectorView< Scalar > > & | targ_sub_vecs, | ||
const Ptr< const ReductTarget > & | reduct_obj | ||
) |
Validate the input to an apply_op(...) function.
op | [in] The RTOpT object we are validating apply_op(...) input for. |
allowed_num_sub_vecs | [in] The allowed number of subvectors for sub_vecs.size(). If allowed_num_sub_vecs < 0 then this number is not valided. |
allowed_num_targ_sub_vecs | [in] The allowed number of subvectors for targ_sub_vecs.size(). If allowed_num_targ_sub_vecs < 0 then this number is not valided. |
expect_reduct_obj | [in] Determines if reduct_obj must be present or not and the type will be validated as well. |
sub_vecs | [in] Input to apply_op(...) being validated |
targ_sub_vecs | [in] Input to apply_op(...) being validated, not modified here |
reduct_obj_in | [in] Input to apply_op(...) being validated |
Definition at line 42 of file RTOpPack_RTOpTHelpers_def.hpp.
void RTOpPack::set_SPMD_apply_op_dump_out | ( | const RCP< FancyOStream > & | dumpOut | ) |
Set up to show a dump of RTOps applied through SPMD_apply_op().
dumOut | [in] RCP to output stream. If non-null, output will be dumped to this stream. If null, then no output will be dumped. |
Definition at line 22 of file RTOpPack_SPMD_apply_op.cpp.
int RTOpPack::serializedSize | ( | int | num_values, |
int | num_indexes, | ||
int | num_chars | ||
) |
Return the size in bytes of an external representation of a ReductTarget
object.
Definition at line 58 of file RTOpPack_SPMD_apply_op_def.hpp.
void RTOpPack::serialize | ( | const RTOpT< Scalar > & | op, |
Ordinal | num_values, | ||
Ordinal | num_indexes, | ||
Ordinal | num_chars, | ||
const ReductTarget & | reduct_obj, | ||
char | reduct_obj_ext[] | ||
) |
Serialize a ReductTarget
object.
Definition at line 72 of file RTOpPack_SPMD_apply_op_def.hpp.
void RTOpPack::deserialize | ( | const RTOpT< Scalar > & | op, |
int | num_values, | ||
int | num_indexes, | ||
int | num_chars, | ||
const char | reduct_obj_ext[], | ||
ReductTarget * | reduct_obj | ||
) |
Deserialize a ReductTarget
object.
Definition at line 111 of file RTOpPack_SPMD_apply_op_def.hpp.
void RTOpPack::SPMD_all_reduce | ( | const Teuchos::Comm< index_type > * | comm, |
const RTOpT< Scalar > & | op, | ||
const int | num_cols, | ||
const ReductTarget *const | i_reduct_objs[], | ||
ReductTarget *const | reduct_objs[] | ||
) |
Reduce a set of reduction objects.
ToDo: Finish documentation!
Definition at line 281 of file RTOpPack_SPMD_apply_op_def.hpp.
void RTOpPack::SPMD_apply_op | ( | const Teuchos::Comm< index_type > * | comm, |
const RTOpT< Scalar > & | op, | ||
const int | num_vecs, | ||
const ConstSubVectorView< Scalar > | sub_vecs[], | ||
const int | num_targ_vecs, | ||
const SubVectorView< Scalar > | targ_sub_vecs[], | ||
ReductTarget * | reduct_obj | ||
) |
Apply an RTOp in SMPD mode to a set of vectors with contiguous storage per process.
ToDo: Finish documentation!
Definition at line 314 of file RTOpPack_SPMD_apply_op_def.hpp.
void RTOpPack::SPMD_apply_op | ( | const Teuchos::Comm< index_type > * | comm, |
const RTOpT< Scalar > & | op, | ||
const int | num_cols, | ||
const int | num_multi_vecs, | ||
const ConstSubMultiVectorView< Scalar > | sub_multi_vecs[], | ||
const int | num_targ_multi_vecs, | ||
const SubMultiVectorView< Scalar > | targ_sub_multi_vecs[], | ||
RTOpPack::ReductTarget *const | reduct_objs[] | ||
) |
Apply an RTOp in SMPD mode to a set of columns to a set of multi-vectors with contiguous storage per process.
ToDo: Finish documentation!
Definition at line 334 of file RTOpPack_SPMD_apply_op_def.hpp.
void RTOpPack::SPMD_apply_op | ( | const Teuchos::Comm< index_type > * | comm, |
const RTOpT< Scalar > & | op, | ||
const int | num_cols, | ||
const int | num_vecs, | ||
const ConstSubVectorView< Scalar > | sub_vecs[], | ||
const int | num_targ_vecs, | ||
const SubVectorView< Scalar > | sub_targ_vecs[], | ||
ReductTarget *const | reduct_objs[] | ||
) |
Apply an RTOp in SMPD mode to a set of columns to a set of multi-vectors with contiguous storage per process.
ToDo: Finish documentation!
Definition at line 386 of file RTOpPack_SPMD_apply_op_def.hpp.
Teuchos::RCP< Teuchos::FancyOStream > & RTOpPack::spmdApplyOpDumpOut | ( | ) |
Definition at line 15 of file RTOpPack_SPMD_apply_op.cpp.
void RTOpPack::print | ( | const ConstSubVectorView< Scalar > & | v, |
Teuchos::FancyOStream & | out_arg | ||
) |
Definition at line 31 of file RTOpPack_SPMD_apply_op_def.hpp.
const int RTOpPack::NUM_ETRANS_ARGS = 3 |
Definition at line 27 of file RTOpPack_LapackWrappers.hpp.
const Teuchos::Tuple< char, RTOpPack::NUM_ETRANS_ARGS > RTOpPack::transpMap = Teuchos::tuple('N', 'T', 'C') |
Definition at line 14 of file RTOpPack_LapackWrappers.cpp.