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

Ifpack_LocalFilter a class for light-weight extraction of the submatrix corresponding to local rows and columns. More...

#include <Ifpack_LocalFilter.h>

Inheritance diagram for Ifpack_LocalFilter:
Inheritance graph
[legend]

Public Member Functions

int SetOwnership (bool)
 Sets ownership. More...
 
int SetUseTranspose (bool UseTranspose_in)
 Sets use transpose (not implemented). More...
 
bool UseTranspose () const
 Returns the current UseTranspose setting. More...
 
bool HasNormInf () const
 Returns true if the this object can provide an approximate Inf-norm, false otherwise. More...
 
const Epetra_CommComm () const
 Returns a pointer to the Epetra_Comm communicator associated with this operator. More...
 
const Epetra_MapOperatorDomainMap () const
 Returns the Epetra_Map object associated with the domain of this operator. More...
 
const Epetra_MapOperatorRangeMap () const
 Returns the Epetra_Map object associated with the range of this operator. More...
 
const Epetra_BlockMapMap () const
 
const char * Label () const
 

Private Attributes

Teuchos::RefCountPtr< const
Epetra_RowMatrix
Matrix_
 Pointer to the matrix to be preconditioned. More...
 
Teuchos::RefCountPtr
< Epetra_SerialComm
SerialComm_
 Communicator containing this process only. More...
 
Teuchos::RefCountPtr< Epetra_MapMap_
 Map based on SerialComm_, containing the local rows only. More...
 
int NumRows_
 Number of rows in the local matrix. More...
 
int NumNonzeros_
 Number of nonzeros in the local matrix. More...
 
int MaxNumEntries_
 Maximum number of nonzero entries in a row for the filtered matrix. More...
 
int MaxNumEntriesA_
 Maximum number of nonzero entries in a row for Matrix_. More...
 
