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::MultiVectorMutable Class Referenceabstract

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

#include <AbstractLinAlgPack_MultiVectorMutable.hpp>

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

Public Types

typedef Teuchos::RCP
< VectorMutable
vec_mut_ptr_t
 
typedef Teuchos::RCP
< MultiVectorMutable
multi_vec_mut_ptr_t
 
- Public Types inherited from AbstractLinAlgPack::MultiVector
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...
 

Clone

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

Provide mutable row, column and/or diagonal access

virtual vec_mut_ptr_t col (index_type j)=0
 Get a mutable column vector. More...
 
virtual vec_mut_ptr_t row (index_type i)=0
 Get a mutable row vector. More...
 
virtual vec_mut_ptr_t diag (int k)=0
 Get a mutable diagonal vector. More...
 

Sub-view methods

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

Overridden from MatrixOp

mat_mut_ptr_t clone ()
 
void zero_out ()
 
void Mt_S (value_type alpha)
 
MatrixOpoperator= (const MatrixOp &mwo_rhs)
 
bool Mp_StM (MatrixOp *mwo_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs) const
 
bool Mp_StM (value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs)
 

Overridden from MultiVector

multi_vec_ptr_t mv_clone () const
 
vec_ptr_t col (index_type j) const
 
vec_ptr_t row (index_type i) const
 
vec_ptr_t diag (int k) const
 
multi_vec_ptr_t mv_sub_view (const Range1D &row_rng, const Range1D &col_rng) const
 

Additional Inherited Members

- Public Member Functions inherited from AbstractLinAlgPack::MultiVector
virtual access_by_t access_by () const =0
 Return a bit field for the types of access that are the most convenient. 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...
 
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...
 
- Public Member Functions inherited from AbstractLinAlgPack::MatrixOp
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::MultiVector
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...
 
- Protected Member Functions inherited from AbstractLinAlgPack::MatrixOp
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 mutable vectors (multi-vector, matrix).

This interface extends the MutiVector interface an allows mutable access to the constituent vectors.

These vectors allow the modification of the matrix row by row, column by column, and/or diagonal by diagonal. Each of the views is transient and should be used and discarded quickly.

Note that the underlying matrix is only guaranteed to be modified after the smart reference counted pointer returned from these methods is destoryed. For example, consider the following code:

void f( MultiVectorMutable* M, index_type i )
{
row_i =M->row(i);
row_i = 0.0;
// The underlying matrix may not be modified at this point.
row_i = NULL;
// Now the underlying matrix is guaranteed to be modified and
// we can assume this in the following code.
...
}

Default implementations of the const access methods row() col() and diag() from MultiVector call the non-const methods defined here and cast the pointers.

Many of the default implementations of the linear algebra operations in MatrixOp and the other matrix interfaces rely on the left hand side matrix objects supporting the MultiVectorMutable interface.

Definition at line 82 of file AbstractLinAlgPack_MultiVectorMutable.hpp.

Member Typedef Documentation

Definition at line 93 of file AbstractLinAlgPack_MultiVectorMutable.hpp.

Definition at line 95 of file AbstractLinAlgPack_MultiVectorMutable.hpp.

Member Function Documentation

MultiVectorMutable::multi_vec_mut_ptr_t AbstractLinAlgPack::MultiVectorMutable::mv_clone ( )
virtual

Clone the non-const multi-vector object.

The default implementation creates a new multi-vector and then copies the values.

Definition at line 89 of file AbstractLinAlgPack_MultiVectorMutable.cpp.

virtual vec_mut_ptr_t AbstractLinAlgPack::MultiVectorMutable::col ( index_type  j)
pure virtual

Get a mutable column vector.

Postconditions:

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

ToDo: Finish documentation!

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

virtual vec_mut_ptr_t AbstractLinAlgPack::MultiVectorMutable::row ( index_type  i)
pure virtual

Get a mutable row vector.

Postconditions:

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

ToDo: Finish documentation!

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

virtual vec_mut_ptr_t AbstractLinAlgPack::MultiVectorMutable::diag ( int  k)
pure virtual

Get a mutable diagonal vector.

Postconditions:

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

ToDo: Finish documentation!

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

MultiVectorMutable::multi_vec_mut_ptr_t AbstractLinAlgPack::MultiVectorMutable::mv_sub_view ( const Range1D row_rng,
const Range1D col_rng 
)
virtual

Returns a mutable sub-view of the multi vector.

ToDo: Finish documentation!

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

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

Definition at line 102 of file AbstractLinAlgPack_MultiVectorMutable.cpp.

MultiVectorMutable::multi_vec_mut_ptr_t AbstractLinAlgPack::MultiVectorMutable::mv_sub_view ( const index_type &  rl,
const index_type &  ru,
const index_type &  cl,
const index_type &  cu 
)
inline

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

Definition at line 211 of file AbstractLinAlgPack_MultiVectorMutable.hpp.

MatrixOp::mat_mut_ptr_t AbstractLinAlgPack::MultiVectorMutable::clone ( )
virtual
void AbstractLinAlgPack::MultiVectorMutable::zero_out ( )
virtual
void AbstractLinAlgPack::MultiVectorMutable::Mt_S ( value_type  alpha)
virtual
MatrixOp & AbstractLinAlgPack::MultiVectorMutable::operator= ( const MatrixOp mwo_rhs)
virtual
bool AbstractLinAlgPack::MultiVectorMutable::Mp_StM ( MatrixOp mwo_lhs,
value_type  alpha,
BLAS_Cpp::Transp  trans_rhs 
) const
virtual
bool AbstractLinAlgPack::MultiVectorMutable::Mp_StM ( value_type  alpha,
const MatrixOp M_rhs,
BLAS_Cpp::Transp  trans_rhs 
)
virtual
MultiVector::multi_vec_ptr_t AbstractLinAlgPack::MultiVectorMutable::mv_clone ( ) const
virtual

Reimplemented from AbstractLinAlgPack::MultiVector.

Definition at line 170 of file AbstractLinAlgPack_MultiVectorMutable.cpp.

MultiVectorMutable::vec_ptr_t AbstractLinAlgPack::MultiVectorMutable::col ( index_type  j) const
virtual
MultiVectorMutable::vec_ptr_t AbstractLinAlgPack::MultiVectorMutable::row ( index_type  i) const
virtual
MultiVectorMutable::vec_ptr_t AbstractLinAlgPack::MultiVectorMutable::diag ( int  k) const
virtual
MultiVectorMutable::multi_vec_ptr_t AbstractLinAlgPack::MultiVectorMutable::mv_sub_view ( const Range1D row_rng,
const Range1D col_rng 
) const
virtual

Reimplemented from AbstractLinAlgPack::MultiVector.

Definition at line 191 of file AbstractLinAlgPack_MultiVectorMutable.cpp.


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