Ifpack Package Browser (Single Doxygen Collection)
Development
|
Ifpack_DenseContainer: a class to define containers for dense matrices. More...
#include <Ifpack_DenseContainer.h>
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... | |
Epetra_SerialDenseMatrix | NonFactoredMatrix_ |
Dense matrix, that contains the non-factored matrix. More... | |
Epetra_SerialDenseMatrix | Matrix_ |
Dense matrix. More... | |
Epetra_SerialDenseMatrix | LHS_ |
Dense vector representing the LHS. More... | |
Epetra_SerialDenseMatrix | RHS_ |
Dense vector representing the RHS. More... | |
Epetra_SerialDenseSolver | Solver_ |
Dense 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_DenseContainer (const int NumRows_in, const int NumVectors_in=1) | |
Default constructor. More... | |
Ifpack_DenseContainer (const Ifpack_DenseContainer &rhs) | |
Copy constructor. More... | |
virtual | ~Ifpack_DenseContainer () |
Destructor. More... | |
Ifpack_DenseContainer & | operator= (const Ifpack_DenseContainer &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 Epetra_SerialDenseMatrix & | Matrix () const |
Returns the dense matrix or its factors. More... | |
virtual const Epetra_SerialDenseMatrix & | 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... | |
Ifpack_DenseContainer: 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_DenseContainer, 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:
A call to Compute() computes the LU factorization of the linear system matrix, using LAPACK (more precisely, by calling the corresponding routines in Epetra_SerialDenseSolver). 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_DenseContainer to maintain in memory a copy of the non-factored matrix.
Definition at line 117 of file Ifpack_DenseContainer.h.
|
inline |
Default constructor.
Definition at line 124 of file Ifpack_DenseContainer.h.
|
inline |
Copy constructor.
Definition at line 136 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
Destructor.
Definition at line 152 of file Ifpack_DenseContainer.h.
|
inline |
Operator=.
Definition at line 159 of file Ifpack_DenseContainer.h.
|
virtual |
Returns the number of rows of the matrix and LHS/RHS.
Implements Ifpack_Container.
Definition at line 48 of file Ifpack_DenseContainer.cpp.
|
inlinevirtual |
Returns the number of vectors in LHS/RHS.
Implements Ifpack_Container.
Definition at line 186 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
Sets the number of vectors for LHS/RHS.
Implements Ifpack_Container.
Definition at line 192 of file Ifpack_DenseContainer.h.
|
virtual |
Returns the i-th component of the vector Vector of LHS.
Implements Ifpack_Container.
Definition at line 91 of file Ifpack_DenseContainer.cpp.
|
virtual |
Returns the i-th component of the vector Vector of RHS.
Implements Ifpack_Container.
Definition at line 97 of file Ifpack_DenseContainer.cpp.
|
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 141 of file Ifpack_DenseContainer.cpp.
|
virtual |
Set the matrix element (row,col) to value
.
Implements Ifpack_Container.
Definition at line 104 of file Ifpack_DenseContainer.cpp.
|
inlinevirtual |
Sets all necessary parameters.
Implements Ifpack_Container.
Definition at line 236 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
Returns true
is the container has been successfully initialized.
Implements Ifpack_Container.
Definition at line 242 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
Returns true
is the container has been successfully computed.
Implements Ifpack_Container.
Definition at line 248 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
Returns the label of this container.
Implements Ifpack_Container.
Definition at line 254 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
If flag
is true
, keeps a copy of the non-factored matrix.
Definition at line 260 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
Returns KeepNonFactoredMatrix_.
Definition at line 267 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
Returns the dense vector containing the LHS.
Definition at line 273 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
Returns the dense vector containing the RHS.
Definition at line 279 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
Returns the dense matrix or its factors.
Definition at line 285 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
Returns the non-factored dense matrix (only if stored).
Definition at line 291 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
Returns the integer dense vector of IDs.
Definition at line 297 of file Ifpack_DenseContainer.h.
|
virtual |
Initialize the container.
Implements Ifpack_Container.
Definition at line 54 of file Ifpack_DenseContainer.cpp.
|
virtual |
Finalizes the linear system matrix and prepares for the application of the inverse.
Implements Ifpack_Container.
Definition at line 204 of file Ifpack_DenseContainer.cpp.
|
virtual |
Apply the matrix to RHS, results are stored in LHS.
Implements Ifpack_Container.
Definition at line 236 of file Ifpack_DenseContainer.cpp.
|
virtual |
Apply the inverse of the matrix to RHS, results are stored in LHS.
Implements Ifpack_Container.
Definition at line 124 of file Ifpack_DenseContainer.cpp.
|
inlinevirtual |
Returns the flops in Initialize().
Implements Ifpack_Container.
Definition at line 319 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
Returns the flops in Compute().
Implements Ifpack_Container.
Definition at line 324 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
Returns the flops in Apply().
Implements Ifpack_Container.
Definition at line 329 of file Ifpack_DenseContainer.h.
|
inlinevirtual |
Returns the flops in ApplyInverse().
Implements Ifpack_Container.
Definition at line 334 of file Ifpack_DenseContainer.h.
|
virtual |
Prints basic information on iostream. This function is used by operator<<.
Implements Ifpack_Container.
Definition at line 254 of file Ifpack_DenseContainer.cpp.
|
privatevirtual |
Extract the submatrices identified by the ID set int ID().
Definition at line 148 of file Ifpack_DenseContainer.cpp.
|
private |
Number of rows in the container.
Definition at line 348 of file Ifpack_DenseContainer.h.
|
private |
Number of vectors in the container.
Definition at line 350 of file Ifpack_DenseContainer.h.
|
private |
Dense matrix, that contains the non-factored matrix.
Definition at line 352 of file Ifpack_DenseContainer.h.
|
private |
Dense matrix.
Definition at line 354 of file Ifpack_DenseContainer.h.
|
private |
Dense vector representing the LHS.
Definition at line 356 of file Ifpack_DenseContainer.h.
|
private |
Dense vector representing the RHS.
Definition at line 358 of file Ifpack_DenseContainer.h.
|
private |
Dense solver (solution will be get using LAPACK).
Definition at line 360 of file Ifpack_DenseContainer.h.
|
private |
Sets of local rows.
Definition at line 362 of file Ifpack_DenseContainer.h.
|
private |
If true
, keeps a copy of the non-factored matrix.
Definition at line 364 of file Ifpack_DenseContainer.h.
|
private |
If true
, the container has been successfully initialized.
Definition at line 366 of file Ifpack_DenseContainer.h.
|
private |
If true
, the container has been successfully computed.
Definition at line 368 of file Ifpack_DenseContainer.h.
|
private |
Label for this
object.
Definition at line 370 of file Ifpack_DenseContainer.h.
|
private |
Flops in Compute().
Definition at line 373 of file Ifpack_DenseContainer.h.
|
private |
Flops in Apply().
Definition at line 375 of file Ifpack_DenseContainer.h.
|
private |
Flops in ApplyInverse().
Definition at line 377 of file Ifpack_DenseContainer.h.