std::vector< int > NumEntries_
 NumEntries_[i] contains the nonzero entries in row `i'. More...
 
std::vector< int > Indices_
 Used in ExtractMyRowCopy, to avoid allocation each time. More...
 
std::vector< double > Values_
 Used in ExtractMyRowCopy, to avoid allocation each time. More...
 
bool UseTranspose_
 If true, the tranpose of the local matrix will be used. More...
 
char Label_ [80]
 Label for this object. More...
 
Teuchos::RefCountPtr
< Epetra_Vector
Diagonal_
 
 Ifpack_LocalFilter (const Teuchos::RefCountPtr< const Epetra_RowMatrix > &Matrix)
 Constructor. More...
 
virtual ~Ifpack_LocalFilter ()
 Destructor. More...
 
virtual int NumMyRowEntries (int MyRow, int &NumEntries) const
 Returns the number of nonzero entries in MyRow. More...
 
virtual int MaxNumEntries () const
 Returns the maximum of NumMyRowEntries() over all rows. More...
 
virtual int ExtractMyRowCopy (int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const
 Returns a copy of the specified local row in user-provided arrays. More...
 
virtual int ExtractDiagonalCopy (Epetra_Vector &Diagonal) const
 Returns a copy of the main diagonal in a user-provided vector. More...
 
virtual int Multiply (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_RowMatrix multiplied by a Epetra_MultiVector X in Y. More...
 
virtual int Solve (bool, bool, bool, const Epetra_MultiVector &, Epetra_MultiVector &) const
 Returns result of a local-only solve using a triangular Epetra_RowMatrix with Epetra_MultiVectors X and Y (NOT IMPLEMENTED). More...
 
virtual int Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 
virtual int ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 
virtual int InvRowSums (Epetra_Vector &) const
 Computes the sum of absolute values of the rows of the Epetra_RowMatrix, results returned in x (NOT IMPLEMENTED). More...
 
virtual int LeftScale (const Epetra_Vector &)
 Scales the Epetra_RowMatrix on the left with a Epetra_Vector x (NOT IMPLEMENTED). More...
 
virtual int InvColSums (Epetra_Vector &) const
 Computes the sum of absolute values of the columns of the Epetra_RowMatrix, results returned in x (NOT IMPLEMENTED). More...
 
virtual int RightScale (const Epetra_Vector &)
 Scales the Epetra_RowMatrix on the right with a Epetra_Vector x (NOT IMPLEMENTED). More...
 
virtual bool Filled () const
 If FillComplete() has been called, this query returns true, otherwise it returns false. More...
 
virtual double NormInf () const
 Returns the infinity norm of the global matrix. More...
 
virtual double NormOne () const
 Returns the one norm of the global matrix. More...
 
virtual int NumGlobalNonzeros () const
 Returns the number of nonzero entries in the global matrix. More...
 
virtual int NumGlobalRows () const
 Returns the number of global matrix rows. More...
 
virtual int NumGlobalCols () const
 Returns the number of global matrix columns. More...
 
virtual int NumGlobalDiagonals () const
 Returns the number of global nonzero diagonal entries, based on global row/column index comparisons. More...
 
virtual long long NumGlobalNonzeros64 () const
 Returns the number of nonzero entries in the global matrix. More...
 
virtual long long NumGlobalRows64 () const
 Returns the number of global matrix rows. More...
 
virtual long long NumGlobalCols64 () const
 Returns the number of global matrix columns. More...
 
virtual long long NumGlobalDiagonals64 () const
 Returns the number of global nonzero diagonal entries, based on global row/column index comparisons. More...
 
virtual int NumMyNonzeros () const
 Returns the number of nonzero entries in the calling processor's portion of the matrix. More...
 
virtual int NumMyRows () const
 Returns the number of matrix rows owned by the calling processor. More...
 
virtual int NumMyCols () const
 Returns the number of matrix columns owned by the calling processor. More...
 
virtual int NumMyDiagonals () const
 Returns the number of local nonzero diagonal entries, based on global row/column index comparisons. More...
 
virtual bool LowerTriangular () const
 If matrix is lower triangular in local index space, this query returns true, otherwise it returns false. More...
 
virtual bool UpperTriangular () const
 If matrix is upper triangular in local index space, this query returns true, otherwise it returns false. More...
 
virtual const Epetra_MapRowMatrixRowMap () const
 Returns the Epetra_Map object associated with the rows of this matrix. More...
 
virtual const Epetra_MapRowMatrixColMap () const
 Returns the Epetra_Map object associated with the columns of this matrix. More...
 
virtual const Epetra_ImportRowMatrixImporter () const
 Returns the Epetra_Import object that contains the import operations for distributed operations. More...
 

Detailed Description

Ifpack_LocalFilter a class for light-weight extraction of the submatrix corresponding to local rows and columns.

Class Ifpack_LocalFilter enables a light-weight contruction of an Epetra_RowMatrix-derived object, containing only the elements of the original, distributed matrix with local row and column ID. The local submatrix is based on a communicator containing the local process only. Each process will have its local object, corresponding to the local submatrix. Submatrices may or may not overlap.

The following instructions can be used to create "localized" matrices:

...
Teuchos::RefCountPtr<Epetra_RowMatrix> A; // fill the elements of A,
A->FillComplete();

Once created, LocalA defined, on each process, the submatrix corresponding to local rows and columns only. The creation and use of LocalA is "cheap", as the elements of the local matrix are obtained through calls to ExtractMyRowCopy on the original, distributed matrix, say A. This means that A must remain in scope every time LocalA is accessed.

A very convenient use of this class is to use Ifpack solvers to compute the LU factorizations of local blocks. If applied to a localized matrix, several Ifpack objects can operator in the same phase in a safe way, without non-required data exchange.

Author
Marzio Sala, SNL 9214
Date
Sep-04

Definition at line 98 of file Ifpack_LocalFilter.h.

Constructor & Destructor Documentation

Ifpack_LocalFilter::Ifpack_LocalFilter ( const Teuchos::RefCountPtr< const Epetra_RowMatrix > &  Matrix)

Constructor.

Definition at line 53 of file Ifpack_LocalFilter.cpp.

virtual Ifpack_LocalFilter::~Ifpack_LocalFilter ( )
inlinevirtual

Destructor.

Definition at line 108 of file Ifpack_LocalFilter.h.

Member Function Documentation

virtual int Ifpack_LocalFilter::NumMyRowEntries ( int  MyRow,
int &  NumEntries 
) const
inlinevirtual

Returns the number of nonzero entries in MyRow.

Parameters
MyRow- (In) Local row.
NumEntries- (Out) Number of nonzero values present.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_RowMatrix.

Definition at line 123 of file Ifpack_LocalFilter.h.

virtual int Ifpack_LocalFilter::MaxNumEntries ( ) const
inlinevirtual

Returns the maximum of NumMyRowEntries() over all rows.

Implements Epetra_RowMatrix.

Definition at line 130 of file Ifpack_LocalFilter.h.

int Ifpack_LocalFilter::ExtractMyRowCopy ( int  MyRow,
int  Length,
int &  NumEntries,
double *  Values,
int *  Indices 
) const
virtual

Returns a copy of the specified local row in user-provided arrays.

Parameters
MyRow- (In) Local row to extract.
Length- (In) Length of Values and Indices.
NumEntries- (Out) Number of nonzero entries extracted.
Values- (Out) Extracted values for this row.
Indices- (Out) Extracted global column indices for the corresponding values.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_RowMatrix.

Definition at line 131 of file Ifpack_LocalFilter.cpp.

int Ifpack_LocalFilter::ExtractDiagonalCopy ( Epetra_Vector Diagonal) const
virtual

Returns a copy of the main diagonal in a user-provided vector.

Parameters
Diagonal- (Out) Extracted main diagonal.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_RowMatrix.

Definition at line 167 of file Ifpack_LocalFilter.cpp.

virtual int Ifpack_LocalFilter::Multiply ( bool  TransA,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
inlinevirtual

Returns the result of a Epetra_RowMatrix multiplied by a Epetra_MultiVector X in Y.

Parameters
TransA-(In) If true, multiply by the transpose of matrix, otherwise just use matrix.
X- (In) A Epetra_MultiVector of dimension NumVectors to multiply with matrix.
Y-(Out) A Epetra_MultiVector of dimension NumVectorscontaining result.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_RowMatrix.

Definition at line 175 of file Ifpack_LocalFilter.h.

virtual int Ifpack_LocalFilter::Solve ( bool  ,
bool  ,
bool  ,
const Epetra_MultiVector ,
Epetra_MultiVector  
) const
inlinevirtual

Returns result of a local-only solve using a triangular Epetra_RowMatrix with Epetra_MultiVectors X and Y (NOT IMPLEMENTED).

Implements Epetra_RowMatrix.

Definition at line 186 of file Ifpack_LocalFilter.h.

int Ifpack_LocalFilter::Apply ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Implements Epetra_Operator.

Definition at line 176 of file Ifpack_LocalFilter.cpp.

int Ifpack_LocalFilter::ApplyInverse ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Implements Epetra_Operator.

Definition at line 209 of file Ifpack_LocalFilter.cpp.

virtual int Ifpack_LocalFilter::InvRowSums ( Epetra_Vector ) const
inlinevirtual

Computes the sum of absolute values of the rows of the Epetra_RowMatrix, results returned in x (NOT IMPLEMENTED).

Implements Epetra_RowMatrix.

Definition at line 198 of file Ifpack_LocalFilter.h.

virtual int Ifpack_LocalFilter::LeftScale ( const Epetra_Vector )
inlinevirtual

Scales the Epetra_RowMatrix on the left with a Epetra_Vector x (NOT IMPLEMENTED).

Implements Epetra_RowMatrix.

Definition at line 204 of file Ifpack_LocalFilter.h.

virtual int Ifpack_LocalFilter::InvColSums ( Epetra_Vector ) const
inlinevirtual

Computes the sum of absolute values of the columns of the Epetra_RowMatrix, results returned in x (NOT IMPLEMENTED).

Implements Epetra_RowMatrix.

Definition at line 210 of file Ifpack_LocalFilter.h.

virtual int Ifpack_LocalFilter::RightScale ( const Epetra_Vector )
inlinevirtual

Scales the Epetra_RowMatrix on the right with a Epetra_Vector x (NOT IMPLEMENTED).

Implements Epetra_RowMatrix.

Definition at line 217 of file Ifpack_LocalFilter.h.

virtual bool Ifpack_LocalFilter::Filled ( ) const
inlinevirtual

If FillComplete() has been called, this query returns true, otherwise it returns false.

Implements Epetra_RowMatrix.

Definition at line 227 of file Ifpack_LocalFilter.h.

virtual double Ifpack_LocalFilter::NormInf ( ) const
inlinevirtual

Returns the infinity norm of the global matrix.

Implements Epetra_RowMatrix.

Definition at line 236 of file Ifpack_LocalFilter.h.

virtual double Ifpack_LocalFilter::NormOne ( ) const
inlinevirtual

Returns the one norm of the global matrix.

Implements Epetra_RowMatrix.

Definition at line 245 of file Ifpack_LocalFilter.h.

virtual int Ifpack_LocalFilter::NumGlobalNonzeros ( ) const
inlinevirtual

Returns the number of nonzero entries in the global matrix.

Implements Epetra_RowMatrix.

Definition at line 252 of file Ifpack_LocalFilter.h.

virtual int Ifpack_LocalFilter::NumGlobalRows ( ) const
inlinevirtual

Returns the number of global matrix rows.

Implements Epetra_RowMatrix.

Definition at line 258 of file Ifpack_LocalFilter.h.

virtual int Ifpack_LocalFilter::NumGlobalCols ( ) const
inlinevirtual

Returns the number of global matrix columns.

Implements Epetra_RowMatrix.

Definition at line 264 of file Ifpack_LocalFilter.h.

virtual int Ifpack_LocalFilter::NumGlobalDiagonals ( ) const
inlinevirtual

Returns the number of global nonzero diagonal entries, based on global row/column index comparisons.

Implements Epetra_RowMatrix.

Definition at line 270 of file Ifpack_LocalFilter.h.

virtual long long Ifpack_LocalFilter::NumGlobalNonzeros64 ( ) const
inlinevirtual

Returns the number of nonzero entries in the global matrix.

Implements Epetra_RowMatrix.

Definition at line 277 of file Ifpack_LocalFilter.h.

virtual long long Ifpack_LocalFilter::NumGlobalRows64 ( ) const
inlinevirtual

Returns the number of global matrix rows.

Implements Epetra_RowMatrix.

Definition at line 283 of file Ifpack_LocalFilter.h.

virtual long long Ifpack_LocalFilter::NumGlobalCols64 ( ) const
inlinevirtual

Returns the number of global matrix columns.

Implements Epetra_RowMatrix.

Definition at line 289 of file Ifpack_LocalFilter.h.

virtual long long Ifpack_LocalFilter::NumGlobalDiagonals64 ( ) const
inlinevirtual

Returns the number of global nonzero diagonal entries, based on global row/column index comparisons.

Implements Epetra_RowMatrix.

Definition at line 295 of file Ifpack_LocalFilter.h.

virtual int Ifpack_LocalFilter::NumMyNonzeros ( ) const
inlinevirtual

Returns the number of nonzero entries in the calling processor's portion of the matrix.

Implements Epetra_RowMatrix.

Definition at line 301 of file Ifpack_LocalFilter.h.

virtual int Ifpack_LocalFilter::NumMyRows ( ) const
inlinevirtual

Returns the number of matrix rows owned by the calling processor.

Implements Epetra_RowMatrix.

Definition at line 307 of file Ifpack_LocalFilter.h.

virtual int Ifpack_LocalFilter::NumMyCols ( ) const
inlinevirtual

Returns the number of matrix columns owned by the calling processor.

Implements Epetra_RowMatrix.

Definition at line 313 of file Ifpack_LocalFilter.h.

virtual int Ifpack_LocalFilter::NumMyDiagonals ( ) const
inlinevirtual

Returns the number of local nonzero diagonal entries, based on global row/column index comparisons.

Implements Epetra_RowMatrix.

Definition at line 319 of file Ifpack_LocalFilter.h.

virtual bool Ifpack_LocalFilter::LowerTriangular ( ) const
inlinevirtual

If matrix is lower triangular in local index space, this query returns true, otherwise it returns false.

Implements Epetra_RowMatrix.

Definition at line 325 of file Ifpack_LocalFilter.h.

virtual bool Ifpack_LocalFilter::UpperTriangular ( ) const
inlinevirtual

If matrix is upper triangular in local index space, this query returns true, otherwise it returns false.

Implements Epetra_RowMatrix.

Definition at line 331 of file Ifpack_LocalFilter.h.

virtual const Epetra_Map& Ifpack_LocalFilter::RowMatrixRowMap ( ) const
inlinevirtual

Returns the Epetra_Map object associated with the rows of this matrix.

Implements Epetra_RowMatrix.

Definition at line 337 of file Ifpack_LocalFilter.h.

virtual const Epetra_Map& Ifpack_LocalFilter::RowMatrixColMap ( ) const
inlinevirtual

Returns the Epetra_Map object associated with the columns of this matrix.

Implements Epetra_RowMatrix.

Definition at line 343 of file Ifpack_LocalFilter.h.

virtual const Epetra_Import* Ifpack_LocalFilter::RowMatrixImporter ( ) const
inlinevirtual

Returns the Epetra_Import object that contains the import operations for distributed operations.

Implements Epetra_RowMatrix.

Definition at line 349 of file Ifpack_LocalFilter.h.

int Ifpack_LocalFilter::SetOwnership ( bool  )
inline

Sets ownership.

Definition at line 358 of file Ifpack_LocalFilter.h.

int Ifpack_LocalFilter::SetUseTranspose ( bool  UseTranspose_in)
inlinevirtual

Sets use transpose (not implemented).

Implements Epetra_Operator.

Definition at line 364 of file Ifpack_LocalFilter.h.

bool Ifpack_LocalFilter::UseTranspose ( ) const
inlinevirtual

Returns the current UseTranspose setting.

Implements Epetra_Operator.

Definition at line 371 of file Ifpack_LocalFilter.h.

bool Ifpack_LocalFilter::HasNormInf ( ) const
inlinevirtual

Returns true if the this object can provide an approximate Inf-norm, false otherwise.

Implements Epetra_Operator.

Definition at line 377 of file Ifpack_LocalFilter.h.

const Epetra_Comm& Ifpack_LocalFilter::Comm ( ) const
inlinevirtual

Returns a pointer to the Epetra_Comm communicator associated with this operator.

Implements Epetra_Operator.

Definition at line 383 of file Ifpack_LocalFilter.h.

const Epetra_Map& Ifpack_LocalFilter::OperatorDomainMap ( ) const
inlinevirtual

Returns the Epetra_Map object associated with the domain of this operator.

Implements Epetra_Operator.

Definition at line 389 of file Ifpack_LocalFilter.h.

const Epetra_Map& Ifpack_LocalFilter::OperatorRangeMap ( ) const
inlinevirtual

Returns the Epetra_Map object associated with the range of this operator.

Implements Epetra_Operator.

Definition at line 395 of file Ifpack_LocalFilter.h.

const Epetra_BlockMap & Ifpack_LocalFilter::Map ( ) const
virtual

Implements Epetra_SrcDistObject.

Definition at line 216 of file Ifpack_LocalFilter.cpp.

const char* Ifpack_LocalFilter::Label ( ) const
inlinevirtual

Implements Epetra_Operator.

Definition at line 403 of file Ifpack_LocalFilter.h.

Member Data Documentation

Teuchos::RefCountPtr<const Epetra_RowMatrix> Ifpack_LocalFilter::Matrix_
private

Pointer to the matrix to be preconditioned.

Definition at line 405 of file Ifpack_LocalFilter.h.

Teuchos::RefCountPtr<Epetra_SerialComm> Ifpack_LocalFilter::SerialComm_
private

Communicator containing this process only.

Definition at line 416 of file Ifpack_LocalFilter.h.

Teuchos::RefCountPtr<Epetra_Map> Ifpack_LocalFilter::Map_
private

Map based on SerialComm_, containing the local rows only.

Definition at line 419 of file Ifpack_LocalFilter.h.

int Ifpack_LocalFilter::NumRows_
private

Number of rows in the local matrix.

Definition at line 421 of file Ifpack_LocalFilter.h.

int Ifpack_LocalFilter::NumNonzeros_
private

Number of nonzeros in the local matrix.

Definition at line 423 of file Ifpack_LocalFilter.h.

int Ifpack_LocalFilter::MaxNumEntries_
private

Maximum number of nonzero entries in a row for the filtered matrix.

Definition at line 425 of file Ifpack_LocalFilter.h.

int Ifpack_LocalFilter::MaxNumEntriesA_
private

Maximum number of nonzero entries in a row for Matrix_.

Definition at line 427 of file Ifpack_LocalFilter.h.

std::vector<int> Ifpack_LocalFilter::NumEntries_
private

NumEntries_[i] contains the nonzero entries in row `i'.

Definition at line 429 of file Ifpack_LocalFilter.h.

std::vector<int> Ifpack_LocalFilter::Indices_
mutableprivate

Used in ExtractMyRowCopy, to avoid allocation each time.

Definition at line 431 of file Ifpack_LocalFilter.h.

std::vector<double> Ifpack_LocalFilter::Values_
mutableprivate

Used in ExtractMyRowCopy, to avoid allocation each time.

Definition at line 433 of file Ifpack_LocalFilter.h.

bool Ifpack_LocalFilter::UseTranspose_
private

If true, the tranpose of the local matrix will be used.

Definition at line 435 of file Ifpack_LocalFilter.h.

char Ifpack_LocalFilter::Label_[80]
private

Label for this object.

Definition at line 437 of file Ifpack_LocalFilter.h.

Teuchos::RefCountPtr<Epetra_Vector> Ifpack_LocalFilter::Diagonal_
private

Definition at line 438 of file Ifpack_LocalFilter.h.


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