AbstractLinAlgPack: C++ Interfaces For Vectors, Matrices And Related Linear Algebra Objects  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Types | List of all members
AbstractLinAlgPack::MultiVector Class Referenceabstract

Interface for a collection of non-mutable vectors (multi-vector, matrix). More...

#include <AbstractLinAlgPack_MultiVector.hpp>

Inheritance diagram for AbstractLinAlgPack::MultiVector:
Inheritance graph
[legend]

Public Types

enum  
 
typedef int access_by_t
 
typedef Teuchos::RCP< const
Vector
vec_ptr_t
 
typedef Teuchos::RCP< const
MultiVector
multi_vec_ptr_t
 
- Public Types inherited from AbstractLinAlgPack::MatrixOp
enum  EMatNormType { MAT_NORM_INF, MAT_NORM_2, MAT_NORM_1, MAT_NORM_FORB }
 Type of matrix norm. More...
 

Friends

void apply_op (EApplyBy apply_by, const RTOpPack::RTOp &primary_op, const size_t num_multi_vecs, const MultiVector *multi_vecs[], const size_t num_targ_multi_vecs, MultiVectorMutable *targ_multi_vecs[], RTOpPack::ReductTarget *reduct_objs[], const index_type primary_first_ele, const index_type primary_sub_dim, const index_type primary_global_offset, const index_type secondary_first_ele, const index_type secondary_sub_dim)
 
void apply_op (EApplyBy apply_by, const RTOpPack::RTOp &primary_op, const RTOpPack::RTOp &secondary_op, const size_t num_multi_vecs, const MultiVector *multi_vecs[], const size_t num_targ_multi_vecs, MultiVectorMutable *targ_multi_vecs[], RTOpPack::ReductTarget *reduct_obj, const index_type primary_first_ele, const index_type primary_sub_dim, const index_type primary_global_offset, const index_type secondary_first_ele, const index_type secondary_sub_dim)
 

Clone

virtual multi_vec_ptr_t mv_clone () const
 Clone the non-const multi-vector object. More...
 

Provide row, column and diagonal access as non-mutable vectors

virtual access_by_t access_by () const =0
 Return a bit field for the types of access that are the most convenient. More...
 
virtual vec_ptr_t col (index_type j) const =0
 Get a non-mutable column vector. More...
 
virtual vec_ptr_t row (index_type i) const =0
 Get a non-mutable row vector. More...
 
virtual vec_ptr_t diag (int k) const =0
 Get a non-mutable diagonal vector. More...
 

Sub-view methods

virtual multi_vec_ptr_t mv_sub_view (const Range1D &row_rng, const Range1D &col_rng) const
 Returns a sub-view of the multi vector. More...
 
multi_vec_ptr_t mv_sub_view (const index_type &rl, const index_type &ru, const index_type &cl, const index_type &cu) const
 Inlined implementation calls this->mv_sub_view(Range1D(rl,ru),Range1D(cl,cu)). More...
 

Collective apply_op() methods

virtual void apply_op (EApplyBy apply_by, const RTOpPack::RTOp &primary_op, const size_t num_multi_vecs, const MultiVector *multi_vecs[], const size_t num_targ_multi_vecs, MultiVectorMutable *targ_multi_vecs[], RTOpPack::ReductTarget *reduct_objs[], const index_type primary_first_ele, const index_type primary_sub_dim, const index_type primary_global_offset, const index_type secondary_first_ele, const index_type secondary_sub_dim) const
 Apply a reduction/transformation operator row by row, or column by column and return an array of the reduction objects. More...
 
virtual void apply_op (EApplyBy apply_by, const RTOpPack::RTOp &primary_op, const RTOpPack::RTOp &secondary_op, const size_t num_multi_vecs, const MultiVector *multi_vecs[], const size_t num_targ_multi_vecs, MultiVectorMutable *targ_multi_vecs[], RTOpPack::ReductTarget *reduct_obj, const index_type primary_first_ele, const index_type primary_sub_dim, const index_type primary_global_offset, const index_type secondary_first_ele, const index_type secondary_sub_dim) const
 Apply a reduction/transformation operator row by row, or column by column and reduce the intermediate reduction objects into one reduction object. More...
 

Overridden from MatrixOp

