Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Member Functions | Private Attributes | List of all members
Ifpack_TriDiContainer Class Reference

Ifpack_TriDiContainer: a class to define containers for dense matrices. More...

#include <Ifpack_TriDiContainer.h>

Inheritance diagram for Ifpack_TriDiContainer:
Inheritance graph
[legend]

Public Member Functions

virtual double InitializeFlops () const
 Returns the flops in Initialize(). More...
 
virtual double ComputeFlops () const
 Returns the flops in Compute(). More...
 
virtual double ApplyFlops () const
 Returns the flops in Apply(). More...
 
virtual double ApplyInverseFlops () const
 Returns the flops in ApplyInverse(). More...
 
virtual std::ostream & Print (std::ostream &os) const
 Prints basic information on iostream. This function is used by operator<<. More...
 
- Public Member Functions inherited from Ifpack_Container
virtual ~Ifpack_Container ()
 Destructor. More...
 

Private Member Functions

virtual int Extract (const Epetra_RowMatrix &Matrix_in)
 Extract the submatrices identified by the ID set int ID(). More...
 

Private Attributes

int NumRows_
 Number of rows in the container. More...
 
int NumVectors_
 Number of vectors in the container. More...
 
Ifpack_SerialTriDiMatrix NonFactoredMatrix_
 TriDi matrix, that contains the non-factored matrix. More...
 
Ifpack_SerialTriDiMatrix Matrix_
 TriDi matrix. More...
 
Epetra_SerialDenseMatrix LHS_
 SerialDense vector representing the LHS. More...
 
Epetra_SerialDenseMatrix RHS_
 SerialDense vector representing the RHS. More...
 
Ifpack_SerialTriDiSolver Solver_
 TriDi solver (solution will be get using LAPACK). More...
 
Epetra_IntSerialDenseVector ID_
 Sets of local rows. More...
 
bool KeepNonFactoredMatrix_
 If true, keeps a copy of the non-factored matrix. More...
 
bool IsInitialized_
 If true, the container has been successfully initialized. More...
 
bool IsComputed_
 If true, the container has been successfully computed. More...
 
std::string Label_
 Label for this object. More...
 
double ComputeFlops_
 Flops in Compute(). More...
 
double ApplyFlops_
 Flops in Apply(). More...
 
double ApplyInverseFlops_
 Flops in ApplyInverse(). More...
 
 Ifpack_TriDiContainer (const int NumRows_in, const int NumVectors_in=1)
 Default constructor. More...
 
 Ifpack_TriDiContainer (const Ifpack_TriDiContainer &rhs)
 Copy constructor. More...
 
virtual ~Ifpack_TriDiContainer ()
 Destructor. More...
 
Ifpack_TriDiContaineroperator= (const Ifpack_TriDiContainer &rhs)
 Operator=. More...
 
virtual int NumRows () const
 Returns the number of rows of the matrix and LHS/RHS. More...
 
virtual int NumVectors () const
 Returns the number of vectors in LHS/RHS. More...
 
virtual int SetNumVectors (const int NumVectors_in)
 Sets the number of vectors for LHS/RHS. More...
 
virtual double & LHS (const int i, const int Vector=0)
 Returns the i-th component of the vector Vector of LHS. More...
 
virtual double & RHS (const int i, const int Vector=0)
 Returns the i-th component of the vector Vector of RHS. More...
 
virtual int & ID (const int i)
 Returns the ID associated to local row i. More...
 
virtual int SetMatrixElement (const int row, const int col, const double value)
 Set the matrix element (row,col) to value. More...
 
virtual int SetParameters (Teuchos::ParameterList &)
 Sets all necessary parameters. More...
 
virtual bool IsInitialized () const
 Returns true is the container has been successfully initialized. More...
 
virtual bool IsComputed () const
 Returns true is the container has been successfully computed. More...
 
virtual const char * Label () const
 Returns the label of this container. More...
 
virtual int SetKeepNonFactoredMatrix (const bool flag)
 If flag is true, keeps a copy of the non-factored matrix. More...
 
virtual bool KeepNonFactoredMatrix () const
 Returns KeepNonFactoredMatrix_. More...
 
virtual const
Epetra_SerialDenseMatrix
LHS () const
 Returns the dense vector containing the LHS. More...
 
virtual const
Epetra_SerialDenseMatrix
RHS () const
 Returns the dense vector containing the RHS. More...
 
virtual const
Ifpack_SerialTriDiMatrix
Matrix () const
 Returns the dense matrix or its factors. More...
 
