MOOCHO (Single Doxygen Collection)
Version of the Day
|
Base class for all matrices that support basic matrix operations. More...
#include <AbstractLinAlgPack_MatrixOp.hpp>
Classes | |
class | IncompatibleMatrices |
Thrown if matrices are not compatible. More... | |
struct | MatNorm |
Returned form calc_norm() . More... | |
class | MethodNotImplemented |
Thrown if a method is not implemented. More... | |
Public types | |
enum | EMatNormType { MAT_NORM_INF, MAT_NORM_2, MAT_NORM_1, MAT_NORM_FORB } |
Type of matrix norm. More... | |
Minimal modifying methods | |
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 MatrixOp & | operator= (const MatrixOp &mwo_rhs) |
M_lhs = mwo_rhs : Virtual assignment operator. More... | |
Clone | |
virtual mat_mut_ptr_t | clone () |
Clone the non-const matrix object (if supported). More... | |
virtual mat_ptr_t | clone () const |
Clone the const matrix object (if supported). More... | |
Output | |
virtual std::ostream & | output (std::ostream &out) const |
Virtual output function. More... | |
Norms | |
const MatNorm | calc_norm (EMatNormType requested_norm_type=MAT_NORM_1, bool allow_replacement=false) const |
Compute a norm of this matrix. More... | |
Sub-matrix views | |
virtual mat_ptr_t | sub_view (const Range1D &row_rng, const Range1D &col_rng) const |
Create a transient constant sub-matrix view of this matrix (if supported). 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... | |
Permuted views | |
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... | |
Level-1 BLAS | |
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... | |
Level-2 BLAS | |
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... | |
Level-3 BLAS | |
virtual 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 |
mwo_lhs = alpha * op(M_rhs1) * op(mwo_rhs2) + beta * mwo_lhs (left) (xGEMM). More... | |
virtual 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 |
mwo_lhs = alpha * op(mwo_rhs1) * op(M_rhs2) + beta * mwo_lhs (right) (xGEMM) 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... | |
Additional Inherited Members | |
Public Member Functions inherited from AbstractLinAlgPack::MatrixBase | |
virtual | ~MatrixBase () |
Virtual destructor. More... | |
virtual const VectorSpace & | space_cols () const =0 |
Vector space for vectors that are compatible with the columns of the matrix. More... | |
virtual const VectorSpace & | space_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... | |
Base class for all matrices that support basic matrix operations.
These basic operations are:
Level-1 BLAS
mwo_lhs += alpha * op(M_rhs)
(BLAS xAXPY)
mwo_lhs += alpha * op(M_rhs) * op(P_rhs)
mwo_lhs += alpha * op(P_rhs) * op(M_rhs)
mwo_lhs += alpha * op(P1_rhs) * op(M_rhs) * op(P2_rhs)
M_lhs += alpha * op(mwo_rhs)
(BLAS xAXPY)
M_lhs += alpha * op(mwo_rhs) * op(P_rhs)
M_lhs += alpha * op(P_rhs) * op(mwo_rhs)
M_lhs += alpha * op(P1_rhs) * op(mwo_rhs) * op(P2_rhs)
Level-2 BLAS
v_lhs = alpha * op(M_rhs1) * v_rhs2 + beta * v_lhs
(BLAS xGEMV)
v_lhs = alpha * op(M_rhs1) * sv_rhs2 + beta * v_lhs
(BLAS xGEMV)
v_lhs = alpha * op(P_rhs1) * op(M_rhs2) * v_rhs3 + beta * v_lhs
v_lhs = alpha * op(P_rhs1) * op(M_rhs2) * sv_rhs3 + beta * v_lhs
result = v_rhs1' * op(M_rhs2) * v_rhs3
result = sv_rhs1' * op(M_rhs2) * sv_rhs3
Level-3 BLAS
mwo_lhs = alpha * op(M_rhs1) * op(mwo_rhs2) + beta * mwo_lhs
(right) (xGEMM)
mwo_lhs = alpha * op(mwo_rhs1) * op(M_rhs2) + beta * mwo_lhs
(left) (xGEMM)
M_lhs = alpha * op(mwo_rhs1) * op(mwo_rhs2) + beta * M_lhs
(xGEMM)
All of the Level-1, Level-2 and Level-3 BLAS operations have default implementations based on the Level-2 BLAS operation:
v_lhs = alpha * op(M_rhs1) * v_rhs2 + beta * v_lhs
(BLAS xGEMV)
The only methods that have to be overridden are space_cols()
, space_rows()
and the single Vp_StMtV()
method shown above. This is to allow fast prototyping of matrix subclasses and for postponement of writting specialized methods of other time critical operations until later if they are needed.
The vector space objects returned by the methods space_cols()
and space_rows()
are specifically bound to this matrix object. The vector space objects returned should only be considered to be transient and may become invalid if this
is modified in some significant way (but not through this MatrixOp
interface obviously).
Most of the Level-1 through Level-3 BLAS methods should not be called directly by clients, but instead be called through the provided non-member functions. The Level-1 and Level-3 matrix methods of this class have a special protocal in order to deal with the multiple dispatch problem. In essence, a poor man's multiple dispatch is used to allow each of the participating matrix objects a chance to implement an operation. In each case, a non-member function must be called by the client which calls the virtual methods on the matrix arguments one at a time. All of the Level-1 and Level-3 matrix methods are implemented for the case where the lhs matrix supports the MultiVectorMutable
interface. These matrix operations are then implemented in terms of AbstractLinAlgPack::Vp_StMtV(...)
which must be implemented for every matrix subclass. Therefore, any combination of rhs matrices can always be used in any matrix operation as long as a compatible (i.e. vector spaces match up) MultiVectorMutable
object is used as the lhs matrix argument.
Note, this behavior is only implemented by the nonmember functions AbstractLinAlgPack::Mp_StM(...)
or AbstractLinAlgPack::Mp_StMtM(...)
. All of the default virtual implementations of Mp_StM(...)
and Mp_StMtM(...)
return false.
This form of multiple dispatach is not ideal in the sense that the first matrix argument that can implement the method will do so instead of perhaps the best implementation that could be provided by another matrix argument. Therefore, a matrix subclass should only override one of these matrix methods if it can provide a significantly better implementation than the default. If a client needs exact control of the implementation of a matrix operation, then they should consider using a ``Strategy'' object.
ToDo: Add more detailed documentation for the default Level-1 and Level-3 BLAS methods.
Definition at line 129 of file AbstractLinAlgPack_MatrixOp.hpp.
Type of matrix norm.
Definition at line 244 of file AbstractLinAlgPack_MatrixOp.hpp.
|
virtual |
M_lhs = 0 : Zero out the matrix.
The default implementation throws an exception. This is not the best design but it meets some needs. Any matrix implementation could implement this method and mimic the behavior (i.e. see the matrix subclass MatrixZero
). However, only matrices that are going to be on the lhs (non-const) of a Level-1 or Level-3 BLAS operation need every implement this method.
Reimplemented in AbstractLinAlgPack::MultiVectorMutableDense, AbstractLinAlgPack::MultiVectorMutable, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MultiVectorMutableCols, and AbstractLinAlgPack::MatrixZero.
Definition at line 64 of file AbstractLinAlgPack_MatrixOp.cpp.
|
virtual |
M_lhs *= alpha : Multiply a matrix by a scalar.
The default implementation throws an exception. This is not the best design but it meets some needs. Any matrix implementation could implement this method and mimic the behavior (i.e. simply implement the matrix M
as (alpha * M
). This method is only called in a few specialized situations.
Reimplemented in AbstractLinAlgPack::MultiVectorMutableDense, AbstractLinAlgPack::MultiVectorMutable, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MultiVectorMutableCols, and AbstractLinAlgPack::MatrixZero.
Definition at line 72 of file AbstractLinAlgPack_MatrixOp.cpp.
M_lhs = mwo_rhs : Virtual assignment operator.
The default implementation just throws a std::logic_error exception if it is not assignment to self. A more specialized implementation could use this to copy the state to this
matrix object from a compatible M
matrix object.
Reimplemented in AbstractLinAlgPack::MatrixSparseCOORSerial, AbstractLinAlgPack::MultiVectorMutableDense, AbstractLinAlgPack::MultiVectorMutable, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MultiVectorMutableCols, AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixSymDiagStd, AbstractLinAlgPack::MatrixOpThyra, AbstractLinAlgPack::COOMatrixPartitionViewSubclass, AbstractLinAlgPack::MultiVectorMutableThyra, AbstractLinAlgPack::MatrixZero, AbstractLinAlgPack::MatrixWithOpConcreteEncap< M >, and AbstractLinAlgPack::MatrixSymDiagSparseStd.
Definition at line 80 of file AbstractLinAlgPack_MatrixOp.cpp.
|
virtual |
Clone the non-const matrix object (if supported).
The primary purpose for this method is to allow a client to capture the current state of a matrix object and be guaranteed that some other client will not alter its behavior. A smart implementation will use reference counting and lazy evaluation internally and will not actually copy any large amount of data unless it has to.
The default implementation returns NULL which is perfectly acceptable. A matrix object is not required to return a non-NULL value but almost every good matrix implementation will.
Reimplemented in AbstractLinAlgPack::MatrixSymOp, AbstractLinAlgPack::MultiVectorMutable, AbstractLinAlgPack::MultiVectorMutableCols, AbstractLinAlgPack::MatrixOpNonsing, AbstractLinAlgPack::MatrixOpThyra, AbstractLinAlgPack::MatrixOpNonsingThyra, AbstractLinAlgPack::MultiVectorMutableThyra, AbstractLinAlgPack::MatrixSymOpNonsing, and AbstractLinAlgPack::MatrixSymOpNonsingSerial.
Definition at line 109 of file AbstractLinAlgPack_MatrixOp.cpp.
|
virtual |
Clone the const matrix object (if supported).
The behavior of this method is the same as for the non-const version above except it returns a smart pointer to a const matrix object.
Reimplemented in AbstractLinAlgPack::MultiVector, AbstractLinAlgPack::MatrixSymOp, AbstractLinAlgPack::MultiVectorMutableCols, AbstractLinAlgPack::MatrixOpNonsing, AbstractLinAlgPack::MatrixSymOpNonsing, and AbstractLinAlgPack::MatrixSymOpNonsingSerial.
Definition at line 115 of file AbstractLinAlgPack_MatrixOp.cpp.
|
virtual |
Virtual output function.
The default implementaion just extracts rows one at a time by calling this->Vp_StMtV()
with EtaVector
objects and then prints the rows.
Reimplemented in ConstrainedOptPack::MatrixSymPosDefLBFGS, AbstractLinAlgPack::MatrixSparseCOORSerial, AbstractLinAlgPack::MatrixOpSerial, ConstrainedOptPack::MatrixSymPosDefLBFGS, ConstrainedOptPack::MatrixVarReductImplicit, AbstractLinAlgPack::MultiVectorMutableDense, AbstractLinAlgPack::MatrixOpSubView, ConstrainedOptPack::MatrixDecompRangeOrthog, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixZero, ConstrainedOptPack::MatrixIdentConcat, AbstractLinAlgPack::MatrixSymIdent, ConstrainedOptPack::MatrixSymIdentitySerial, and AbstractLinAlgPack::MatrixSymDiagSparse.
Definition at line 91 of file AbstractLinAlgPack_MatrixOp.cpp.
const MatrixOp::MatNorm MatrixOp::calc_norm | ( | EMatNormType | requested_norm_type = MAT_NORM_1 , |
bool | allow_replacement = false |
||
) | const |
Compute a norm of this matrix.
requested_norm_type | [in] Determines the requested type of norm to be computed. |
allow_replacement | [in] Determines if the requested norm in specified in norm_type can be replaced with another norm that can be computed by the matrix. |
return.value
gives the value of the norm of type return.type
.Postconditions:
allow_replacement==true
, the matrix object must return a computted norm who's type is given in return.type
. allow_replacement==false
and the underlying matrix object can not compute the norm requested in norm_type
, then a MethodNotImplemented
exception will be thrown. If the matrix object can compute this norm, then return.type
will be equal to requested_norm_type
. The default implementation of this method uses Algorithm 2.5 in "Applied Numerical Linear Algebra" by James Demmel (1997) to estimate ||M||1 or ||M||inf. The algorithm uses some of the refinements in the referenced algorithm by Highman. This algorithm only requires mat-vecs and transposed mat-vecs so every matrix object can implement this method. The main purpose of this default implementation is to allow a default implementation of the estimation of the ||.||1
or ||.||inf
normed condition number in the class MatrixOpNonsing
. The default arguments for this function will compute a norm and will not thrown an exception. The default implementation will throw an exception for any other norm type than requested_norm_type == MAT_NORM_1
or requested_norm_type = MAT_NORM_INF
.
Definition at line 123 of file AbstractLinAlgPack_MatrixOp.cpp.
|
virtual |
Create a transient constant sub-matrix view of this matrix (if supported).
This view is to be used immediatly and then discarded.
This method can only be expected to return return.get() != NULL
if this->space_cols().sub_space(row_rng) != NULL
and this->space_rows().sub_space(col_rng) != NULL
.
It is allows for a matrix implementation to return return.get() == NULL
for any arbitrary subview.
The default implementation uses the matrix subclass MatrixOpSubView
and therefore, can return any arbitrary subview. More specialized implementations may want to restrict the subview that can be created somewhat.
Reimplemented in AbstractLinAlgPack::MatrixComposite, AbstractLinAlgPack::MultiVector, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, and AbstractLinAlgPack::MatrixOpNonsingAggr.
Definition at line 180 of file AbstractLinAlgPack_MatrixOp.cpp.
|
inline |
Inlined implementation calls this->sub_view(Range1D(rl,ru),Range1D(cl,cu))
.
Definition at line 975 of file AbstractLinAlgPack_MatrixOp.hpp.
|
virtual |
Create a permuted view: M_perm = P_row' * M * P_col
.
P_row | [in] Row permutation. If P_row == NULL then the indentity permutation is used. |
row_part | [in] Array (length num_row_part+1 ) storing the row indexes that may be passed to return->sub_view(r1,r2,...) . If row_part == NULL then the assumed array is { 1, this->rows() } . |
num_row_part | [in] Length of the array row_part . If row_part == NULL then this argument is ignored. |
P_col | [in] Column permutation. If P_col == NULL then the indentity permutation is used. |
col_part | [in] Array (length num_col_part+1 ) storing the column indexes that may be passed to return->sub_view(...,c1,c2) . If col_part == NULL then the assumed array is { 1, this->cols() } . |
num_col_part | [in] Length of the array col_part . If col_part == NULL then this argument is ignored. |
Preconditions:
P_row != NULL
] P_row->space().is_compatible(this->space_cols())
(throw VectorSpace::IncompatibleVectorSpaces
) P_col != NULL
] P_col->space().is_compatible(this->space_rows())
(throw VectorSpace::IncompatibleVectorSpaces
) row_part != NULL
] 1 <= row_part[i-1] < row_part[i] <= this->rows(), for i = 1...num_row_part
(throw ???) col_part != NULL
] 1 <= col_part[i-1] < col_part[i] <= this->cols(), for i = 1...num_col_part
(throw ???) Postconditions:
return->sub_view(R,C)
should be able to be created efficiently where R = [row_part[kr-1],row_part[kr]-1], for kr = 1...num_row_part
and C = [col_part[kc-1],col_part[kc]-1], for kc = 1...num_col_part
. The default implementation returns a MatrixPermAggr
object.
Definition at line 203 of file AbstractLinAlgPack_MatrixOp.cpp.
|
virtual |
Reinitialize a permuted view: M_perm = P_row' * M * P_col
.
P_row | [in] Same as input to perm_view() . |
row_part | [in] Same as input to perm_view() . |
num_row_part | [in] Same as input to perm_view() . |
P_col | [in] Same as input to perm_view() . |
col_part | [in] Same as input to perm_view() . |
num_col_part | [in] Same as input to perm_view() . |
perm_view | [in] Smart pointer to a permuted view returned from this->perm_view() . |
Preconditions:
P_row != NULL
] P_row->space().is_compatible(this->space_cols())
(throw VectorSpace::IncompatibleVectorSpaces
) P_col != NULL
] P_col->space().is_compatible(this->space_rows())
(throw VectorSpace::IncompatibleVectorSpaces
) row_part != NULL
] 1 <= row_part[i-1] < row_part[i] <= this->rows(), for i = 1...num_row_part
(throw ???) col_part != NULL
] 1 <= col_part[i-1] < col_part[i] <= this->cols(), for i = 1...num_col_part
(throw ???) The default implementation simply returns this->perm_view()
Definition at line 223 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
mwo_lhs += alpha * op(M_rhs) (BLAS xAXPY).
The default implementation does nothing returns false.
A client can not call this method call this method directly. Instead, use AbstractLinAlgPack::Mp_StM()
.
Reimplemented in AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MultiVectorMutableDense, AbstractLinAlgPack::MultiVectorMutable, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixSymDiagStd, and AbstractLinAlgPack::MatrixZero.
Definition at line 240 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
M_lhs += alpha * op(mwo_rhs) (BLAS xAXPY).
The default implementation does nothing and returns false.
A client can not call this method call this method directly. Instead, use AbstractLinAlgPack::Mp_StM()
.
Reimplemented in AbstractLinAlgPack::MultiVectorMutableDense, AbstractLinAlgPack::MatrixOpSubView, and AbstractLinAlgPack::MultiVectorMutable.
Definition at line 247 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
mwo_lhs += alpha * op(M_rhs) * op(P_rhs).
The default implementation does nothing and returns false.
A client can not call this method call this method directly. Instead, use AbstractLinAlgPack::Mp_StMtP()
.
Reimplemented in AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixOpNonsingAggr, and AbstractLinAlgPack::MatrixZero.
Definition at line 253 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
M_lhs += alpha * op(mwo_rhs) * op(P_rhs).
The default implementation does nothing and returns false.
A client can not call this method call this method directly. Instead, use AbstractLinAlgPack::Mp_StMtP()
.
Reimplemented in AbstractLinAlgPack::MatrixOpSubView.
Definition at line 262 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
mwo_lhs += alpha * op(P_rhs) * op(M_rhs).
The default implementation does nothing and returns false.
A client can not call this method call this method directly. Instead, use AbstractLinAlgPack::Mp_StPtM()
.
Reimplemented in AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixOpNonsingAggr, and AbstractLinAlgPack::MatrixZero.
Definition at line 271 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
M_lhs += alpha * op(P_rhs) * op(mwo_rhs).
The default implementation does nothing and returns false.
A client can not call this method call this method directly. Instead, use AbstractLinAlgPack::Mp_StPtM()
.
Reimplemented in AbstractLinAlgPack::MatrixOpSubView.
Definition at line 280 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
mwo_lhs += alpha * op(P_rhs1) * op(M_rhs) * op(P_rhs2).
The default implementation does nothing and returns false.
A client can not call this method call this method directly. Instead, use AbstractLinAlgPack::Mp_StPtMtP()
.
Reimplemented in AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixOpNonsingAggr, and AbstractLinAlgPack::MatrixZero.
Definition at line 289 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
M_lhs += alpha * op(P_rhs1) * op(mwo_rhs) * op(P_rhs2).
The default implementation does nothing and returns false.
A client can not call this method call this method directly. Instead, use AbstractLinAlgPack::Mp_StPtMtP()
.
Reimplemented in AbstractLinAlgPack::MatrixOpSubView.
Definition at line 299 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedpure virtual |
v_lhs = alpha * op(M_rhs1) * v_rhs2 + beta * v_lhs (BLAS xGEMV)
Implemented in AbstractLinAlgPack::MatrixComposite, AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MatrixSparseCOORSerial, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MultiVectorMutableCols, AbstractLinAlgPack::MatrixSymDiagStd, AbstractLinAlgPack::MatrixZero, AbstractLinAlgPack::MatrixOpThyra, AbstractLinAlgPack::MultiVectorMutableThyra, and AbstractLinAlgPack::MatrixSymIdent.
|
protectedvirtual |
v_lhs = alpha * op(M_rhs1) * sv_rhs2 + beta * v_lhs (BLAS xGEMV)
Reimplemented in AbstractLinAlgPack::MatrixComposite, AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MultiVectorMutableCols, AbstractLinAlgPack::MatrixSymDiagStd, and AbstractLinAlgPack::MatrixZero.
Definition at line 311 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
v_lhs = alpha * op(P_rhs1) * op(M_rhs2) * v_rhs3 + beta * v_rhs
Reimplemented in AbstractLinAlgPack::MatrixComposite, AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixOpNonsingAggr, and AbstractLinAlgPack::MatrixZero.
Definition at line 330 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
v_lhs = alpha * op(P_rhs1) * op(M_rhs2) * sv_rhs3 + beta * v_rhs
Reimplemented in AbstractLinAlgPack::MatrixComposite, AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixOpNonsingAggr, and AbstractLinAlgPack::MatrixZero.
Definition at line 345 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
result = v_rhs1' * op(M_rhs2) * v_rhs3
Reimplemented in AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixOpNonsingAggr, and AbstractLinAlgPack::MatrixZero.
Definition at line 360 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
result = sv_rhs1' * op(M_rhs2) * sv_rhs3
Reimplemented in AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixZero, and AbstractLinAlgPack::COOMatrixPartitionViewSubclass.
Definition at line 368 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
Perform a specialized rank-2k update of a dense symmetric matrix of the form:
symwo_lhs += alpha*op(P1')*op(M)*op(P2) + alpha*op(P2')*op(M')*op(P1) + beta*symwo_lhs
The reason that this operation is being classified as a level-2 operation is that the total flops should be of O(n^2)
and not O(n^2*k)
.
The default implementation is based on Mp_StMtP(...)
and Mp_StPtM(...)
. Of course in situations where this default implemention is inefficient the subclass should override this method.
Reimplemented in AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixOpNonsingAggr, and AbstractLinAlgPack::MatrixZero.
Definition at line 376 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
mwo_lhs = alpha * op(M_rhs1) * op(mwo_rhs2) + beta * mwo_lhs (left) (xGEMM).
The default implementation does nothing and returns false.
Warning! A client should never call this method call this method directly. Instead, use AbstractLinAlgPack::Mp_StMtM()
.
Reimplemented in AbstractLinAlgPack::MultiVector, AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MultiVectorMutableDense, AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixZero, AbstractLinAlgPack::MatrixOpThyra, and AbstractLinAlgPack::MultiVectorMutableThyra.
Definition at line 387 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
mwo_lhs = alpha * op(mwo_rhs1) * op(M_rhs2) + beta * mwo_lhs (right) (xGEMM)
The default implementation does nothing and returns false.
Warning! A client should never call this method call this method directly. Instead, use AbstractLinAlgPack::Mp_StMtM()
.
Reimplemented in AbstractLinAlgPack::MultiVector, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MultiVectorMutableDense, AbstractLinAlgPack::MatrixOpNonsingAggr, and AbstractLinAlgPack::MatrixZero.
Definition at line 395 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
M_lhs = alpha * op(mwo_rhs1) * op(mwo_rhs2) + beta * mwo_lhs (left) (xGEMM)
The default implementation does nothing and returns false.
Warning! A client should never call this method call this method directly. Instead, use AbstractLinAlgPack::Mp_StMtM()
.
Reimplemented in AbstractLinAlgPack::MatrixOpSubView.
Definition at line 403 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
Perform a rank-k update of a symmetric matrix of the form:
symwo_lhs += alpha*op(M)*op(M') + beta*symwo_lhs
where this
is the rhs matrix argument.
Never call this method directly. Instead use the nonmember function AbstractLinAlgPack::syrk()
.
The default implementation returns false
and does nothing.
Reimplemented in AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MultiVectorMutableDense, AbstractLinAlgPack::MatrixZero, AbstractLinAlgPack::MultiVectorMutableCols, and AbstractLinAlgPack::MatrixSymDiagStd.
Definition at line 412 of file AbstractLinAlgPack_MatrixOp.cpp.
|
protectedvirtual |
Perform a rank-k update of a symmetric matrix of the form:
M += alpha*op(mwo_rhs)*op(mwo_rhs') + beta*M
where this
is the lhs matrix argument.
Never call this method directly. Instead use the nonmember function AbstractLinAlgPack::syrk()
.
The default implementation returns false
and does nothing.
Definition at line 422 of file AbstractLinAlgPack_MatrixOp.cpp.
|
friend |
mwo_lhs *= alpha.
If alpha == 0.0
then mwo_lhs->zero_out()
will be called, otherwise mwo_lhs->Mt_S(alpha)
will be called. If alpha == 1.0
then nothing is done.
|
friend |
mwo_lhs += alpha * op(M_rhs) (BLAS xAXPY).
Entry point for (poor man's) multiple dispatch.
This method first calls M_rhs->Mp_StM(mwo_lhs,alpha,trans_rhs)
to give the rhs argument a chance to implement the operation. If M_rhs->Mp_StM(...)
returns false, then mwo_lhs->Mp_StM(alpha,*this,trans_rhs)
is called to give the lhs matrix argument a chance to implement the method. If mwo_lhs->Mp_StM(...)
returns false, then an attempt to perform a dynamic cast the lhs matrix argument to MultiVectorMutable
is attempted. If this cast failes, then an exception is thrown.
|
friend |
Entry point for (poor man's) multiple dispatch.
ToDo: Finish documentation!
|
friend |
Entry point for (poor man's) multiple dispatch.
ToDo: Finish documentation!
|
friend |
Entry point for (poor man's) multiple dispatch.
ToDo: Finish documentation!
|
friend |
Definition at line 857 of file AbstractLinAlgPack_MatrixOp.hpp.
|
friend |
Definition at line 866 of file AbstractLinAlgPack_MatrixOp.hpp.
|
friend |
Definition at line 875 of file AbstractLinAlgPack_MatrixOp.hpp.
|
friend |
Definition at line 886 of file AbstractLinAlgPack_MatrixOp.hpp.
|
friend |
Definition at line 897 of file AbstractLinAlgPack_MatrixOp.hpp.
|
friend |
Definition at line 906 of file AbstractLinAlgPack_MatrixOp.hpp.
|
friend |
Definition at line 915 of file AbstractLinAlgPack_MatrixOp.hpp.
|
friend |
This method first calls mwo_rhs1.Mp_StMtM(...)
to perform the opeation. If mwo_rhs1.Mp_StMtM(...)
returns false, then mwo_rhs2.Mp_StMtM(...)
is called. If mwo_rhs2.Mp_StMtM(...)
returns false, then mwo_lhs.Mp_StMtM(...)
is called.
As a last resort, the function attempts to cast dynamic_cast<MultiVectorMutable*>(mwo_lhs)
. If this dynamic cast fails, the this function throws an exception. Otherwise, the operation is implemented in terms of Vp_StMtV()
.
|
friend |
symwo_lhs += alpha*op(mwo_rhs)*op(mwo_rhs') + beta*symwo_lhs
The default implementation returns false
and does nothing.