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
Friends | List of all members
AbstractLinAlgPack::DirectSparseSolverImp::BasisMatrixImp Class Referenceabstract

Implementation node subclass that combines factorization structure and factorization nonzeros into a single basis matrix object. More...

#include <AbstractLinAlgPack_DirectSparseSolverImp.hpp>

Inheritance diagram for AbstractLinAlgPack::DirectSparseSolverImp::BasisMatrixImp:
Inheritance graph
[legend]

Friends

class DirectSparseSolverImp
 Allow only DirectSparseSolverImp to initialize objects of this type. Important !!!!! Even though DirectSparseSolverImp has access to these private members it is strictly not to access them directly !!!! More...
 

Public types

typedef Teuchos::RCP
< FactorizationNonzeros
fact_nonzeros_ptr_t
 

Access

virtual const fact_nonzeros_ptr_tget_fact_nonzeros () const
 Return a reference to a smart pointer to the object that represents the factorization nonzeros. More...
 

Overridden from MatrixBase

const VectorSpacespace_cols () const
 
const VectorSpacespace_rows () const
 
size_type rows () const
 
size_type cols () const
 

Overridden from MatrixNonsinguar

mat_mns_mut_ptr_t clone_mns ()
 

Overridden from BasisMatrix

virtual const fact_struc_ptr_tget_fact_struc () const
 

Constructors/initializers

 BasisMatrixImp ()
 Default initializers to uninitialized. More...
 
 BasisMatrixImp (size_type dim, const fact_struc_ptr_t &fact_struc, const fact_nonzeros_ptr_t &fact_nonzeros)
 Calls this->initialize() More...
 
virtual void initialize (size_type dim, const fact_struc_ptr_t &fact_struc, const fact_nonzeros_ptr_t &fact_nonzeros)
 Initialize given initialized factorization structure and factorization nonzeros objects. More...
 
void set_uninitialized ()
 Make uninitialized. More...
 

Factory methods to be overridden by subclasses

virtual Teuchos::RCP
< BasisMatrixImp
create_matrix () const =0
 Called by this->clone-mns(). More...
 

Additional Inherited Members

- Public Types inherited from AbstractLinAlgPack::DirectSparseSolver::BasisMatrix
typedef Teuchos::RCP
< FactorizationStructure
fact_struc_ptr_t
 
- Public Member Functions inherited from AbstractLinAlgPack::MatrixNonsing
virtual mat_mns_ptr_t clone_mns () const
 Clone the const matrix object (if supported). More...
 
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 size_type nz () const
 Return the number of nonzero elements in the matrix. More...
 

Detailed Description

Implementation node subclass that combines factorization structure and factorization nonzeros into a single basis matrix object.

The only methods that subclasses must override for a complete basis matrix object are:

AbstractLinAlgPack::MatrixNonsing::V_InvMtV(), and create_matrix().

Note that is if very important that subclasses do not maintain their own copies of smart reference counting pointer objects to the factorization structure or factorization nonzeros. Instead, the methods get_fact_struc() and get_fact_nonzeros() should be called inside of the implementation of the the V_InvMtV() method and then dynamic_cast<> used to get at the concreate objects.

Definition at line 98 of file AbstractLinAlgPack_DirectSparseSolverImp.hpp.

Member Typedef Documentation

Constructor & Destructor Documentation

AbstractLinAlgPack::DirectSparseSolverImp::BasisMatrixImp::BasisMatrixImp ( )
protected

Default initializers to uninitialized.

Definition at line 54 of file AbstractLinAlgPack_DirectSparseSolverImp.cpp.

AbstractLinAlgPack::DirectSparseSolverImp::BasisMatrixImp::BasisMatrixImp ( size_type  dim,
const fact_struc_ptr_t fact_struc,
const fact_nonzeros_ptr_t fact_nonzeros 
)
protected

Calls this->initialize()

Definition at line 58 of file AbstractLinAlgPack_DirectSparseSolverImp.cpp.

Member Function Documentation

const DirectSparseSolverImp::BasisMatrixImp::fact_nonzeros_ptr_t & AbstractLinAlgPack::DirectSparseSolverImp::BasisMatrixImp::get_fact_nonzeros ( ) const
virtual

Return a reference to a smart pointer to the object that represents the factorization nonzeros.

Returning a reference to a RCP<> object verses returning a RCP<> object itself is critical so that we can rely on RCP<>::count() to tell us how many clients have a reference to this object.

Definition at line 94 of file AbstractLinAlgPack_DirectSparseSolverImp.cpp.

const VectorSpace & AbstractLinAlgPack::DirectSparseSolverImp::BasisMatrixImp::space_cols ( ) const
virtual
const VectorSpace & AbstractLinAlgPack::DirectSparseSolverImp::BasisMatrixImp::space_rows ( ) const
virtual
size_type AbstractLinAlgPack::DirectSparseSolverImp::BasisMatrixImp::rows ( ) const
virtual

Reimplemented from AbstractLinAlgPack::MatrixBase.

Definition at line 111 of file AbstractLinAlgPack_DirectSparseSolverImp.cpp.

size_type AbstractLinAlgPack::DirectSparseSolverImp::BasisMatrixImp::cols ( ) const
virtual

Reimplemented from AbstractLinAlgPack::MatrixBase.

Definition at line 116 of file AbstractLinAlgPack_DirectSparseSolverImp.cpp.

MatrixNonsing::mat_mns_mut_ptr_t AbstractLinAlgPack::DirectSparseSolverImp::BasisMatrixImp::clone_mns ( )
virtual
const DirectSparseSolver::BasisMatrix::fact_struc_ptr_t & AbstractLinAlgPack::DirectSparseSolverImp::BasisMatrixImp::get_fact_struc ( ) const
virtual
void AbstractLinAlgPack::DirectSparseSolverImp::BasisMatrixImp::initialize ( size_type  dim,
const fact_struc_ptr_t fact_struc,
const fact_nonzeros_ptr_t fact_nonzeros 
)
protectedvirtual

Initialize given initialized factorization structure and factorization nonzeros objects.

Definition at line 67 of file AbstractLinAlgPack_DirectSparseSolverImp.cpp.

void AbstractLinAlgPack::DirectSparseSolverImp::BasisMatrixImp::set_uninitialized ( )
protected

Make uninitialized.

Postconditions:

Definition at line 85 of file AbstractLinAlgPack_DirectSparseSolverImp.cpp.

virtual Teuchos::RCP<BasisMatrixImp> AbstractLinAlgPack::DirectSparseSolverImp::BasisMatrixImp::create_matrix ( ) const
protectedpure virtual

Called by this->clone-mns().

Friends And Related Function Documentation

friend class DirectSparseSolverImp
friend

Allow only DirectSparseSolverImp to initialize objects of this type. Important !!!!! Even though DirectSparseSolverImp has access to these private members it is strictly not to access them directly !!!!

Definition at line 205 of file AbstractLinAlgPack_DirectSparseSolverImp.hpp.


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