virtual const
Ifpack_SerialTriDiMatrix
NonFactoredMatrix () const
 Returns the non-factored dense matrix (only if stored). More...
 
virtual const
Epetra_IntSerialDenseVector
ID () const
 Returns the integer dense vector of IDs. More...
 
virtual int Initialize ()
 Initialize the container. More...
 
virtual int Compute (const Epetra_RowMatrix &Matrix_in)
 Finalizes the linear system matrix and prepares for the application of the inverse. More...
 
virtual int Apply ()
 Apply the matrix to RHS, results are stored in LHS. More...
 
virtual int ApplyInverse ()
 Apply the inverse of the matrix to RHS, results are stored in LHS. More...
 

Detailed Description

Ifpack_TriDiContainer: a class to define containers for dense matrices.

To understand what an IFPACK container is, please refer to the documentation of the pure virtual class Ifpack_Container. Currently, containers are used by class Ifpack_BlockRelaxation.

Using block methods, one needs to store all diagonal blocks and to be also to apply the inverse of each diagonal block. Using class Ifpack_TriDiContainer, one can store the blocks as dense matrices, which can be advantageous when the blocks are small. Otherwise, class Ifpack_SparseContainer is probably more appropriate.

A typical use of a container is as follows:

...
// local matrix of (5,5), with two vectors for solution and rhs.
Ifpack_Container* Container = new
// assign local rows 1, 5, 12, 13, 16 to this container
Container(0) = 1;
Container(1) = 5;
Container(2) = 12;
Container(3) = 13;
Container(4) = 16;
// Now extract the submatrix corresponding to rows and columns:
// 1. initialize the container.
Container.Initialize();
// 2. extract matrix values from an Epetra_RowMatrix A,
// and compute LU factors of the submatrix identified by rows
// and columns 1, 5, 12, 13 and 16 using LAPACK
Container.Compute(A);
// We can set the RHS as follows:
Container.RHS(0) = 1.0;
Container.RHS(1) = 2.0;
Container.RHS(2) = 3.0;
Container.RHS(3) = 4.0;
Container.RHS(4) = 5.0;
// The linear system with the submatrix is solved as follows:
Container.ApplyInverse().

A call to Compute() computes the LU factorization of the linear system matrix, using LAPACK (more precisely, by calling the corresponding routines in Ifpack_SerialTriDiSolver). The default behavior is to store the matrix factors by overwriting the linear system matrix itself. This way, method Apply() fails, as the original matrix does no longer exists. An alternative is to call KeepNonFactoredMatrix(true), which forces Ifpack_TriDiContainer to maintain in memory a copy of the non-factored matrix.

Author
Marzio Sala, SNL 9214.
Date
Last update Nov-04.

Definition at line 118 of file Ifpack_TriDiContainer.h.

Constructor & Destructor Documentation

Ifpack_TriDiContainer::Ifpack_TriDiContainer ( const int  NumRows_in,
const int  NumVectors_in = 1 
)
inline

Default constructor.

Definition at line 125 of file Ifpack_TriDiContainer.h.

Ifpack_TriDiContainer::Ifpack_TriDiContainer ( const Ifpack_TriDiContainer rhs)
inline

Copy constructor.

Definition at line 137 of file Ifpack_TriDiContainer.h.

virtual Ifpack_TriDiContainer::~Ifpack_TriDiContainer ( )
inlinevirtual

Destructor.

Definition at line 153 of file Ifpack_TriDiContainer.h.

Member Function Documentation

Ifpack_TriDiContainer& Ifpack_TriDiContainer::operator= ( const Ifpack_TriDiContainer rhs)
inline

Operator=.

Definition at line 160 of file Ifpack_TriDiContainer.h.

int Ifpack_TriDiContainer::NumRows ( ) const
virtual

Returns the number of rows of the matrix and LHS/RHS.

Implements Ifpack_Container.

Definition at line 48 of file Ifpack_TriDiContainer.cpp.

virtual int Ifpack_TriDiContainer::NumVectors ( ) const
inlinevirtual

Returns the number of vectors in LHS/RHS.

Implements Ifpack_Container.

Definition at line 187 of file Ifpack_TriDiContainer.h.

virtual int Ifpack_TriDiContainer::SetNumVectors ( const int  NumVectors_in)
inlinevirtual

Sets the number of vectors for LHS/RHS.

Implements Ifpack_Container.

Definition at line 193 of file Ifpack_TriDiContainer.h.

