#include <AbstractLinAlgPack_MatrixOpSerial.hpp>
Inheritance diagram for AbstractLinAlgPack::MatrixOpSerial:

Level-1 BLAS | |
| virtual void | Mp_StM (DMatrixSlice *gms_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs) const |
| gms_lhs += alpha * op(M_rhs) (BLAS xAXPY) | |
| virtual void | Mp_StMtP (DMatrixSlice *gms_lhs, value_type alpha, BLAS_Cpp::Transp M_trans, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans) const |
| gms_lhs += alpha * op(M_rhs) * op(P_rhs) | |
| virtual void | Mp_StPtM (DMatrixSlice *gms_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs, BLAS_Cpp::Transp P_rhs_trans, BLAS_Cpp::Transp M_trans) const |
| gms_lhs += alpha * op(P) * op(M_rhs) | |
| virtual void | Mp_StPtMtP (DMatrixSlice *gms_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 |
| gms_lhs += alpha * op(P_rhs1) * op(M_rhs) * op(P_rhs2) | |
Level-2 BLAS | |
| virtual void | Vp_StMtV (DVectorSlice *vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const DVectorSlice &vs_rhs2, value_type beta) const =0 |
| vs_lhs = alpha * op(M_rhs1) * vs_rhs2 + beta * vs_lhs (BLAS xGEMV) | |
| virtual void | Vp_StMtV (DVectorSlice *vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice &sv_rhs2, value_type beta) const |
| vs_lhs = alpha * op(M_rhs1) * sv_rhs2 + beta * vs_lhs (BLAS xGEMV) | |
| virtual void | Vp_StPtMtV (DVectorSlice *vs_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, BLAS_Cpp::Transp M_rhs2_trans, const DVectorSlice &vs_rhs3, value_type beta) const |
| vs_lhs = alpha * op(P_rhs1) * op(M_rhs2) * vs_rhs3 + beta * vs_rhs | |
| virtual void | Vp_StPtMtV (DVectorSlice *vs_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 |
| vs_lhs = alpha * op(P_rhs1) * op(M_rhs2) * sv_rhs3 + beta * vs_rhs | |
| virtual value_type | transVtMtV (const DVectorSlice &vs_rhs1, BLAS_Cpp::Transp trans_rhs2, const DVectorSlice &vs_rhs3) const |
| result = vs_rhs1' * op(M_rhs2) * vs_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, DMatrixSliceSym *sym_lhs) const |
| Perform a specialized rank-2k update of a dense symmetric matrix of the form:. | |
Level-3 BLAS | |
| virtual void | Mp_StMtM (DMatrixSlice *gms_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const DMatrixSlice &gms_rhs2, BLAS_Cpp::Transp trans_rhs2, value_type beta) const |
| gms_lhs = alpha * op(M_rhs1) * op(gms_rhs2) + beta * gms_lhs (right) (xGEMM) | |
| virtual void | Mp_StMtM (DMatrixSlice *gms_lhs, value_type alpha, const DMatrixSlice &gms_rhs1, BLAS_Cpp::Transp trans_rhs1, BLAS_Cpp::Transp trans_rhs2, value_type beta) const |
| gms_lhs = alpha * op(gms_rhs1) * op(M_rhs2) + beta * gms_lhs (left) (xGEMM) | |
| virtual void | Mp_StMtM (DMatrixSlice *gms_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const MatrixOpSerial &mwo_rhs2, BLAS_Cpp::Transp trans_rhs2, value_type beta) const |
| gms_lhs = alpha * op(M_rhs1) * op(mwo_rhs2) + beta * gms_lhs (right) (xGEMM) | |
| virtual void | Mp_StMtM (DMatrixSlice *gms_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const DMatrixSliceSym &sym_rhs2, BLAS_Cpp::Transp trans_rhs2, value_type beta) const |
| gms_lhs = alpha * op(M_rhs1) * op(sym_rhs2) + beta * gms_lhs (right) (xSYMM) | |
| virtual void | Mp_StMtM (DMatrixSlice *gms_lhs, value_type alpha, const DMatrixSliceSym &sym_rhs1, BLAS_Cpp::Transp trans_rhs1, BLAS_Cpp::Transp trans_rhs2, value_type beta) const |
| gms_lhs = alpha * op(sym_rhs1) * op(M_rhs2) + beta * gms_lhs (left) (xSYMM) | |
| virtual void | Mp_StMtM (DMatrixSlice *gms_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const DMatrixSliceTri &tri_rhs2, BLAS_Cpp::Transp trans_rhs2, value_type beta) const |
| gms_lhs = alpha * op(M_rhs1) * op(tri_rhs2) + beta * gms_lhs (right) (xTRMM) | |
| virtual void | Mp_StMtM (DMatrixSlice *gms_lhs, value_type alpha, const DMatrixSliceTri &tri_rhs1, BLAS_Cpp::Transp trans_rhs1, BLAS_Cpp::Transp trans_rhs2, value_type beta) const |
| gms_lhs = alpha * op(tri_rhs1) * op(M_rhs2) + beta * gms_lhs (left) (xTRMM) | |
| virtual void | syrk (BLAS_Cpp::Transp M_trans, value_type alpha, value_type beta, DMatrixSliceSym *sym_lhs) const |
| Perform a rank-k update of a dense symmetric matrix of the form:. | |
Overridden from MatrixOp | |
| const VectorSpace & | space_cols () const |
| | |
| const VectorSpace & | space_rows () const |
| | |
| std::ostream & | output (std::ostream &out) const |
| | |
| bool | Mp_StM (MatrixOp *mwo_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs) const |
| | |
| 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 |
| | |
| 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 |
| | |
| 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 |
| | |
| void | Vp_StMtV (VectorMutable *v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta) const |
| | |
| void | Vp_StMtV (VectorMutable *v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice &sv_rhs2, value_type beta) const |
| | |
| 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 |
| | |
| 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 |
| | |
| value_type | transVtMtV (const Vector &v_rhs1, BLAS_Cpp::Transp trans_rhs2, const Vector &v_rhs3) const |
| | |
| value_type | transVtMtV (const SpVectorSlice &sv_rhs1, BLAS_Cpp::Transp trans_rhs2, const SpVectorSlice &sv_rhs3) const |
| | |
| 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 |
| | |
| 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 |
| | |
| bool | syrk (BLAS_Cpp::Transp M_trans, value_type alpha, value_type beta, MatrixSymOp *sym_lhs) const |
| | |
This base class inherites from AbstractLinAlgPack::MatrixOp and provides a means of transition from purely abstract (i.e. parallel, out-of-core) linear algebra to serial linear algebra (with explicit dense vectors and matrices.
This matrix subclass simply uses dynamic_cast<> to map from the abstract interfaces AbstractLinAlgPack::Vector, AbstractLinAlgPack::VectorMutable, AbstractLinAlgPack::MatrixOp to a concreate types using a few standarized serial interfaces VectorGetSparse, VectorMutableDense ... If these interfaces are not supported by the vector and matrix arguments, then exceptions may be thrown.
Note that every method from MatrixOp is not overridden here. Specifically, the methods MatrixOp::zero_out() and MatrixOp::sub_view() are not overridden. These methods have default implementations that should be sufficient for most uses but subclasses may want to provide specialized implementations anyway.
These methods should not be called directly but instead through a set of provided nonmember functions.
ToDo: Finish documentation!
Definition at line 64 of file AbstractLinAlgPack_MatrixOpSerial.hpp.
| void AbstractLinAlgPack::MatrixOpSerial::Mp_StM | ( | DMatrixSlice * | gms_lhs, | |
| value_type | alpha, | |||
| BLAS_Cpp::Transp | trans_rhs | |||
| ) | const [virtual] |
gms_lhs += alpha * op(M_rhs) (BLAS xAXPY)
Definition at line 59 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::Mp_StMtP | ( | DMatrixSlice * | gms_lhs, | |
| value_type | alpha, | |||
| BLAS_Cpp::Transp | M_trans, | |||
| const GenPermMatrixSlice & | P_rhs, | |||
| BLAS_Cpp::Transp | P_rhs_trans | |||
| ) | const [virtual] |
gms_lhs += alpha * op(M_rhs) * op(P_rhs)
Definition at line 79 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::Mp_StPtM | ( | DMatrixSlice * | gms_lhs, | |
| value_type | alpha, | |||
| const GenPermMatrixSlice & | P_rhs, | |||
| BLAS_Cpp::Transp | P_rhs_trans, | |||
| BLAS_Cpp::Transp | M_trans | |||
| ) | const [virtual] |
gms_lhs += alpha * op(P) * op(M_rhs)
Definition at line 88 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::Mp_StPtMtP | ( | DMatrixSlice * | gms_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 [virtual] |
gms_lhs += alpha * op(P_rhs1) * op(M_rhs) * op(P_rhs2)
Definition at line 97 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| virtual void AbstractLinAlgPack::MatrixOpSerial::Vp_StMtV | ( | DVectorSlice * | vs_lhs, | |
| value_type | alpha, | |||
| BLAS_Cpp::Transp | trans_rhs1, | |||
| const DVectorSlice & | vs_rhs2, | |||
| value_type | beta | |||
| ) | const [pure virtual] |
vs_lhs = alpha * op(M_rhs1) * vs_rhs2 + beta * vs_lhs (BLAS xGEMV)
Implemented in AbstractLinAlgPack::MatrixSymDiagSparse, and AbstractLinAlgPack::MultiVectorMutableDense.
| void AbstractLinAlgPack::MatrixOpSerial::Vp_StMtV | ( | DVectorSlice * | vs_lhs, | |
| value_type | alpha, | |||
| BLAS_Cpp::Transp | trans_rhs1, | |||
| const SpVectorSlice & | sv_rhs2, | |||
| value_type | beta | |||
| ) | const [virtual] |
vs_lhs = alpha * op(M_rhs1) * sv_rhs2 + beta * vs_lhs (BLAS xGEMV)
Reimplemented in AbstractLinAlgPack::MultiVectorMutableDense.
Definition at line 109 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::Vp_StPtMtV | ( | DVectorSlice * | vs_lhs, | |
| value_type | alpha, | |||
| const GenPermMatrixSlice & | P_rhs1, | |||
| BLAS_Cpp::Transp | P_rhs1_trans, | |||
| BLAS_Cpp::Transp | M_rhs2_trans, | |||
| const DVectorSlice & | vs_rhs3, | |||
| value_type | beta | |||
| ) | const [virtual] |
vs_lhs = alpha * op(P_rhs1) * op(M_rhs2) * vs_rhs3 + beta * vs_rhs
Definition at line 132 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::Vp_StPtMtV | ( | DVectorSlice * | vs_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 [virtual] |
vs_lhs = alpha * op(P_rhs1) * op(M_rhs2) * sv_rhs3 + beta * vs_rhs
Definition at line 146 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| value_type AbstractLinAlgPack::MatrixOpSerial::transVtMtV | ( | const DVectorSlice & | vs_rhs1, | |
| BLAS_Cpp::Transp | trans_rhs2, | |||
| const DVectorSlice & | vs_rhs3 | |||
| ) | const [virtual] |
result = vs_rhs1' * op(M_rhs2) * vs_rhs3
Definition at line 160 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::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, | |||
| DMatrixSliceSym * | sym_lhs | |||
| ) | const [virtual] |
Perform a specialized rank-2k update of a dense symmetric matrix of the form:.
sym_lhs += alpha*op(P1')*op(M)*op(P2) + alpha*op(P2')*op(M')*op(P1) + beta*sym_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.
Definition at line 187 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::Mp_StMtM | ( | DMatrixSlice * | gms_lhs, | |
| value_type | alpha, | |||
| BLAS_Cpp::Transp | trans_rhs1, | |||
| const DMatrixSlice & | gms_rhs2, | |||
| BLAS_Cpp::Transp | trans_rhs2, | |||
| value_type | beta | |||
| ) | const [virtual] |
gms_lhs = alpha * op(M_rhs1) * op(gms_rhs2) + beta * gms_lhs (right) (xGEMM)
Definition at line 284 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::Mp_StMtM | ( | DMatrixSlice * | gms_lhs, | |
| value_type | alpha, | |||
| const DMatrixSlice & | gms_rhs1, | |||
| BLAS_Cpp::Transp | trans_rhs1, | |||
| BLAS_Cpp::Transp | trans_rhs2, | |||
| value_type | beta | |||
| ) | const [virtual] |
gms_lhs = alpha * op(gms_rhs1) * op(M_rhs2) + beta * gms_lhs (left) (xGEMM)
Definition at line 302 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::Mp_StMtM | ( | DMatrixSlice * | gms_lhs, | |
| value_type | alpha, | |||
| BLAS_Cpp::Transp | trans_rhs1, | |||
| const MatrixOpSerial & | mwo_rhs2, | |||
| BLAS_Cpp::Transp | trans_rhs2, | |||
| value_type | beta | |||
| ) | const [virtual] |
gms_lhs = alpha * op(M_rhs1) * op(mwo_rhs2) + beta * gms_lhs (right) (xGEMM)
Definition at line 320 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::Mp_StMtM | ( | DMatrixSlice * | gms_lhs, | |
| value_type | alpha, | |||
| BLAS_Cpp::Transp | trans_rhs1, | |||
| const DMatrixSliceSym & | sym_rhs2, | |||
| BLAS_Cpp::Transp | trans_rhs2, | |||
| value_type | beta | |||
| ) | const [virtual] |
gms_lhs = alpha * op(M_rhs1) * op(sym_rhs2) + beta * gms_lhs (right) (xSYMM)
Definition at line 344 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::Mp_StMtM | ( | DMatrixSlice * | gms_lhs, | |
| value_type | alpha, | |||
| const DMatrixSliceSym & | sym_rhs1, | |||
| BLAS_Cpp::Transp | trans_rhs1, | |||
| BLAS_Cpp::Transp | trans_rhs2, | |||
| value_type | beta | |||
| ) | const [virtual] |
gms_lhs = alpha * op(sym_rhs1) * op(M_rhs2) + beta * gms_lhs (left) (xSYMM)
Definition at line 351 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::Mp_StMtM | ( | DMatrixSlice * | gms_lhs, | |
| value_type | alpha, | |||
| BLAS_Cpp::Transp | trans_rhs1, | |||
| const DMatrixSliceTri & | tri_rhs2, | |||
| BLAS_Cpp::Transp | trans_rhs2, | |||
| value_type | beta | |||
| ) | const [virtual] |
gms_lhs = alpha * op(M_rhs1) * op(tri_rhs2) + beta * gms_lhs (right) (xTRMM)
Definition at line 357 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::Mp_StMtM | ( | DMatrixSlice * | gms_lhs, | |
| value_type | alpha, | |||
| const DMatrixSliceTri & | tri_rhs1, | |||
| BLAS_Cpp::Transp | trans_rhs1, | |||
| BLAS_Cpp::Transp | trans_rhs2, | |||
| value_type | beta | |||
| ) | const [virtual] |
gms_lhs = alpha * op(tri_rhs1) * op(M_rhs2) + beta * gms_lhs (left) (xTRMM)
Definition at line 364 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::syrk | ( | BLAS_Cpp::Transp | M_trans, | |
| value_type | alpha, | |||
| value_type | beta, | |||
| DMatrixSliceSym * | sym_lhs | |||
| ) | const [virtual] |
Perform a rank-k update of a dense symmetric matrix of the form:.
sym_lhs += op(M)*op(M') + beta*sym_lhs#
The default implementation is based on Vp_StMtV(...)#. Of course in situations where this default implemention is inefficient the subclass should override this method.
Definition at line 370 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| const VectorSpace & AbstractLinAlgPack::MatrixOpSerial::space_cols | ( | ) | const [virtual] |
Implements AbstractLinAlgPack::MatrixBase.
Definition at line 426 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| const VectorSpace & AbstractLinAlgPack::MatrixOpSerial::space_rows | ( | ) | const [virtual] |
Implements AbstractLinAlgPack::MatrixBase.
Reimplemented in AbstractLinAlgPack::MatrixSymOpSerial.
Definition at line 434 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| std::ostream & AbstractLinAlgPack::MatrixOpSerial::output | ( | std::ostream & | out | ) | const [virtual] |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Reimplemented in AbstractLinAlgPack::MatrixSymDiagSparse, and AbstractLinAlgPack::MultiVectorMutableDense.
Definition at line 442 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| bool AbstractLinAlgPack::MatrixOpSerial::Mp_StM | ( | MatrixOp * | mwo_lhs, | |
| value_type | alpha, | |||
| BLAS_Cpp::Transp | trans_rhs | |||
| ) | const [virtual] |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Reimplemented in AbstractLinAlgPack::MultiVectorMutableDense.
Definition at line 448 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| bool AbstractLinAlgPack::MatrixOpSerial::Mp_StMtP | ( | MatrixOp * | mwo_lhs, | |
| value_type | alpha, | |||
| BLAS_Cpp::Transp | M_trans, | |||
| const GenPermMatrixSlice & | P_rhs, | |||
| BLAS_Cpp::Transp | P_rhs_trans | |||
| ) | const [virtual] |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 461 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| bool AbstractLinAlgPack::MatrixOpSerial::Mp_StPtM | ( | MatrixOp * | mwo_lhs, | |
| value_type | alpha, | |||
| const GenPermMatrixSlice & | P_rhs, | |||
| BLAS_Cpp::Transp | P_rhs_trans, | |||
| BLAS_Cpp::Transp | M_trans | |||
| ) | const [virtual] |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 475 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| bool AbstractLinAlgPack::MatrixOpSerial::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 [virtual] |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 489 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::Vp_StMtV | ( | VectorMutable * | v_lhs, | |
| value_type | alpha, | |||
| BLAS_Cpp::Transp | trans_rhs1, | |||
| const Vector & | v_rhs2, | |||
| value_type | beta | |||
| ) | const [virtual] |
Implements AbstractLinAlgPack::MatrixOp.
Definition at line 504 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::Vp_StMtV | ( | VectorMutable * | v_lhs, | |
| value_type | alpha, | |||
| BLAS_Cpp::Transp | trans_rhs1, | |||
| const SpVectorSlice & | sv_rhs2, | |||
| value_type | beta | |||
| ) | const [virtual] |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 513 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::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 [virtual] |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 521 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::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 [virtual] |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 532 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| value_type AbstractLinAlgPack::MatrixOpSerial::transVtMtV | ( | const Vector & | v_rhs1, | |
| BLAS_Cpp::Transp | trans_rhs2, | |||
| const Vector & | v_rhs3 | |||
| ) | const [virtual] |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 542 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| value_type AbstractLinAlgPack::MatrixOpSerial::transVtMtV | ( | const SpVectorSlice & | sv_rhs1, | |
| BLAS_Cpp::Transp | trans_rhs2, | |||
| const SpVectorSlice & | sv_rhs3 | |||
| ) | const [virtual] |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 169 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| void AbstractLinAlgPack::MatrixOpSerial::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 [virtual] |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 551 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| bool AbstractLinAlgPack::MatrixOpSerial::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] |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Reimplemented in AbstractLinAlgPack::MultiVectorMutableDense.
Definition at line 569 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
| bool AbstractLinAlgPack::MatrixOpSerial::syrk | ( | BLAS_Cpp::Transp | M_trans, | |
| value_type | alpha, | |||
| value_type | beta, | |||
| MatrixSymOp * | sym_lhs | |||
| ) | const [virtual] |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Reimplemented in AbstractLinAlgPack::MultiVectorMutableDense.
Definition at line 602 of file AbstractLinAlgPack_MatrixOpSerial.cpp.
1.4.7