mat_ptr_t clone () const
 Returns this->mv_clone(). More...
 
mat_ptr_t sub_view (const Range1D &row_rng, const Range1D &col_rng) const
 Returns this->mv_sub_view(row_rng,col_rng) casted to a MatrixOp. More...
 
bool Mp_StMtM (MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &mwo_rhs1, BLAS_Cpp::Transp trans_rhs1, BLAS_Cpp::Transp trans_rhs2, value_type beta) const
 Provides a specialized implementation for mwo_rhs1 of type MatrixSymDiag. More...
 
bool Mp_StMtM (MatrixOp *mwo_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const MatrixOp &mwo_rhs2, BLAS_Cpp::Transp trans_rhs2, value_type beta) const
 Provides a specialized implementation for mwo_rhs2 of type MatrixSymDiag. More...
 

Additional Inherited Members

- Public Member Functions inherited from AbstractLinAlgPack::MatrixOp
virtual void zero_out ()
 M_lhs = 0 : Zero out the matrix. More...
 
virtual void Mt_S (value_type alpha)
 M_lhs *= alpha : Multiply a matrix by a scalar. More...
 
virtual MatrixOpoperator= (const MatrixOp &mwo_rhs)
 M_lhs = mwo_rhs : Virtual assignment operator. More...
 
virtual mat_mut_ptr_t clone ()
 Clone the non-const matrix object (if supported). More...
 
virtual std::ostream & output (std::ostream &out) const
 Virtual output function. More...
 
const MatNorm calc_norm (EMatNormType requested_norm_type=MAT_NORM_1, bool allow_replacement=false) const
 Compute a norm of this matrix. More...
 
mat_ptr_t sub_view (const index_type &rl, const index_type &ru, const index_type &cl, const index_type &cu) const
 Inlined implementation calls this->sub_view(Range1D(rl,ru),Range1D(cl,cu)). More...
 
virtual mat_ptr_t perm_view (const Permutation *P_row, const index_type row_part[], int num_row_part, const Permutation *P_col, const index_type col_part[], int num_col_part) const
 Create a permuted view: M_perm = P_row' * M * P_col. More...
 
virtual mat_ptr_t perm_view_update (const Permutation *P_row, const index_type row_part[], int num_row_part, const Permutation *P_col, const index_type col_part[], int num_col_part, const mat_ptr_t &perm_view) const
 Reinitialize a permuted view: M_perm = P_row' * M * P_col. More...
 
- Public Member Functions inherited from AbstractLinAlgPack::MatrixBase
virtual ~MatrixBase ()
 Virtual destructor. More...
 
virtual const VectorSpacespace_cols () const =0
 Vector space for vectors that are compatible with the columns of the matrix. More...
 
virtual const VectorSpacespace_rows () const =0
 Vector space for vectors that are compatible with the rows of the matrix. More...
 
virtual size_type rows () const
 Return the number of rows in the matrix. More...
 
virtual size_type cols () const
 Return the number of columns in the matrix. More...
 
virtual size_type nz () const
 Return the number of nonzero elements in the matrix. More...
 
- Protected Member Functions inherited from AbstractLinAlgPack::MatrixOp
virtual bool Mp_StM (MatrixOp *mwo_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs) const
 mwo_lhs += alpha * op(M_rhs) (BLAS xAXPY). More...
 
virtual bool Mp_StM (value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs)
 M_lhs += alpha * op(mwo_rhs) (BLAS xAXPY). More...
 
virtual bool Mp_StMtP (MatrixOp *mwo_lhs, value_type alpha, BLAS_Cpp::Transp M_trans, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans) const
 mwo_lhs += alpha * op(M_rhs) * op(P_rhs). More...
 
virtual bool Mp_StMtP (value_type alpha, const MatrixOp &mwo_rhs, BLAS_Cpp::Transp M_trans, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans)
 M_lhs += alpha * op(mwo_rhs) * op(P_rhs). More...
 
virtual bool Mp_StPtM (MatrixOp *mwo_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans, BLAS_Cpp::Transp M_trans) const
 mwo_lhs += alpha * op(P_rhs) * op(M_rhs). More...
 
virtual bool Mp_StPtM (value_type alpha, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans, const MatrixOp &mwo_rhs, BLAS_Cpp::Transp M_trans)
 M_lhs += alpha * op(P_rhs) * op(mwo_rhs). More...
 