double & Ifpack_TriDiContainer::LHS ( const int  i,
const int  Vector = 0 
)
virtual

Returns the i-th component of the vector Vector of LHS.

Implements Ifpack_Container.

Definition at line 95 of file Ifpack_TriDiContainer.cpp.

double & Ifpack_TriDiContainer::RHS ( const int  i,
const int  Vector = 0 
)
virtual

Returns the i-th component of the vector Vector of RHS.

Implements Ifpack_Container.

Definition at line 101 of file Ifpack_TriDiContainer.cpp.

int & Ifpack_TriDiContainer::ID ( const int  i)
virtual

Returns the ID associated to local row i.

The set of (local) rows assigned to this container is defined by calling ID(i) = j, where i (from 0 to NumRows()) indicates the container-row, and j indicates the local row in the calling process.

This is usually used to recorder the local row ID (on calling process) of the i-th row in the container.

Implements Ifpack_Container.

Definition at line 145 of file Ifpack_TriDiContainer.cpp.

int Ifpack_TriDiContainer::SetMatrixElement ( const int  row,
const int  col,
const double  value 
)
virtual

Set the matrix element (row,col) to value.

Implements Ifpack_Container.

Definition at line 108 of file Ifpack_TriDiContainer.cpp.

virtual int Ifpack_TriDiContainer::SetParameters ( Teuchos::ParameterList )
inlinevirtual

Sets all necessary parameters.

Implements Ifpack_Container.

Definition at line 237 of file Ifpack_TriDiContainer.h.

virtual bool Ifpack_TriDiContainer::IsInitialized ( ) const
inlinevirtual

Returns true is the container has been successfully initialized.

Implements Ifpack_Container.

Definition at line 243 of file Ifpack_TriDiContainer.h.

virtual bool Ifpack_TriDiContainer::IsComputed ( ) const
inlinevirtual

Returns true is the container has been successfully computed.

Implements Ifpack_Container.

Definition at line 249 of file Ifpack_TriDiContainer.h.

virtual const char* Ifpack_TriDiContainer::Label ( ) const
inlinevirtual

Returns the label of this container.

Implements Ifpack_Container.

Definition at line 255 of file Ifpack_TriDiContainer.h.

virtual int Ifpack_TriDiContainer::SetKeepNonFactoredMatrix ( const bool  flag)
inlinevirtual

If flag is true, keeps a copy of the non-factored matrix.

Definition at line 261 of file Ifpack_TriDiContainer.h.

virtual bool Ifpack_TriDiContainer::KeepNonFactoredMatrix ( ) const
inlinevirtual

Returns KeepNonFactoredMatrix_.

Definition at line 268 of file Ifpack_TriDiContainer.h.

virtual const Epetra_SerialDenseMatrix& Ifpack_TriDiContainer::LHS ( ) const
inlinevirtual

Returns the dense vector containing the LHS.

Definition at line 274 of file Ifpack_TriDiContainer.h.

virtual const Epetra_SerialDenseMatrix& Ifpack_TriDiContainer::RHS ( ) const
inlinevirtual

Returns the dense vector containing the RHS.

Definition at line 280 of file Ifpack_TriDiContainer.h.

virtual const Ifpack_SerialTriDiMatrix& Ifpack_TriDiContainer::Matrix ( ) const
inlinevirtual

Returns the dense matrix or its factors.

Definition at line 286 of file Ifpack_TriDiContainer.h.

virtual const Ifpack_SerialTriDiMatrix& Ifpack_TriDiContainer::NonFactoredMatrix ( ) const
inlinevirtual

Returns the non-factored dense matrix (only if stored).

Definition at line 292 of file Ifpack_TriDiContainer.h.

virtual const Epetra_IntSerialDenseVector& Ifpack_TriDiContainer::ID ( ) const
inlinevirtual

Returns the integer dense vector of IDs.

Definition at line 298 of file Ifpack_TriDiContainer.h.

int Ifpack_TriDiContainer::Initialize ( )
virtual

Initialize the container.

Implements Ifpack_Container.

Definition at line 54 of file Ifpack_TriDiContainer.cpp.

int Ifpack_TriDiContainer::Compute ( const Epetra_RowMatrix Matrix_in)
virtual

Finalizes the linear system matrix and prepares for the application of the inverse.

Implements Ifpack_Container.

Definition at line 208 of file Ifpack_TriDiContainer.cpp.

int Ifpack_TriDiContainer::Apply ( )
virtual

Apply the matrix to RHS, results are stored in LHS.

