42 #ifndef ABSTRACT_LINALG_PACK_MATRIX_WITH_OP_H
43 #define ABSTRACT_LINALG_PACK_MATRIX_WITH_OP_H
47 #include "AbstractLinAlgPack_MatrixBase.hpp"
50 namespace AbstractLinAlgPack {
183 ,
const Vector& v_rhs3, value_type beta
236 #ifndef DOXYGEN_COMPILE
290 virtual void Mt_S( value_type alpha );
318 virtual mat_mut_ptr_t
clone();
325 virtual mat_ptr_t
clone()
const;
338 virtual std::ostream&
output(std::ostream& out)
const;
377 ,
bool allow_replacement =
false
405 const index_type& rl,
const index_type& ru
406 ,
const index_type& cl,
const index_type& cu
456 ,
const index_type row_part[]
459 ,
const index_type col_part[]
494 ,
const index_type row_part[]
497 ,
const index_type col_part[]
504 #ifdef TEMPLATE_FRIENDS_NOT_SUPPORTED
625 ,
const Vector& v_rhs2, value_type beta
639 ,
const Vector& v_rhs3, value_type beta
862 M_rhs1.
Vp_StMtV(v_lhs,alpha,trans_rhs1,v_rhs2,beta);
871 M_rhs1.
Vp_StMtV(v_lhs,alpha,trans_rhs1,sv_rhs2,beta);
879 ,
const Vector& v_rhs3, value_type beta = 1.0
882 M_rhs2.
Vp_StPtMtV(v_lhs,alpha,P_rhs1,P_rhs1_trans,M_rhs2_trans,v_rhs3,beta);
893 M_rhs2.
Vp_StPtMtV(v_lhs,alpha,P_rhs1,P_rhs1_trans,M_rhs2_trans,sv_rhs3,beta);
902 return M_rhs2.
transVtMtV(v_rhs1,trans_rhs2,v_rhs3);
911 return M_rhs2.
transVtMtV(sv_rhs1,trans_rhs2,sv_rhs3);
922 M.
syr2k(M_trans,alpha,P1,P1_trans,P2,P2_trans,beta,symwo_lhs);
944 MatrixOp* mwo_lhs, value_type alpha
947 ,value_type beta = 1.0
957 const MatrixOp &mwo_rhs
961 ,MatrixSymOp *sym_lhs
976 const index_type& rl,
const index_type& ru
977 ,
const index_type& cl,
const index_type& cu
985 #endif // ABSTRACT_LINALG_PACK_MATRIX_WITH_OP_H
Base class for all polymorphic matrices.
friend 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)
Thrown if matrices are not compatible.
Abstract interface for immutable, finite dimensional, coordinate vectors {abstract}.
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=1.0)
mwo_lhs = alpha * op(mwo_rhs1) * op(mwo_rhs2) + beta * mwo_lhs (right) (xGEMM).
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=1.0)
v_lhs = alpha * op(M_rhs1) * sv_rhs2 + beta * v_lhs (BLAS xGEMV)
Sparse Vector Slice class template.
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).
virtual mat_mut_ptr_t clone()
Clone the non-const matrix object (if supported).
friend void Mt_S(MatrixOp *mwo_lhs, value_type alpha)
The induced infinity norm ||M||inf, i.e. max abs row sum.
virtual MatrixOp & operator=(const MatrixOp &mwo_rhs)
M_lhs = mwo_rhs : Virtual assignment operator.
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.
friend value_type transVtMtV(const Vector &v_rhs1, const MatrixOp &M_rhs2, BLAS_Cpp::Transp trans_rhs2, const Vector &v_rhs3)
friend 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)
Interface adding operations specific for a symmetric matrix {abstract}.
virtual void zero_out()
M_lhs = 0 : Zero out the matrix.
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)
symwo_lhs += alpha*op(P1')*op(M)*op(P2) + alpha*op(P2')*op(M')*op(P1) + beta*symwo_lhs ...
virtual std::ostream & output(std::ostream &out) const
Virtual output function.
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)
mwo_lhs += alpha * op(P_rhs1) * op(M_rhs) * op(P_rhs2).
The induced two (i.e. Euclidean norm) norm ||M||2.
The Forbenious norm, i.e. max abs element.
friend 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)
Returned form calc_norm().
void syrk(const MatrixOp &mwo_rhs, BLAS_Cpp::Transp M_trans, value_type alpha, value_type beta, MatrixSymOp *sym_lhs)
Perform a rank-k update of a symmetric matrix of the form:
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)
mwo_lhs += alpha * op(M_rhs) * op(P_rhs).
Thrown if a method is not implemented.
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)
mwo_lhs += alpha * op(P) * op(M_rhs).
void Mp_StM(MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs)
friend void Mp_StM(MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs)
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.
Base class for all matrices that support basic matrix operations.
const MatNorm calc_norm(EMatNormType requested_norm_type=MAT_NORM_1, bool allow_replacement=false) const
Compute a norm of this matrix.
void Mt_S(MatrixOp *mwo_lhs, value_type alpha)
The induced one norm ||M||1, i.e. max abs col sum.
Abstract interface to permutation matrices.
EMatNormType
Type of matrix norm.
friend 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)
Abstract interface for mutable coordinate vectors {abstract}.
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=1.0)
v_lhs = alpha * op(P_rhs1) * op(M_rhs2) * sv_rhs3 + beta * v_rhs
value_type transVtMtV(const SpVectorSlice &sv_rhs1, const MatrixOp &M_rhs2, BLAS_Cpp::Transp trans_rhs2, const SpVectorSlice &sv_rhs3)
result = sv_rhs1' * op(M_rhs2) * sv_rhs3
friend void syrk(const MatrixOp &mwo_rhs, BLAS_Cpp::Transp M_trans, value_type alpha, value_type beta, MatrixSymOp *sym_lhs)
friend 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)
Concrete matrix type to represent general permutation (mapping) matrices.
friend 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)
friend 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)