virtual bool Mp_StPtMtP (MatrixOp *mwo_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, BLAS_Cpp::Transp M_trans, const GenPermMatrixSlice &P_rhs2, BLAS_Cpp::Transp P_rhs2_trans) const
 mwo_lhs += alpha * op(P_rhs1) * op(M_rhs) * op(P_rhs2). More...
 
virtual bool Mp_StPtMtP (value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, const MatrixOp &mwo_rhs, BLAS_Cpp::Transp M_trans, const GenPermMatrixSlice &P_rhs2, BLAS_Cpp::Transp P_rhs2_trans)
 M_lhs += alpha * op(P_rhs1) * op(mwo_rhs) * op(P_rhs2). More...
 
virtual void Vp_StMtV (VectorMutable *v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta) const =0
 v_lhs = alpha * op(M_rhs1) * v_rhs2 + beta * v_lhs (BLAS xGEMV) More...
 
virtual void Vp_StMtV (VectorMutable *v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice &sv_rhs2, value_type beta) const
 v_lhs = alpha * op(M_rhs1) * sv_rhs2 + beta * v_lhs (BLAS xGEMV) More...
 
virtual void Vp_StPtMtV (VectorMutable *v_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, BLAS_Cpp::Transp M_rhs2_trans, const Vector &v_rhs3, value_type beta) const
 v_lhs = alpha * op(P_rhs1) * op(M_rhs2) * v_rhs3 + beta * v_rhs More...
 
virtual void Vp_StPtMtV (VectorMutable *v_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, BLAS_Cpp::Transp M_rhs2_trans, const SpVectorSlice &sv_rhs3, value_type beta) const
 v_lhs = alpha * op(P_rhs1) * op(M_rhs2) * sv_rhs3 + beta * v_rhs More...
 
virtual value_type transVtMtV (const Vector &v_rhs1, BLAS_Cpp::Transp trans_rhs2, const Vector &v_rhs3) const
 result = v_rhs1' * op(M_rhs2) * v_rhs3 More...
 
virtual value_type transVtMtV (const SpVectorSlice &sv_rhs1, BLAS_Cpp::Transp trans_rhs2, const SpVectorSlice &sv_rhs3) const
 result = sv_rhs1' * op(M_rhs2) * sv_rhs3 More...
 
virtual void syr2k (BLAS_Cpp::Transp M_trans, value_type alpha, const GenPermMatrixSlice &P1, BLAS_Cpp::Transp P1_trans, const GenPermMatrixSlice &P2, BLAS_Cpp::Transp P2_trans, value_type beta, MatrixSymOp *symwo_lhs) const
 Perform a specialized rank-2k update of a dense symmetric matrix of the form: More...
 
virtual bool Mp_StMtM (value_type alpha, const MatrixOp &mwo_rhs1, BLAS_Cpp::Transp trans_rhs1, const MatrixOp &mwo_rhs2, BLAS_Cpp::Transp trans_rhs2, value_type beta)
 M_lhs = alpha * op(mwo_rhs1) * op(mwo_rhs2) + beta * mwo_lhs (left) (xGEMM) More...
 
virtual bool syrk (BLAS_Cpp::Transp M_trans, value_type alpha, value_type beta, MatrixSymOp *sym_lhs) const
 Perform a rank-k update of a symmetric matrix of the form: More...
 
virtual bool syrk (const MatrixOp &mwo_rhs, BLAS_Cpp::Transp M_trans, value_type alpha, value_type beta)
 Perform a rank-k update of a symmetric matrix of the form: More...
 

Detailed Description

Interface for a collection of non-mutable vectors (multi-vector, matrix).

This interface is quite restrictive in that it allows a client to access a matrix by accessing rows, columns and/or diagonals. The vector objects returned from the provided access methods row(), col() and diag() are abstract vectors so there is still good implementation flexibility but many MatrixOp implementations will not be able to support this interface.

The primary purpose for this interface is to allow for convienent aggregations of column vectors. Such an orderly arrangement allows for better optimized linear algebra operations such as matrix-matrix multiplication and the solution of linear systems for multiple right hand sides. Every application area (serial parallel, out-of-core etc.) should be able to define at least one reasonbly efficient implementation of a MultiVector (or a MultiVectorMutable) subclass.