Implements Ifpack_Container.

Definition at line 240 of file Ifpack_TriDiContainer.cpp.

int Ifpack_TriDiContainer::ApplyInverse ( )
virtual

Apply the inverse of the matrix to RHS, results are stored in LHS.

Implements Ifpack_Container.

Definition at line 128 of file Ifpack_TriDiContainer.cpp.

virtual double Ifpack_TriDiContainer::InitializeFlops ( ) const
inlinevirtual

Returns the flops in Initialize().

Implements Ifpack_Container.

Definition at line 320 of file Ifpack_TriDiContainer.h.

virtual double Ifpack_TriDiContainer::ComputeFlops ( ) const
inlinevirtual

Returns the flops in Compute().

Implements Ifpack_Container.

Definition at line 325 of file Ifpack_TriDiContainer.h.

virtual double Ifpack_TriDiContainer::ApplyFlops ( ) const
inlinevirtual

Returns the flops in Apply().

Implements Ifpack_Container.

Definition at line 330 of file Ifpack_TriDiContainer.h.

virtual double Ifpack_TriDiContainer::ApplyInverseFlops ( ) const
inlinevirtual

Returns the flops in ApplyInverse().

Implements Ifpack_Container.

Definition at line 335 of file Ifpack_TriDiContainer.h.

std::ostream & Ifpack_TriDiContainer::Print ( std::ostream &  os) const
virtual

Prints basic information on iostream. This function is used by operator<<.

Implements Ifpack_Container.

Definition at line 259 of file Ifpack_TriDiContainer.cpp.

int Ifpack_TriDiContainer::Extract ( const Epetra_RowMatrix Matrix_in)
privatevirtual

Extract the submatrices identified by the ID set int ID().

Definition at line 152 of file Ifpack_TriDiContainer.cpp.

Member Data Documentation

int Ifpack_TriDiContainer::NumRows_
private

Number of rows in the container.

Definition at line 349 of file Ifpack_TriDiContainer.h.

int Ifpack_TriDiContainer::NumVectors_
private

Number of vectors in the container.

Definition at line 351 of file Ifpack_TriDiContainer.h.

Ifpack_SerialTriDiMatrix Ifpack_TriDiContainer::NonFactoredMatrix_
private

TriDi matrix, that contains the non-factored matrix.

Definition at line 353 of file Ifpack_TriDiContainer.h.

Ifpack_SerialTriDiMatrix Ifpack_TriDiContainer::Matrix_
private

TriDi matrix.

Definition at line 355 of file Ifpack_TriDiContainer.h.

Epetra_SerialDenseMatrix Ifpack_TriDiContainer::LHS_
private

SerialDense vector representing the LHS.

Definition at line 357 of file Ifpack_TriDiContainer.h.

Epetra_SerialDenseMatrix Ifpack_TriDiContainer::RHS_
private

SerialDense vector representing the RHS.

Definition at line 359 of file Ifpack_TriDiContainer.h.

Ifpack_SerialTriDiSolver Ifpack_TriDiContainer::Solver_
private

TriDi solver (solution will be get using LAPACK).

Definition at line 361 of file Ifpack_TriDiContainer.h.

Epetra_IntSerialDenseVector Ifpack_TriDiContainer::ID_
private

Sets of local rows.

Definition at line 363 of file Ifpack_TriDiContainer.h.

bool Ifpack_TriDiContainer::KeepNonFactoredMatrix_
private

If true, keeps a copy of the non-factored matrix.

Definition at line 365 of file Ifpack_TriDiContainer.h.

bool Ifpack_TriDiContainer::IsInitialized_
private

If true, the container has been successfully initialized.

Definition at line 367 of file Ifpack_TriDiContainer.h.

bool Ifpack_TriDiContainer::IsComputed_
private

If true, the container has been successfully computed.

Definition at line 369 of file Ifpack_TriDiContainer.h.

std::string Ifpack_TriDiContainer::Label_
private

Label for this object.

Definition at line 371 of file Ifpack_TriDiContainer.h.

double Ifpack_TriDiContainer::ComputeFlops_
private

Flops in Compute().

Definition at line 374 of file Ifpack_TriDiContainer.h.

double Ifpack_TriDiContainer::ApplyFlops_
private

Flops in Apply().

Definition at line 376 of file Ifpack_TriDiContainer.h.

double Ifpack_TriDiContainer::ApplyInverseFlops_
private

Flops in ApplyInverse().

Definition at line 378 of file Ifpack_TriDiContainer.h.


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