MOOCHO (Single Doxygen Collection)
Version of the Day
|
Abstract base class for all nonsingular polymorphic matrices that can be used to compute matrix-vector products and solve for linear systems efficiently. More...
#include <AbstractLinAlgPack_MatrixOpNonsing.hpp>
Public Member Functions | |
MatrixOpNonsing & | operator= (const MatrixOpNonsing &M) |
Calls operator=(MatrixOp&) More... | |
Public Member Functions inherited from AbstractLinAlgPack::MatrixOp | |
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... | |
virtual std::ostream & | output (std::ostream &out) const |
Virtual output function. 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... | |
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... | |
Public Member Functions inherited from AbstractLinAlgPack::MatrixNonsing | |
virtual void | V_InvMtV (VectorMutable *v_lhs, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2) const =0 |
v_lhs = inv(op(M_rhs1)) * vs_rhs2 More... | |
virtual void | V_InvMtV (VectorMutable *v_lhs, BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice &sv_rhs2) const |
v_lhs = inv(op(M_rhs1)) * sv_rhs2 More... | |
virtual value_type | transVtInvMtV (const Vector &v_rhs1, BLAS_Cpp::Transp trans_rhs2, const Vector &v_rhs3) const |
result = vs_rhs1' * inv(op(M_rhs2)) * vs_rhs3 More... | |
virtual value_type | transVtInvMtV (const SpVectorSlice &sv_rhs1, BLAS_Cpp::Transp trans_rhs2, const SpVectorSlice &sv_rhs3) const |
result = sv_rhs1' * inv(op(M_rhs2)) * sv_rhs3 More... | |
virtual void | M_StInvMtM (MatrixOp *m_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const MatrixOp &mwo_rhs2, BLAS_Cpp::Transp trans_rhs2) const |
m_lhs = alpha * inv(op(M_rhs1)) * op(mwo_rhs2) (right). More... | |
virtual void | M_StMtInvM (MatrixOp *m_lhs, value_type alpha, const MatrixOp &mwo_rhs1, BLAS_Cpp::Transp trans_rhs1, BLAS_Cpp::Transp trans_rhs2) const |
m_lhs = alpha * op(mwo_rhs1) * inv(op(M_rhs2)) (left). More... | |
Clone | |
virtual mat_mwons_mut_ptr_t | clone_mwons () |
Clone the non-const matrix object (if supported). More... | |
virtual mat_mwons_ptr_t | clone_mwons () const |
Clone the const matrix object (if supported). More... | |
Condition number estimation | |
const MatNorm | calc_cond_num (EMatNormType requested_norm_type=MAT_NORM_1, bool allow_replacement=false) const |
Compute an estimate of the condition number of this matrix. More... | |
Overridden from MatrixOp | |
mat_mut_ptr_t | clone () |
Returns this->clone_mwons() . More... | |
mat_ptr_t | clone () const |
Returns this->clone_mwons() . More... | |
Overridden from MatrixNonsing | |
mat_mns_mut_ptr_t | clone_mns () |
Returns this->clone_mwons() . More... | |
mat_mns_ptr_t | clone_mns () const |
Returns this->clone_mwons() . More... | |
Abstract base class for all nonsingular polymorphic matrices that can be used to compute matrix-vector products and solve for linear systems efficiently.
Definition at line 54 of file AbstractLinAlgPack_MatrixOpNonsing.hpp.
|
virtual |
Clone the non-const matrix object (if supported).
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::MatrixSymOpNonsing.
Definition at line 52 of file AbstractLinAlgPack_MatrixOpNonsing.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::MatrixOpNonsingThyra, and AbstractLinAlgPack::MatrixSymOpNonsing.
Definition at line 58 of file AbstractLinAlgPack_MatrixOpNonsing.cpp.
const MatrixOp::MatNorm AbstractLinAlgPack::MatrixOpNonsing::calc_cond_num | ( | EMatNormType | requested_norm_type = MAT_NORM_1 , |
bool | allow_replacement = false |
||
) | const |
Compute an estimate of the condition number of this matrix.
requested_norm_type | [in] Determines the requested type of norm for the condition number. |
allow_replacement | [in] Determines if the requested norm in specified in norm_type can be replaced with another norm that can be computde by the matrix. |
return.value
gives the value of the condition number in the norm of type return.type
.Postconditions:
allow_replacement==true
, the matrix object must return a computed condition number who's type is given in return.type
. allow_replacement==false
and the underlying matrix object can not compute condition number for the norm requested in norm_type
, then a MethodNotImplemented
exception will be thrown. If the matrix object can an estimate of the condition number for 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 ||inv(M)||1 or ||inv(M)||inf. The algorithm uses some of the refinements in the referenced algorithm by Highman. This algorithm only requires solves and transposed solves so every nonsingular matrix object can implement this method. The default arguments for this function will compute an estimate of the condition number 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 64 of file AbstractLinAlgPack_MatrixOpNonsing.cpp.
|
virtual |
Returns this->clone_mwons()
.
Reimplemented from AbstractLinAlgPack::MatrixOp.
Reimplemented in AbstractLinAlgPack::MatrixOpNonsingThyra, AbstractLinAlgPack::MatrixSymOpNonsing, and AbstractLinAlgPack::MatrixSymOpNonsingSerial.
Definition at line 121 of file AbstractLinAlgPack_MatrixOpNonsing.cpp.
|
virtual |
Returns this->clone_mwons()
.
Reimplemented from AbstractLinAlgPack::MatrixOp.
Reimplemented in AbstractLinAlgPack::MatrixSymOpNonsing, and AbstractLinAlgPack::MatrixSymOpNonsingSerial.
Definition at line 127 of file AbstractLinAlgPack_MatrixOpNonsing.cpp.
|
virtual |
Returns this->clone_mwons()
.
Reimplemented from AbstractLinAlgPack::MatrixNonsing.
Reimplemented in AbstractLinAlgPack::MatrixSymOpNonsing, and AbstractLinAlgPack::MatrixSymOpNonsingSerial.
Definition at line 135 of file AbstractLinAlgPack_MatrixOpNonsing.cpp.
|
virtual |
Returns this->clone_mwons()
.
Reimplemented from AbstractLinAlgPack::MatrixNonsing.
Reimplemented in AbstractLinAlgPack::MatrixSymOpNonsing, and AbstractLinAlgPack::MatrixSymOpNonsingSerial.
Definition at line 141 of file AbstractLinAlgPack_MatrixOpNonsing.cpp.
|
inline |
Calls operator=(MatrixOp&)
Definition at line 152 of file AbstractLinAlgPack_MatrixOpNonsing.hpp.