IFPACK  Development
 All Classes Namespaces Files Functions Variables Enumerations Friends Pages
Public Member Functions | List of all members
Ifpack_SparseContainer< T > Class Template Reference

Ifpack_SparseContainer: a class for storing and solving linear systems using sparse matrices. More...

#include <Ifpack_SparseContainer.h>

Inheritance diagram for Ifpack_SparseContainer< T >:
Inheritance graph
[legend]
Collaboration diagram for Ifpack_SparseContainer< T >:
Collaboration graph
[legend]

Public Member Functions

virtual double InitializeFlops () const
 Returns the flops in Compute().
 
virtual double ComputeFlops () const
 Returns the flops in Compute().
 
virtual double ApplyFlops () const
 Returns the flops in Apply().
 
virtual double ApplyInverseFlops () const
 Returns the flops in ApplyInverse().
 
virtual std::ostream & Print (std::ostream &os) const
 Prints basic information on iostream. This function is used by operator<<.
 
 Ifpack_SparseContainer (const int NumRows, const int NumVectors=1)
 Constructor.
 
 Ifpack_SparseContainer (const Ifpack_SparseContainer< T > &rhs)
 Copy constructor.
 
virtual ~Ifpack_SparseContainer ()
 Destructor.
 
Ifpack_SparseContaineroperator= (const Ifpack_SparseContainer< T > &rhs)
 Operator =.
 
virtual int NumRows () const
 Returns the number of rows of the matrix and LHS/RHS.
 
virtual int NumVectors () const
 Returns the number of vectors in LHS/RHS.
 
virtual int SetNumVectors (const int NumVectors_in)
 Sets the number of vectors for LHS/RHS.
 
virtual double & LHS (const int i, const int Vector=0)
 Returns the i-th component of the vector Vector of LHS.
 
virtual double & RHS (const int i, const int Vector=0)
 Returns the i-th component of the vector Vector of RHS.
 
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.
 
virtual bool IsInitialized () const
 Returns true is the container has been successfully initialized.
 
virtual bool IsComputed () const
 Returns true is the container has been successfully computed.
 
virtual int SetParameters (Teuchos::ParameterList &List)
 Sets all necessary parameters.
 
virtual const char * Label () const
 Returns the label of this container.
 
Teuchos::RCP< const Epetra_MapMap () const
 Returns a pointer to the internally stored map.
 
Teuchos::RCP< const
Epetra_MultiVector
LHS () const
 Returns a pointer to the internally stored solution multi-vector.
 
Teuchos::RCP< const
Epetra_MultiVector
RHS () const
 Returns a pointer to the internally stored rhs multi-vector.
 
Teuchos::RCP< const
Epetra_CrsMatrix
Matrix () const
 Returns a pointer to the internally stored matrix.
 
const Epetra_IntSerialDenseVectorID () const
 Returns a pointer to the internally stored ID's.
 
Teuchos::RCP< const T > Inverse () const
 Returns a pointer to the internally stored inverse operator.
 
virtual int Initialize ()
 Initializes the container, by completing all the operations based on matrix structure. More...
 
virtual int Compute (const Epetra_RowMatrix &Matrix_in)
 Finalizes the linear system matrix and prepares for the application of the inverse.
 
virtual int Apply ()
 Apply the matrix to RHS, result is stored in LHS.
 
virtual int ApplyInverse ()
 Apply the inverse of the matrix to RHS, result is stored in LHS.
 
virtual int Destroy ()
 Destroys all data.
 
- Public Member Functions inherited from Ifpack_Container
virtual ~Ifpack_Container ()
 Destructor.
 

Detailed Description

template<typename T>
class Ifpack_SparseContainer< T >

Ifpack_SparseContainer: a class for storing and solving linear systems using sparse 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_DenseContainer, one can store the blocks as sparse matrices (Epetra_CrsMatrix), which can be advantageous when the blocks are large. Otherwise, class Ifpack_DenseContainer is probably more appropriate.

Sparse containers are templated with a type T, which represent the class to use in the application of the inverse. (T is not used in Ifpack_DenseContainer). In SparseContainer, T must be an Ifpack_Preconditioner derived class. The container will allocate a T object, use SetParameters() and Compute(), then use T every time the linear system as to be solved (using the ApplyInverse() method of T).

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

Definition at line 93 of file Ifpack_SparseContainer.h.

Member Function Documentation

template<typename T >
int & Ifpack_SparseContainer< T >::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 542 of file Ifpack_SparseContainer.h.

template<typename T >
int Ifpack_SparseContainer< T >::Initialize ( )
virtual

Initializes the container, by completing all the operations based on matrix structure.

Note
After a call to Initialize(), no new matrix entries can be added.

Implements Ifpack_Container.

Definition at line 378 of file Ifpack_SparseContainer.h.

References Copy.


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