#include <AbstractLinAlgPack_MatrixOp.hpp>
Inheritance diagram for AbstractLinAlgPack::MatrixOp:
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. | |
| virtual void | Mt_S (value_type alpha) |
| M_lhs *= alpha : Multiply a matrix by a scalar. | |
| virtual MatrixOp & | operator= (const MatrixOp &mwo_rhs) |
| M_lhs = mwo_rhs : Virtual assignment operator. | |
Clone | |
| virtual mat_mut_ptr_t | clone () |
| Clone the non-const matrix object (if supported). | |
| virtual mat_ptr_t | clone () const |
| Clone the const matrix object (if supported). | |
Output | |
| virtual std::ostream & | output (std::ostream &out) const |
| Virtual output function. | |
Norms | |
| const MatNorm | calc_norm (EMatNormType requested_norm_type=MAT_NORM_1, bool allow_replacement=false) const |
| Compute a norm of this matrix. | |
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). | |
| 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)). | |
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. | |
| 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. | |
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). | |
| virtual bool | Mp_StM (value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs) |
| M_lhs += alpha * op(mwo_rhs) (BLAS xAXPY). | |
| 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). | |
| 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). | |
| 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). | |
| 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). | |
| 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). | |
| 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). | |
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) | |
| 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) | |
| 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 | |
| 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 | |
| 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 | |
| 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 | |
| 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:. | |
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). | |
| 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) | |
| 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). | |
| 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:. | |
| 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:. | |
Friends | |
| void | Mt_S (MatrixOp *mwo_lhs, value_type alpha) |
| | |
| void | Mp_StM (MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs) |
| | |
| void | Mp_StMtP (MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp M_trans, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans) |
| | |
| void | Mp_StPtM (MatrixOp *mwo_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans, const MatrixOp &M_rhs, BLAS_Cpp::Transp M_trans) |
| | |
| void | Mp_StPtMtP (MatrixOp *mwo_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs, const GenPermMatrixSlice &P_rhs2, BLAS_Cpp::Transp P_rhs2_trans) |
| | |
| void | Vp_StMtV (VectorMutable *v_lhs, value_type alpha, const MatrixOp &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta) |
| | |
| void | Vp_StMtV (VectorMutable *v_lhs, value_type alpha, const MatrixOp &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice &sv_rhs2, value_type beta) |
| | |
| void | Vp_StPtMtV (VectorMutable *v_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, const MatrixOp &M_rhs2, BLAS_Cpp::Transp M_rhs2_trans, const Vector &v_rhs3, value_type beta) |
| | |
| void | Vp_StPtMtV (VectorMutable *v_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, const MatrixOp &M_rhs2, BLAS_Cpp::Transp M_rhs2_trans, const SpVectorSlice &sv_rhs3, value_type beta) |
| | |
| value_type | transVtMtV (const Vector &v_rhs1, const MatrixOp &M_rhs2, BLAS_Cpp::Transp trans_rhs2, const Vector &v_rhs3) |
| | |
| value_type | transVtMtV (const SpVectorSlice &sv_rhs1, const MatrixOp &M_rhs2, BLAS_Cpp::Transp trans_rhs2, const SpVectorSlice &sv_rhs3) |
| | |
| void | syr2k (const MatrixOp &M, 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) |
| | |
| void | Mp_StMtM (MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &mwo_rhs1, BLAS_Cpp::Transp trans_rhs1, const MatrixOp &mwo_rhs2, BLAS_Cpp::Transp trans_rhs2, value_type beta) |
| | |
| void | syrk (const MatrixOp &mwo_rhs, BLAS_Cpp::Transp M_trans, value_type alpha, value_type beta, MatrixSymOp *sym_lhs) |
| | |
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 116 of file AbstractLinAlgPack_MatrixOp.hpp.
|
|
Type of matrix norm.
Definition at line 231 of file AbstractLinAlgPack_MatrixOp.hpp. |
|
|
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 Reimplemented in AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MultiVectorMutable, AbstractLinAlgPack::MatrixZero, AbstractLinAlgPack::MultiVectorMutableCols, and AbstractLinAlgPack::MultiVectorMutableDense. Definition at line 51 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
|
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 Reimplemented in AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MultiVectorMutable, AbstractLinAlgPack::MatrixZero, AbstractLinAlgPack::MultiVectorMutableCols, and AbstractLinAlgPack::MultiVectorMutableDense. Definition at line 59 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 Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MultiVectorMutable, AbstractLinAlgPack::MatrixSymDiagStd, AbstractLinAlgPack::MatrixZero, AbstractLinAlgPack::COOMatrixPartitionViewSubclass, AbstractLinAlgPack::MatrixSparseCOORSerial, AbstractLinAlgPack::MatrixSymDiagSparseStd, AbstractLinAlgPack::MatrixWithOpConcreteEncap< M >, AbstractLinAlgPack::MultiVectorMutableCols, AbstractLinAlgPack::MultiVectorMutableDense, and AbstractLinAlgPack::MatrixWithOpConcreteEncap< COOMatrixWithPartitionedView >. Definition at line 67 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
|
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::MatrixOpNonsing, AbstractLinAlgPack::MatrixSymOp, AbstractLinAlgPack::MatrixSymOpNonsing, AbstractLinAlgPack::MultiVectorMutable, and AbstractLinAlgPack::MultiVectorMutableCols. Definition at line 96 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
|
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::MatrixOpNonsing, AbstractLinAlgPack::MatrixSymOp, AbstractLinAlgPack::MatrixSymOpNonsing, AbstractLinAlgPack::MultiVector, and AbstractLinAlgPack::MultiVectorMutableCols. Definition at line 102 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
|
Virtual output function.
The default implementaion just extracts rows one at a time by calling Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixSymIdent, AbstractLinAlgPack::MatrixZero, AbstractLinAlgPack::COOMatrixWithPartitionedViewSubclass, AbstractLinAlgPack::MatrixSparseCOORSerial, AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MatrixSymDiagSparse, and AbstractLinAlgPack::MultiVectorMutableDense. Definition at line 78 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||
|
Compute a norm of this matrix.
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 Definition at line 110 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||
|
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
It is allows for a matrix implementation to return
The default implementation uses the matrix subclass Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MultiVector, and AbstractLinAlgPack::MatrixComposite. Definition at line 167 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||
|
Inlined implementation calls
Definition at line 958 of file AbstractLinAlgPack_MatrixOp.hpp. |
|
||||||||||||||||||||||||||||
|
Create a permuted view:
Postconditions:
The default implementation returns a Definition at line 190 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||||||||||||||
|
Reinitialize a permuted view:
The default implementation simply returns Definition at line 210 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||
|
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 Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MultiVectorMutable, AbstractLinAlgPack::MatrixSymDiagStd, AbstractLinAlgPack::MatrixZero, AbstractLinAlgPack::MatrixOpSerial, and AbstractLinAlgPack::MultiVectorMutableDense. Definition at line 227 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||
|
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 Reimplemented in AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MultiVectorMutable, and AbstractLinAlgPack::MultiVectorMutableDense. Definition at line 234 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||||||
|
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 Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixZero, and AbstractLinAlgPack::MatrixOpSerial. Definition at line 240 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||||||
|
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 Reimplemented in AbstractLinAlgPack::MatrixOpSubView. Definition at line 249 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||||||
|
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 Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixZero, and AbstractLinAlgPack::MatrixOpSerial. Definition at line 258 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||||||
|
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 Reimplemented in AbstractLinAlgPack::MatrixOpSubView. Definition at line 267 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||||||||||||||
|
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 Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixZero, and AbstractLinAlgPack::MatrixOpSerial. Definition at line 276 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||||||||||||||
|
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 Reimplemented in AbstractLinAlgPack::MatrixOpSubView. Definition at line 286 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||||||
|
v_lhs = alpha * op(M_rhs1) * v_rhs2 + beta * v_lhs (BLAS xGEMV)
Implemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixComposite, AbstractLinAlgPack::MatrixSymDiagStd, AbstractLinAlgPack::MatrixSymIdent, AbstractLinAlgPack::MatrixZero, AbstractLinAlgPack::MatrixSparseCOORSerial, AbstractLinAlgPack::MatrixOpSerial, and AbstractLinAlgPack::MultiVectorMutableCols. |
|
||||||||||||||||||||||||
|
v_lhs = alpha * op(M_rhs1) * sv_rhs2 + beta * v_lhs (BLAS xGEMV)
Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixComposite, AbstractLinAlgPack::MatrixSymDiagStd, AbstractLinAlgPack::MatrixZero, AbstractLinAlgPack::MatrixOpSerial, and AbstractLinAlgPack::MultiVectorMutableCols. Definition at line 298 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||||||||||||||
|
v_lhs = alpha * op(P_rhs1) * op(M_rhs2) * v_rhs3 + beta * v_rhs
Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixComposite, AbstractLinAlgPack::MatrixZero, and AbstractLinAlgPack::MatrixOpSerial. Definition at line 317 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||||||||||||||
|
v_lhs = alpha * op(P_rhs1) * op(M_rhs2) * sv_rhs3 + beta * v_rhs
Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixComposite, AbstractLinAlgPack::MatrixZero, and AbstractLinAlgPack::MatrixOpSerial. Definition at line 332 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||
|
result = v_rhs1' * op(M_rhs2) * v_rhs3
Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixZero, and AbstractLinAlgPack::MatrixOpSerial. Definition at line 347 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||
|
result = sv_rhs1' * op(M_rhs2) * sv_rhs3
Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixZero, AbstractLinAlgPack::COOMatrixPartitionViewSubclass, AbstractLinAlgPack::COOMatrixWithPartitionedViewSubclass, and AbstractLinAlgPack::MatrixOpSerial. Definition at line 355 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||||||||||||||||||
|
Perform a specialized rank-2k update of a dense symmetric matrix of the form:.
The reason that this operation is being classified as a level-2 operation is that the total flops should be of
The default implementation is based on Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixZero, and AbstractLinAlgPack::MatrixOpSerial. Definition at line 363 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||||||||||
|
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 Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MultiVector, AbstractLinAlgPack::MatrixZero, AbstractLinAlgPack::MatrixOpSerial, and AbstractLinAlgPack::MultiVectorMutableDense. Definition at line 374 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||||||||||
|
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 Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MultiVector, AbstractLinAlgPack::MatrixZero, and AbstractLinAlgPack::MultiVectorMutableDense. Definition at line 382 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||||||||||
|
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 Reimplemented in AbstractLinAlgPack::MatrixOpSubView. Definition at line 390 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||
|
Perform a rank-k update of a symmetric matrix of the form:.
where
Never call this method directly. Instead use the nonmember function
The default implementation returns Reimplemented in AbstractLinAlgPack::MatrixOpNonsingAggr, AbstractLinAlgPack::MatrixOpSubView, AbstractLinAlgPack::MatrixPermAggr, AbstractLinAlgPack::MatrixSymDiagStd, AbstractLinAlgPack::MatrixZero, AbstractLinAlgPack::MatrixOpSerial, AbstractLinAlgPack::MultiVectorMutableCols, and AbstractLinAlgPack::MultiVectorMutableDense. Definition at line 399 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||||||||||
|
Perform a rank-k update of a symmetric matrix of the form:.
where
Never call this method directly. Instead use the nonmember function
The default implementation returns Definition at line 409 of file AbstractLinAlgPack_MatrixOp.cpp. |
|
||||||||||||
|
mwo_lhs *= alpha.
If |
|
||||||||||||||||||||
|
mwo_lhs += alpha * op(M_rhs) (BLAS xAXPY). Entry point for (poor man's) multiple dispatch.
This method first calls |
|
||||||||||||||||||||||||||||
|
Entry point for (poor man's) multiple dispatch. ToDo: Finish documentation! |
|
||||||||||||||||||||||||||||
|
Entry point for (poor man's) multiple dispatch. ToDo: Finish documentation! |
|
||||||||||||||||||||||||||||||||||||
|
Entry point for (poor man's) multiple dispatch. ToDo: Finish documentation! |
|
||||||||||||||||||||||||||||
|
Definition at line 840 of file AbstractLinAlgPack_MatrixOp.hpp. |
|
||||||||||||||||||||||||||||
|
Definition at line 849 of file AbstractLinAlgPack_MatrixOp.hpp. |
|
||||||||||||||||||||||||||||||||||||
|
Definition at line 858 of file AbstractLinAlgPack_MatrixOp.hpp. |
|
||||||||||||||||||||||||||||||||||||
|
Definition at line 869 of file AbstractLinAlgPack_MatrixOp.hpp. |
|
||||||||||||||||||||
|
Definition at line 880 of file AbstractLinAlgPack_MatrixOp.hpp. |
|
||||||||||||||||||||
|
Definition at line 889 of file AbstractLinAlgPack_MatrixOp.hpp. |
|
||||||||||||||||||||||||||||||||||||||||
|
Definition at line 898 of file AbstractLinAlgPack_MatrixOp.hpp. |
|
||||||||||||||||||||||||||||||||
|
This method first calls
As a last resort, the function attempts to cast |
|
||||||||||||||||||||||||
|
The default implementation returns |
1.3.9.1