38#ifndef IFPACK_SUBDOMAINFILTER_H
39#define IFPACK_SUBDOMAINFILTER_H
41#if defined(Ifpack_SHOW_DEPRECATED_WARNINGS)
43#warning "The Ifpack package is deprecated"
47#include "Ifpack_ConfigDefs.h"
49#ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
52#include "Epetra_MpiComm.h"
54#include "Epetra_SerialComm.h"
56#include "Epetra_RowMatrix.h"
57#include "Epetra_CrsMatrix.h"
58#include "Epetra_IntVector.h"
59#include "Teuchos_RCP.hpp"
60#include "Ifpack_OverlappingRowMatrix.h"
89 Ifpack_SubdomainFilter(
const Teuchos::RCP<const Epetra_RowMatrix>& Matrix,
int subdomainId);
94 ~Ifpack_SubdomainFilter();
110 NumEntries = NumEntries_[MyRow];
117 return(MaxNumEntries_);
135 virtual inline int ExtractMyRowCopy(
int MyRow,
int Length,
int & NumEntries,
double *Values,
int * Indices)
const;
162 if (TransA ==
true) {
166 IFPACK_CHK_ERR(
Apply(X,Y));
212 virtual bool Filled()
const
235#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
239 return(NumGlobalNonzeros_);
246 return(NumGlobalRows_);
252 return(NumGlobalCols_);
258 return(NumGlobalRows_);
263 virtual long long NumGlobalNonzeros64()
const
265 return(NumGlobalNonzeros_);
269 virtual long long NumGlobalRows64()
const
271 return(NumGlobalRows_);
275 virtual long long NumGlobalCols64()
const
277 return(NumGlobalRows_);
281 virtual long long NumGlobalDiagonals64()
const
283 return(NumGlobalRows_);
289 return(NumMyNonzeros_);
313 return(Matrix_->LowerTriangular());
319 return(Matrix_->UpperTriangular());
341 virtual const Epetra_Import* Importer()
const {
return(&*Importer_);}
343 virtual const Epetra_Export* Exporter()
const {
return(&*Exporter_);}
348 int SetOwnership(
bool ownership)
356 UseTranspose_ = UseTranspose_in;
363 return(UseTranspose_);
393const char*
Label()
const{
398 void UpdateImportVector(
int NumVectors)
const;
399 void UpdateExportVector(
int NumVectors)
const;
402 Teuchos::RCP<const Epetra_RowMatrix> Matrix_;
406 Teuchos::RCP<Epetra_MpiComm> SubComm_;
407 MPI_Comm subdomainMPIComm_;
410 Teuchos::RCP<Epetra_SerialComm> SubComm_;
413 Teuchos::RCP<Epetra_Map> Map_;
415 Teuchos::RCP<Epetra_Map> colMap_;
427 int NumGlobalNonzeros_;
433 std::vector<int> NumEntries_;
435 mutable std::vector<int> Indices_;
437 mutable std::vector<double> Values_;
442 Teuchos::RCP<Epetra_Vector> Diagonal_;
463 Teuchos::RCP<Epetra_Import> Importer_;
464 Teuchos::RCP<Epetra_Export> Exporter_;
virtual int SetUseTranspose(bool UseTranspose)=0
virtual int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const=0
virtual const Epetra_Comm & Comm() const=0
virtual const char * Label() const=0
virtual int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const=0
virtual const Epetra_Map & OperatorDomainMap() const=0
virtual bool HasNormInf() const=0
virtual const Epetra_Map & OperatorRangeMap() const=0
virtual bool UseTranspose() const=0
virtual int NumMyRows() const=0
virtual int NumMyCols() const=0
virtual int Multiply(bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const=0
virtual int ExtractDiagonalCopy(Epetra_Vector &Diagonal) const=0
virtual int NumGlobalCols() const=0
virtual int NumGlobalNonzeros() const=0
virtual const Epetra_Map & RowMatrixColMap() const=0
virtual const Epetra_Map & RowMatrixRowMap() const=0
virtual int NumMyNonzeros() const=0
virtual int NumMyDiagonals() const=0
virtual int NumMyRowEntries(int MyRow, int &NumEntries) const=0
virtual int Solve(bool Upper, bool Trans, bool UnitDiagonal, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const=0
virtual int NumGlobalRows() const=0
virtual bool LowerTriangular() const=0
virtual int NumGlobalDiagonals() const=0
virtual int MaxNumEntries() const=0
virtual int ExtractMyRowCopy(int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const=0
virtual int InvRowSums(Epetra_Vector &x) const=0
virtual const Epetra_Import * RowMatrixImporter() const=0
virtual bool UpperTriangular() const=0
virtual int RightScale(const Epetra_Vector &x)=0
virtual double NormOne() const=0
virtual double NormInf() const=0
virtual bool Filled() const=0
virtual int InvColSums(Epetra_Vector &x) const=0
virtual int LeftScale(const Epetra_Vector &x)=0
virtual const Epetra_BlockMap & Map() const=0
Ifpack_OverlappingRowMatrix: matrix with ghost rows, based on Epetra_RowMatrix.