The MultiVector interface is derived from the MatrixOp interface and therefore a MultiVector can be considered as a matrix which has some interesting implications. As an extended matrix interface, this is somewhat of a "last resort" interface that allows many matrix operations to have default implementations based on vector operations. None of the linear algebra methods in MatrixOp or any of the other matrix interfaces have methods that directly accept MultiVector objects. However, since MultiVector is derived from MatrixOp, a MultiVector object can be used anywere a MatrixOp object is accepted. In fact, many of the default implementations for the linear algebra methods in MatrixOp test to see if the matrix arguments support MultiVector (or MultiVectorMutable) and will fail if these interfaces are not supported.

Note that only certain kinds of access may be preferred and it is allowed for subclasses to return NULL vector objects for some types of access. For example, a matrix may be naturally oriented by column (for the primary role of as a multi-vector) but row or diagonal access may be very inefficient. For this reason, the client should call the access_by() method which returns a bit field that the client can compare to the constants ROW_ACCESS, COL_ACCESS and DIAG_ACCESS. The method access_by() only returns the types of access that are guarrentted to be efficient, but does not necessarily imply that a type of access is not supported. For example, (this->access_by() & ROW_ACCESS) == false this does not mean that this->row(i) will return NULL, it only means that row access will be inefficient. To determine if a certain type of access is even possible, check the return for row(), col() and/or diag(). For example, if this->rows(1) returns NULL, then this is a flag that row access for every row is not supported. Diagonal access may be a different story. For some matrix subclasses, only the center diagonal my be easily accessable in which case diag(0) may return != NULL but diag(k) for all k != 0 may return NULL.

Note that since, this interface is derived from MatrixOp that it must support the methods space_rows() and space_cols(). This does not imply however that either of the access methods row() or col() must return non-NULL.

Examples of matrix implementations that can support this interface are a dense BLAS compatible matrix (ROW_ACCESS, COL_ACCESS and DIAG_ACCESS), a compressed column sparse matrix (COL_ACCESS only), a compressed row sparse matrix (ROW_ACCESS only), a unordered sparse matrix with the diagonal explicity stored (diag(0) != NULL) etc.

Another very powerfull feature of this interface is the ability to apply reduction/transformation operators over a sub-set of rows and columns in a set of multi-vector objects. The behavior is identical as if the client extracted the rows or columns in a set of multi-vectors and called Vector::apply_op() itself. However, the advantage of using the multi-vector methods is that there may be greater opertunity to exploit parallelism. Also, the intermediate reduction objects over a set of rows or columns can be reduced by a secondary reduction object.

ToDo: Finish documentation!

Definition at line 170 of file AbstractLinAlgPack_MultiVector.hpp.

Member Typedef Documentation

Definition at line 179 of file AbstractLinAlgPack_MultiVector.hpp.

Definition at line 187 of file AbstractLinAlgPack_MultiVector.hpp.

Definition at line 189 of file AbstractLinAlgPack_MultiVector.hpp.

Member Enumeration Documentation

anonymous enum

Definition at line 181 of file AbstractLinAlgPack_MultiVector.hpp.

Member Function Documentation

MultiVector::multi_vec_ptr_t AbstractLinAlgPack::MultiVector::mv_clone ( ) const
virtual

Clone the non-const multi-vector object.

The default implementation returns return.get()==NULL.

Reimplemented in AbstractLinAlgPack::MultiVectorMutable.

Definition at line 162 of file AbstractLinAlgPack_MultiVector.cpp.

virtual access_by_t AbstractLinAlgPack::MultiVector::access_by ( ) const
pure virtual

Return a bit field for the types of access that are the most convenient.

Postconditions:

  • return & COL_ACCESS || return & ROW_ACCESS || return & DIAG_ACCESS

Implemented in AbstractLinAlgPack::MultiVectorMutableDense, and AbstractLinAlgPack::MultiVectorMutableCols.

virtual vec_ptr_t AbstractLinAlgPack::MultiVector::col ( index_type  j) const
pure virtual

Get a non-mutable column vector.

Postconditions:

  • [this->access_by() & COL_ACCESS] return.get() != NULL
  • [return.get() != NULL] space_cols().is_compatible(return->space()) == true

