AbstractLinAlgPack: C++ Interfaces For Vectors, Matrices And Related Linear Algebra Objects  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
List of all members
AbstractLinAlgPack::MatrixSymNonsing Class Reference

Abstract base class for all polymorphic symmetrix nonsingular matrices that can be used to solve for linear systems relatively efficently. More...

#include <AbstractLinAlgPack_MatrixSymNonsing.hpp>

Inheritance diagram for AbstractLinAlgPack::MatrixSymNonsing:
Inheritance graph
[legend]

Public types

enum  EMatrixDummyArg
 

Friends

void M_StMtInvMtM (MatrixSymOp *sym_gms_lhs, value_type alpha, const MatrixOp &mwo, BLAS_Cpp::Transp mwo_trans, const MatrixSymNonsing &mswof, EMatrixDummyArg mwo_rhs)
 

Clone

virtual mat_msns_mut_ptr_t clone_msns ()
 Clone the non-const matrix object (if supported). More...
 
virtual mat_msns_ptr_t clone_msns () const
 Clone the const matrix object (if supported). More...
 

Level-3

virtual void M_StMtInvMtM (MatrixSymOp *symwo_lhs, value_type alpha, const MatrixOp &mwo, BLAS_Cpp::Transp mwo_trans, EMatrixDummyArg) const
 symwo_lhs = alpha * op(mwo) * inv(M) * op(mwo)'. More...
 
mat_mns_mut_ptr_t clone_mns ()
 Returns this->clone_msns(). More...
 
mat_mns_ptr_t clone_mns () const
 Returns this->clone_msns(). More...
 

Additional Inherited Members

- 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...
 
- Public Member Functions inherited from AbstractLinAlgPack::MatrixBase
virtual ~MatrixBase ()
 Virtual destructor. More...
 
virtual const VectorSpacespace_cols () const =0
 Vector space for vectors that are compatible with the columns of the matrix. More...
 
virtual const VectorSpacespace_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...
 

Detailed Description

Abstract base class for all polymorphic symmetrix nonsingular matrices that can be used to solve for linear systems relatively efficently.

This interface defines a single addition method to those found in MatrixNonsing:

symwo_lhs = alpha * op(mwo) * inv(M) * op(mwo)'

The reason that this method could not be defined in the MatrixNonsing interface is that the lhs matrix matrix argument symwo_lhs is only guaranteed to be symmetric if the rhs matrix argument M (which is this matrix) is guaranteed to be symmetric. Since a MatrixNonsing matrix object may be unsymmetric, it can not implement this operation, only a symmetric nonsingular matrix can.

Any symmetric nonsingular matrix abstraction that can be used to solve for nonlinear systems should also be able to support the MatrixSymOp interface. Therefore, this interface is more of an implementation artifact than an a legitimate domain abstraction. However, some symmetric linear solvers that can implement this interface, can not easily implement the MatrixSymOp interface and therefore this interface is justified. A general client should never use this interface directly. Instead, the combined interface MatrixSymOpNonsing should be used with fully formed symmetric matrix abstractions.

Clients should use the provided non-member functions to call the methods and not the methods themselves.

Definition at line 74 of file AbstractLinAlgPack_MatrixSymNonsing.hpp.

Member Enumeration Documentation

Definition at line 89 of file AbstractLinAlgPack_MatrixSymNonsing.hpp.

Member Function Documentation

MatrixSymNonsing::mat_msns_mut_ptr_t AbstractLinAlgPack::MatrixSymNonsing::clone_msns ( )
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 50 of file AbstractLinAlgPack_MatrixSymNonsing.cpp.

MatrixSymNonsing::mat_msns_ptr_t AbstractLinAlgPack::MatrixSymNonsing::clone_msns ( ) const
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::MatrixSymOpNonsing.

Definition at line 56 of file AbstractLinAlgPack_MatrixSymNonsing.cpp.

void AbstractLinAlgPack::MatrixSymNonsing::M_StMtInvMtM ( MatrixSymOp symwo_lhs,
value_type  alpha,
const MatrixOp mwo,
BLAS_Cpp::Transp  mwo_trans,
EMatrixDummyArg   
) const
protectedvirtual

symwo_lhs = alpha * op(mwo) * inv(M) * op(mwo)'.

The default implementation is based on the operation M_StInvMtM(...) assuming that this #M# is a symmetric matrix. For an efficient implementation (for this = L*L' for instance) the subclass may want to override this function.

Reimplemented in AbstractLinAlgPack::MatrixSymNonsingSerial.

Definition at line 61 of file AbstractLinAlgPack_MatrixSymNonsing.cpp.

MatrixSymNonsing::mat_mns_mut_ptr_t AbstractLinAlgPack::MatrixSymNonsing::clone_mns ( )
virtual

Returns this->clone_msns().

Overridden from MatrixNonsing

Reimplemented from AbstractLinAlgPack::MatrixNonsing.

Reimplemented in AbstractLinAlgPack::MatrixSymOpNonsing.

Definition at line 71 of file AbstractLinAlgPack_MatrixSymNonsing.cpp.

MatrixSymNonsing::mat_mns_ptr_t AbstractLinAlgPack::MatrixSymNonsing::clone_mns ( ) const
virtual

Returns this->clone_msns().

Reimplemented from AbstractLinAlgPack::MatrixNonsing.

Reimplemented in AbstractLinAlgPack::MatrixSymOpNonsing.

Definition at line 77 of file AbstractLinAlgPack_MatrixSymNonsing.cpp.

Friends And Related Function Documentation

void M_StMtInvMtM ( MatrixSymOp sym_gms_lhs,
value_type  alpha,
const MatrixOp mwo,
BLAS_Cpp::Transp  mwo_trans,
const MatrixSymNonsing mswof,
MatrixSymNonsing::EMatrixDummyArg  mwo_rhs 
)
friend

Definition at line 167 of file AbstractLinAlgPack_MatrixSymNonsing.hpp.


The documentation for this class was generated from the following files: