AbstractLinAlgPack: C++ Interfaces For Vectors, Matrices And Related Linear Algebra Objects
Version of the Day
|
Implementation of a matrix with all zeros. More...
#include <AbstractLinAlgPack_MatrixZero.hpp>
Constructors/initializers | |
MatrixZero (const VectorSpace::space_ptr_t &space_cols=Teuchos::null, const VectorSpace::space_ptr_t &space_rows=Teuchos::null) | |
Calls this->initalize() More... | |
void | initialize (const VectorSpace::space_ptr_t &space_cols, const VectorSpace::space_ptr_t &space_rows) |
Initialize (or initialize) given the columns and rows vector spaces. More... | |
Overridden from MatrixBase | |
size_type | rows () const |
size_type | cols () const |
size_type | nz () const |
Overridden from MatrixOp | |
const VectorSpace & | space_cols () const |
const VectorSpace & | space_rows () const |
void | zero_out () |
void | Mt_S (value_type alpha) |
MatrixOp & | operator= (const MatrixOp &M) |
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 *vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta) const |
void | Vp_StMtV (VectorMutable *vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice &sv_rhs2, value_type beta) const |
void | Vp_StPtMtV (VectorMutable *vs_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 *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 |
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 | 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 |
bool | syrk (BLAS_Cpp::Transp M_trans, value_type alpha, value_type beta, MatrixSymOp *sym_lhs) const |
Additional Inherited Members | |
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... | |
Public Member Functions inherited from AbstractLinAlgPack::MatrixOp | |
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... | |
const MatNorm | calc_norm (EMatNormType requested_norm_type=MAT_NORM_1, bool allow_replacement=false) const |
Compute a norm of this matrix. More... | |
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... | |
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... | |
Protected Member Functions inherited from AbstractLinAlgPack::MatrixOp | |
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 (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 (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 (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 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 (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... | |
Implementation of a matrix with all zeros.
This may seem like a silly class but it is helpful in some circumstances. This class needs to be updated whenever methods are added or removed from MatrixOp
.
Definition at line 56 of file AbstractLinAlgPack_MatrixZero.hpp.
AbstractLinAlgPack::MatrixZero::MatrixZero | ( | const VectorSpace::space_ptr_t & | space_cols = Teuchos::null , |
const VectorSpace::space_ptr_t & | space_rows = Teuchos::null |
||
) |
Calls this->initalize()
Definition at line 53 of file AbstractLinAlgPack_MatrixZero.cpp.
void AbstractLinAlgPack::MatrixZero::initialize | ( | const VectorSpace::space_ptr_t & | space_cols, |
const VectorSpace::space_ptr_t & | space_rows | ||
) |
Initialize (or initialize) given the columns and rows vector spaces.
Preconditions:
space_cols.get() == NULL
] space_rows.get() == NULL
(throw std::invalid_argument
)
[space_cols.get() != NULL
] space_rows.get() != NULL
(throw std::invalid_argument
)
Postconditions:
space_cols.get() == NULL
] this->rows() == 0
space_cols.get() == NULL
] this->cols() == 0
space_cols.get() != NULL
] &this->space_cols() == space_cols.get()
space_cols.get() != NULL
] &this->space_rows() == space_rows.get()
Definition at line 61 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixBase.
Definition at line 78 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixBase.
Definition at line 83 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixBase.
Definition at line 88 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Implements AbstractLinAlgPack::MatrixBase.
Definition at line 95 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Implements AbstractLinAlgPack::MatrixBase.
Definition at line 101 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 107 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 113 of file AbstractLinAlgPack_MatrixZero.cpp.
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 119 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 126 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 134 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 142 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 152 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 162 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Implements AbstractLinAlgPack::MatrixOp.
Definition at line 175 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 184 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 192 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 202 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 212 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 220 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 228 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 240 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 250 of file AbstractLinAlgPack_MatrixZero.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 260 of file AbstractLinAlgPack_MatrixZero.cpp.