MOOCHO (Single Doxygen Collection)
Version of the Day
|
Abstract base class for all nonsingular polymorphic matrices that can solve for linear system with but it may not be convienent to compute matrix vector products {abstract}. More...
#include <AbstractLinAlgPack_MatrixNonsing.hpp>
Classes | |
class | SingularMatrix |
This exception will be thrown if it turns out at runtime that the matrix is numerically singular. More... | |
Friends | |
void | V_InvMtV (VectorMutable *v_lhs, const MatrixNonsing &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2) |
void | V_InvMtV (VectorMutable *v_lhs, const MatrixNonsing &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice &sv_rhs2) |
value_type | transVtInvMtV (const Vector &v_rhs1, const MatrixNonsing &M_rhs2, BLAS_Cpp::Transp trans_rhs2, const Vector &v_rhs3) |
value_type | transVtInvMtV (const SpVectorSlice &sv_rhs1, const MatrixNonsing &M_rhs2, BLAS_Cpp::Transp trans_rhs2, const SpVectorSlice &sv_rhs3) |
void | M_StInvMtM (MatrixOp *m_lhs, value_type alpha, const MatrixNonsing &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const MatrixOp &mwo_rhs2, BLAS_Cpp::Transp trans_rhs2) |
void | M_StMtInvM (MatrixOp *m_lhs, value_type alpha, const MatrixOp &mwo_rhs1, BLAS_Cpp::Transp trans_rhs1, const MatrixNonsing &M_rhs2, BLAS_Cpp::Transp trans_rhs2) |
Clone | |
virtual mat_mns_mut_ptr_t | clone_mns () |
Clone the non-const matrix object (if supported). More... | |
virtual mat_mns_ptr_t | clone_mns () const |
Clone the const matrix object (if supported). More... | |
Level-2 BLAS | |
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... | |
Level-3 BLAS | |
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... | |
Additional Inherited Members | |
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... | |
Abstract base class for all nonsingular polymorphic matrices that can solve for linear system with but it may not be convienent to compute matrix vector products {abstract}.
The operations supported are:
Level-2 BLAS
v_lhs = inv(op(M_rhs1)) * vs_rhs2
v_lhs = inv(op(M_rhs1)) * sv_rhs2
result = v_rhs1' * inv(op(M_rhs2)) * v_rhs3
result = sv_rhs1' * inv(op(M_rhs2)) * sv_rhs3
Level-3 BLAS
m_lhs = alpha * inv(op(M_rhs1)) * op(mwo_rhs2) (right)
m_lhs = alpha * op(mwo_rhs1) * inv(op(M_rhs2)) (left)
For the solve operations, the lhs and rhs arguments may not be the same in general so don't assume that you can alias the lhs with the rhs and get correct results.
Any nonsingular matrix abstraction that can be used to solve for nonlinear systems should also be able to support the MatrixOp
interface. Therefore, this interface is more of an implementation artifact than an a legitimate domain abstraction. However, some linear solvers that can implement this interface, can not easily implement the MatrixOp
interface and therefore this interface is justified. A general client should never use this interface directly. Instead, the combined interface MatrixOpNonsing
should be used with fully formed matrix abstractions.
All these Level-2 and Level-3 BLAS operations have default implementations based on the Level-2 BLAS operations:
v_lhs = inv(op(M_rhs1)) * vs_rhs2
which allows for fast prototyping of new matrix subclasses.
The member functions should not be called directly but instead through the provided non-member functions.
The multiple dispatch approach taken in MatrixOp
is not taken in this interface. This is because it is considered here that the nonsingular matrix takes procedence of a general matrix arguemnt and we can not expect a general matrix to know how to solve for a linear system with some other nonsigular matrix.
Definition at line 97 of file AbstractLinAlgPack_MatrixNonsing.hpp.