Implemented in AbstractLinAlgPack::MultiVectorMutable.

virtual vec_ptr_t AbstractLinAlgPack::MultiVector::row ( index_type  i) const
pure virtual

Get a non-mutable row vector.

Postconditions:

  • [this->access_by() & ROW_ACCESS] return.get() != NULL
  • [return.get() != NULL] space_rows().is_compatible(return->space()) == true

Implemented in AbstractLinAlgPack::MultiVectorMutable.

virtual vec_ptr_t AbstractLinAlgPack::MultiVector::diag ( int  k) const
pure virtual

Get a non-mutable diagonal vector.

Postconditions:

  • [this->access_by() & DIAG_ACCESS] return.get() != NULL

Implemented in AbstractLinAlgPack::MultiVectorMutable.

MultiVector::multi_vec_ptr_t AbstractLinAlgPack::MultiVector::mv_sub_view ( const Range1D row_rng,
const Range1D col_rng 
) const
virtual

Returns a sub-view of the multi vector.

ToDo: Finish documentation!

The default implementation returns a MultiVectorSubView object for any valid arbitary sub-view.

Reimplemented in AbstractLinAlgPack::MultiVectorMutable.

Definition at line 168 of file AbstractLinAlgPack_MultiVector.cpp.

MultiVector::multi_vec_ptr_t AbstractLinAlgPack::MultiVector::mv_sub_view ( const index_type &  rl,
const index_type &  ru,
const index_type &  cl,
const index_type &  cu 
) const
inline

Inlined implementation calls this->mv_sub_view(Range1D(rl,ru),Range1D(cl,cu)).

Definition at line 418 of file AbstractLinAlgPack_MultiVector.hpp.

void AbstractLinAlgPack::MultiVector::apply_op ( EApplyBy  apply_by,
const RTOpPack::RTOp &  primary_op,
const size_t  num_multi_vecs,
const MultiVector multi_vecs[],
const size_t  num_targ_multi_vecs,
MultiVectorMutable targ_multi_vecs[],
RTOpPack::ReductTarget reduct_objs[],
const index_type  primary_first_ele,
const index_type  primary_sub_dim,
const index_type  primary_global_offset,
const index_type  secondary_first_ele,
const index_type  secondary_sub_dim 
) const
protectedvirtual

Apply a reduction/transformation operator row by row, or column by column and return an array of the reduction objects.

Preconditions:

  • [apply_by == APPLY_BY_COL] (this->access_by() & COL_ACCESS) == true) (throw ???)
  • [apply_by == APPLY_BY_ROW] (this->access_by() & ROW_ACCESS) == true) (throw ???)
  • ToDo: Finish!

The default implementation calls this->apply_op().

ToDo: Finish documentation!

Definition at line 176 of file AbstractLinAlgPack_MultiVector.cpp.

void AbstractLinAlgPack::MultiVector::apply_op ( EApplyBy  apply_by,
const RTOpPack::RTOp &  primary_op,
const RTOpPack::RTOp &  secondary_op,
const size_t  num_multi_vecs,
const MultiVector multi_vecs[],
const size_t  num_targ_multi_vecs,
MultiVectorMutable targ_multi_vecs[],
RTOpPack::ReductTarget reduct_obj,
const index_type  primary_first_ele,
const index_type  primary_sub_dim,
const index_type  primary_global_offset,
const index_type  secondary_first_ele,
const index_type  secondary_sub_dim 
) const
protectedvirtual

Apply a reduction/transformation operator row by row, or column by column and reduce the intermediate reduction objects into one reduction object.

Preconditions:

  • [apply_by == APPLY_BY_COL] (this->access_by() & COL_ACCESS) == true) (throw ???)
  • [apply_by == APPLY_BY_ROW] (this->access_by() & ROW_ACCESS) == true) (throw ???)
  • ToDo: Finish!

The default implementation calls this->apply_op().

ToDo: Finish documentation!

Definition at line 235 of file AbstractLinAlgPack_MultiVector.cpp.

MatrixOp::mat_ptr_t AbstractLinAlgPack::MultiVector::clone ( ) const
virtual

Returns this->mv_clone().

Reimplemented from AbstractLinAlgPack::MatrixOp.

Reimplemented in AbstractLinAlgPack::MultiVectorMutableCols.

Definition at line 285 of file AbstractLinAlgPack_MultiVector.cpp.

MatrixOp::mat_ptr_t AbstractLinAlgPack::MultiVector::sub_view ( const Range1D row_rng,
const Range1D col_rng 
) const
virtual

Returns this->mv_sub_view(row_rng,col_rng) casted to a MatrixOp.

Reimplemented from AbstractLinAlgPack::MatrixOp.

Definition at line 291 of file AbstractLinAlgPack_MultiVector.cpp.

bool AbstractLinAlgPack::MultiVector::Mp_StMtM ( MatrixOp mwo_lhs,
value_type  alpha,
const MatrixOp mwo_rhs1,
BLAS_Cpp::Transp  trans_rhs1,
BLAS_Cpp::Transp  trans_rhs2,
value_type  beta 
) const
virtual

Provides a specialized implementation for mwo_rhs1 of type MatrixSymDiag.

Returns
Returns true and implements the operation if dynamic_cast<MatrixSymDiag>(&mwo_rhs1) != NULL && op(<em>this).access_by() =& MultiVector::COL_ACCESS && (mvm_lhs = dynamic_cast<MultiVectorMutable>(&mwo_lhs)) != NULL && mvm_lhs->access_by() & MultiVector::COL_ACCESS. Otherwise, this function returns false and does not implement the operation. or dynamic_cast<const MatrixSymDiag>(&mwo_rhs1) != NULL.

The default implementation relies on column access of op(*this) and mwo_lhs to implement this method.

Reimplemented from AbstractLinAlgPack::MatrixOp.

Reimplemented in AbstractLinAlgPack::MultiVectorMutableDense.

Definition at line 296 of file AbstractLinAlgPack_MultiVector.cpp.

bool AbstractLinAlgPack::MultiVector::Mp_StMtM ( MatrixOp mwo_lhs,
value_type  alpha,
BLAS_Cpp::Transp  trans_rhs1,
const MatrixOp mwo_rhs2,
BLAS_Cpp::Transp  trans_rhs2,
value_type  beta 
) const
virtual

Provides a specialized implementation for mwo_rhs2 of type MatrixSymDiag.

Returns
Returns true and implements the operation if dynamic_cast<MatrixSymDiag>(&mwo_rhs1) != NULL && op(<em>this).access_by() =& MultiVector::ROW_ACCESS && (mvm_lhs = dynamic_cast<MultiVectorMutable>(&mwo_lhs)) != NULL && mvm_lhs->access_by() & MultiVector::ROW_ACCESS. Otherwise, this function returns false and does not implement the operation. or dynamic_cast<const MatrixSymDiag>(&mwo_rhs1) != NULL.

The default implementation relies on row access of op(*this) and mwo_lhs to implement this method.

Reimplemented from AbstractLinAlgPack::MatrixOp.

Reimplemented in AbstractLinAlgPack::MultiVectorMutableDense.

Definition at line 311 of file AbstractLinAlgPack_MultiVector.cpp.

Friends And Related Function Documentation

void apply_op ( EApplyBy  apply_by,
const RTOpPack::RTOp &  primary_op,
const size_t  num_multi_vecs,
const MultiVector multi_vecs[],
const size_t  num_targ_multi_vecs,
MultiVectorMutable targ_multi_vecs[],
RTOpPack::ReductTarget reduct_objs[],
const index_type  primary_first_ele,
const index_type  primary_sub_dim,
const index_type  primary_global_offset,
const index_type  secondary_first_ele,
const index_type  secondary_sub_dim 
)
friend

ToDo: Finish documentation!

void apply_op ( EApplyBy  apply_by,
const RTOpPack::RTOp &  primary_op,
const RTOpPack::RTOp &  secondary_op,
const size_t  num_multi_vecs,
const MultiVector multi_vecs[],
const size_t  num_targ_multi_vecs,
MultiVectorMutable targ_multi_vecs[],
RTOpPack::ReductTarget reduct_obj,
const index_type  primary_first_ele,
const index_type  primary_sub_dim,
const index_type  primary_global_offset,
const index_type  secondary_first_ele,
const index_type  secondary_sub_dim 
)
friend

ToDo: Finish documentation!


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