44 #ifndef IFPACK_ADDITIVESCHWARZ_H
45 #define IFPACK_ADDITIVESCHWARZ_H
47 #if defined(Ifpack_SHOW_DEPRECATED_WARNINGS)
49 #warning "The Ifpack package is deprecated"
66 #include "Epetra_CombineMode.h"
67 #include "Epetra_MultiVector.h"
68 #include "Epetra_Map.h"
69 #include "Epetra_Comm.h"
70 #include "Epetra_Time.h"
71 #include "Epetra_LinearProblem.h"
72 #include "Epetra_RowMatrix.h"
73 #include "Epetra_CrsMatrix.h"
74 #include "Teuchos_ParameterList.hpp"
75 #include "Teuchos_RefCountPtr.hpp"
77 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
79 #include "EpetraExt_Reindex_CrsMatrix.h"
80 #include "EpetraExt_Reindex_MultiVector.h"
82 #ifdef IFPACK_NODE_AWARE_CODE
83 #include "EpetraExt_OperatorOut.h"
84 #include "EpetraExt_RowMatrixOut.h"
85 #include "EpetraExt_BlockMapOut.h"
88 #ifdef HAVE_IFPACK_AMESOS
166 int OverlapLevel_in = 0);
217 virtual double NormInf()
const;
223 virtual const char *
Label()
const;
287 const int MaxIters = 1550,
288 const double Tol = 1e-9,
310 virtual std::ostream&
Print(std::ostream&)
const;
399 Teuchos::RefCountPtr<const Epetra_RowMatrix>
Matrix_;
408 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
427 # ifdef IFPACK_NODE_AWARE_CODE
433 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
440 int NumMpiProcsPerSubdomain_;
497 Teuchos::RefCountPtr<Epetra_Time>
Time_;
501 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
502 mutable Teuchos::RefCountPtr<Epetra_MultiVector> OverlappingX;
503 mutable Teuchos::RefCountPtr<Epetra_MultiVector> OverlappingY;
505 # ifdef IFPACK_NODE_AWARE_CODE
506 mutable Teuchos::RefCountPtr<Epetra_MultiVector> OverlappingX;
507 mutable Teuchos::RefCountPtr<Epetra_MultiVector> OverlappingY;
516 int OverlapLevel_in) :
517 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
520 NumMpiProcsPerSubdomain_(1),
524 IsInitialized_(
false),
526 UseTranspose_(
false),
527 IsOverlapping_(
false),
528 OverlapLevel_(OverlapLevel_in),
531 ComputeCondest_(
true),
532 UseReordering_(
false),
533 ReorderingType_(
"none"),
534 FilterSingletons_(
false),
538 InitializeTime_(0.0),
540 ApplyInverseTime_(0.0),
541 InitializeFlops_(0.0),
543 ApplyInverseFlops_(0.0)
548 if (
Matrix_->Comm().NumProc() == 1)
558 #ifdef IFPACK_NODE_AWARE_CODE
559 extern int ML_NODE_ID;
571 #ifndef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
572 # ifdef IFPACK_NODE_AWARE_CODE
590 try{ nodeID = List_.get(
"ML node id",0);}
591 catch(...){fprintf(stderr,
"%s",
"Ifpack_AdditiveSchwarz<T>::Setup(): no parameter \"ML node id\"\n\n");
592 std::cout << List_ << std::endl;}
599 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
600 if (NumMpiProcsPerSubdomain_ > 1) {
601 LocalizedMatrix_ =
Teuchos::rcp(
new Ifpack_SubdomainFilter(OverlappingMatrix_, SubdomainId_));
606 # ifdef IFPACK_NODE_AWARE_CODE
607 Ifpack_NodeFilter *tt =
new Ifpack_NodeFilter(OverlappingMatrix_,nodeID);
617 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
619 if (NumMpiProcsPerSubdomain_ > 1) {
620 LocalizedMatrix_ =
Teuchos::rcp(
new Ifpack_SubdomainFilter(Matrix_, SubdomainId_));
625 # ifdef IFPACK_NODE_AWARE_CODE
626 Ifpack_NodeFilter *tt =
new Ifpack_NodeFilter(Matrix_,nodeID);
636 fprintf(stderr,
"%s",
"AdditiveSchwarz Setup: problem creating local filter matrix.\n");
643 if (FilterSingletons_) {
645 MatrixPtr = &*SingletonFilter_;
648 MatrixPtr = &*LocalizedMatrix_;
650 if (UseReordering_) {
653 if (ReorderingType_ ==
"rcm")
655 else if (ReorderingType_ ==
"metis")
657 #ifdef HAVE_IFPACK_AMESOS
658 else if (ReorderingType_ ==
"amd" )
662 cerr <<
"reordering type not correct (" << ReorderingType_ <<
")" << endl;
671 ReorderedLocalizedMatrix_ =
676 MatrixPtr = &*ReorderedLocalizedMatrix_;
679 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
687 SubdomainCrsMatrix_->Import(*MatrixPtr, *tempImporter,
Insert);
688 SubdomainCrsMatrix_->FillComplete();
690 if (NumMpiProcsPerSubdomain_ > 1) {
691 tempMap_.reset(
new Epetra_Map(SubdomainCrsMatrix_->RowMap().NumGlobalElements(),
692 SubdomainCrsMatrix_->RowMap().NumMyElements(),
693 0, SubdomainCrsMatrix_->Comm()));
694 tempRangeMap_.reset(
new Epetra_Map(SubdomainCrsMatrix_->OperatorRangeMap().NumGlobalElements(),
695 SubdomainCrsMatrix_->OperatorRangeMap().NumMyElements(),
696 0, SubdomainCrsMatrix_->Comm()));
697 tempDomainMap_.reset(
new Epetra_Map(SubdomainCrsMatrix_->OperatorDomainMap().NumGlobalElements(),
698 SubdomainCrsMatrix_->OperatorDomainMap().NumMyElements(),
699 0, SubdomainCrsMatrix_->Comm()));
705 ReindexedCrsMatrix_.reset(&((*SubdomainMatrixReindexer_)(*SubdomainCrsMatrix_)),
false);
707 MatrixPtr = &*ReindexedCrsMatrix_;
709 Inverse_ =
Teuchos::rcp(
new T(&*ReindexedCrsMatrix_) );
711 Inverse_ =
Teuchos::rcp(
new T(&*SubdomainCrsMatrix_) );
727 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
728 MpiRank_ = Matrix_->Comm().MyPID();
729 NumMpiProcs_ = Matrix_->Comm().NumProc();
730 NumMpiProcsPerSubdomain_ = List_in.
get(
"subdomain: number-of-processors", 1);
731 NumSubdomains_ = NumMpiProcs_ / NumMpiProcsPerSubdomain_;
732 SubdomainId_ = MpiRank_ / NumMpiProcsPerSubdomain_;
734 if (NumSubdomains_ == 1) {
735 IsOverlapping_ =
false;
740 ComputeCondest_ = List_in.
get(
"schwarz: compute condest", ComputeCondest_);
744 if(
typeid(std::string) == combineModeEntry->getAny().type() )
746 std::string mode = List_in.
get(
"schwarz: combine mode",
"Add");
749 else if (mode ==
"Zero")
751 else if (mode ==
"Insert")
753 else if (mode ==
"InsertAdd")
755 else if (mode ==
"Average")
757 else if (mode ==
"AbsMax")
762 true,std::logic_error
763 ,
"Error, The (Epetra) combine mode of \""<<mode<<
"\" is not valid! Only the values"
764 " \"Add\", \"Zero\", \"Insert\", \"InsertAdd\", \"Average\", and \"AbsMax\" are accepted!"
775 Teuchos::getParameter<std::string>(List_in,
"schwarz: combine mode");
781 List_in.
get(
"schwarz: combine mode",
"Zero");
784 ReorderingType_ = List_in.
get(
"schwarz: reordering type", ReorderingType_);
785 if (ReorderingType_ ==
"none")
786 UseReordering_ =
false;
788 UseReordering_ =
true;
793 FilterSingletons_ = List_in.
get(
"schwarz: filter singletons", FilterSingletons_);
805 IsInitialized_ =
false;
812 Time_->ResetStartTime();
815 if (IsOverlapping_) {
816 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
817 if (NumMpiProcsPerSubdomain_ > 1) {
823 # ifdef IFPACK_NODE_AWARE_CODE
825 try{ myNodeID = List_.get(
"ML node id",-1);}
826 catch(...){fprintf(stderr,
"pid %d: no such entry (returned %d)\n",Comm().MyPID(),myNodeID);}
843 # ifdef IFPACK_NODE_AWARE_CODE
852 # ifdef IFPACK_NODE_AWARE_CODE
870 Label_ =
"Ifpack_AdditiveSchwarz, ";
872 Label_ +=
", transp";
874 +
", local solver = \n\t\t***** `" + std::string(Inverse_->Label()) +
"'";
876 IsInitialized_ =
true;
878 InitializeTime_ += Time_->ElapsedTime();
880 #ifdef IFPACK_FLOPCOUNTERS
884 double partial = Inverse_->InitializeFlops();
886 Comm().SumAll(&partial, &total, 1);
887 InitializeFlops_ += total;
897 if (IsInitialized() ==
false)
900 Time_->ResetStartTime();
908 ComputeTime_ += Time_->ElapsedTime();
910 #ifdef IFPACK_FLOPCOUNTERS
912 double partial = Inverse_->ComputeFlops();
914 Comm().SumAll(&partial, &total, 1);
915 ComputeFlops_ += total;
921 R = ReorderingType_ +
" reord, ";
927 Label_ =
"Ifpack_AdditiveSchwarz, ov = " +
Ifpack_toString(OverlapLevel_)
928 +
", local solver = \n\t\t***** `" + std::string(Inverse_->Label()) +
"'"
929 +
"\n\t\t***** " + R +
"Condition number estimate = "
941 UseTranspose_ = UseTranspose_in;
969 return(Label_.c_str());
976 return(UseTranspose_);
990 return(Matrix_->Comm());
997 return(Matrix_->OperatorDomainMap());
1001 template<
typename T>
1004 return(Matrix_->OperatorRangeMap());
1008 template<
typename T>
1018 if (NumVectors != Y.NumVectors())
1021 Time_->ResetStartTime();
1023 Teuchos::RefCountPtr<Epetra_MultiVector> OverlappingX;
1024 Teuchos::RefCountPtr<Epetra_MultiVector> OverlappingY;
1025 Teuchos::RefCountPtr<Epetra_MultiVector> Xtmp;
1028 #ifdef IFPACK_FLOPCOUNTERS
1029 double pre_partial = Inverse_->ApplyInverseFlops();
1031 Comm().SumAll(&pre_partial, &pre_total, 1);
1035 if (IsOverlapping()) {
1036 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
1040 }
else assert(OverlappingX->NumVectors() == X.NumVectors());
1044 }
else assert(OverlappingY->NumVectors() == Y.NumVectors());
1046 # ifdef IFPACK_NODE_AWARE_CODE
1051 }
else assert(OverlappingX->NumVectors() == X.NumVectors());
1056 }
else assert(OverlappingY->NumVectors() == Y.NumVectors());
1065 OverlappingY->PutScalar(0.0);
1066 OverlappingX->PutScalar(0.0);
1074 OverlappingX = Xtmp;
1078 if (FilterSingletons_) {
1082 IFPACK_CHK_ERR(SingletonFilter_->SolveSingletons(*OverlappingX,*OverlappingY));
1083 IFPACK_CHK_ERR(SingletonFilter_->CreateReducedRHS(*OverlappingY,*OverlappingX,ReducedX));
1086 if (!UseReordering_) {
1098 IFPACK_CHK_ERR(SingletonFilter_->UpdateLHS(ReducedY,*OverlappingY));
1102 if (!UseReordering_) {
1103 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
1104 if (NumMpiProcsPerSubdomain_ > 1) {
1105 tempX_.reset(&((*RangeVectorReindexer_)(*OverlappingX)),
false);
1106 tempY_.reset(&((*DomainVectorReindexer_)(*OverlappingY)),
false);
1109 IFPACK_CHK_ERR(Inverse_->ApplyInverse(*OverlappingX, *OverlappingY));
1112 IFPACK_CHK_ERR(Inverse_->ApplyInverse(*OverlappingX,*OverlappingY));
1124 if (IsOverlapping()) {
1125 IFPACK_CHK_ERR(OverlappingMatrix_->ExportMultiVector(*OverlappingY,Y,
1129 #ifdef IFPACK_FLOPCOUNTERS
1132 double partial = Inverse_->ApplyInverseFlops();
1134 Comm().SumAll(&partial, &total, 1);
1135 ApplyInverseFlops_ += total - pre_total;
1140 ApplyInverseTime_ += Time_->ElapsedTime();
1147 template<
typename T>
1153 #ifdef IFPACK_FLOPCOUNTERS
1154 double IF = InitializeFlops();
1155 double CF = ComputeFlops();
1156 double AF = ApplyInverseFlops();
1158 double IFT = 0.0, CFT = 0.0, AFT = 0.0;
1159 if (InitializeTime() != 0.0)
1160 IFT = IF / InitializeTime();
1161 if (ComputeTime() != 0.0)
1162 CFT = CF / ComputeTime();
1163 if (ApplyInverseTime() != 0.0)
1164 AFT = AF / ApplyInverseTime();
1167 if (Matrix().Comm().MyPID())
1171 os <<
"================================================================================" << endl;
1172 os <<
"Ifpack_AdditiveSchwarz, overlap level = " << OverlapLevel_ << endl;
1173 if (CombineMode_ ==
Insert)
1174 os <<
"Combine mode = Insert" << endl;
1175 else if (CombineMode_ ==
Add)
1176 os <<
"Combine mode = Add" << endl;
1177 else if (CombineMode_ ==
Zero)
1178 os <<
"Combine mode = Zero" << endl;
1179 else if (CombineMode_ ==
Average)
1180 os <<
"Combine mode = Average" << endl;
1181 else if (CombineMode_ ==
AbsMax)
1182 os <<
"Combine mode = AbsMax" << endl;
1184 os <<
"Condition number estimate = " << Condest_ << endl;
1185 os <<
"Global number of rows = " << Matrix_->NumGlobalRows64() << endl;
1187 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
1189 os <<
"================================================================================" << endl;
1190 os <<
"Subcommunicator stats" << endl;
1191 os <<
"Number of MPI processes in simulation: " << NumMpiProcs_ << endl;
1192 os <<
"Number of subdomains: " << NumSubdomains_ << endl;
1193 os <<
"Number of MPI processes per subdomain: " << NumMpiProcsPerSubdomain_ << endl;
1197 os <<
"Phase # calls Total Time (s) Total MFlops MFlops/s" << endl;
1198 os <<
"----- ------- -------------- ------------ --------" << endl;
1199 os <<
"Initialize() " << std::setw(5) << NumInitialize()
1200 <<
" " << std::setw(15) << InitializeTime()
1201 #ifdef IFPACK_FLOPCOUNTERS
1202 <<
" " << std::setw(15) << 1.0e-6 * IF
1203 <<
" " << std::setw(15) << 1.0e-6 * IFT
1206 os <<
"Compute() " << std::setw(5) << NumCompute()
1207 <<
" " << std::setw(15) << ComputeTime()
1208 #ifdef IFPACK_FLOPCOUNTERS
1209 <<
" " << std::setw(15) << 1.0e-6 * CF
1210 <<
" " << std::setw(15) << 1.0e-6 * CFT
1213 os <<
"ApplyInverse() " << std::setw(5) << NumApplyInverse()
1214 <<
" " << std::setw(15) << ApplyInverseTime()
1215 #ifdef IFPACK_FLOPCOUNTERS
1216 <<
" " << std::setw(15) << 1.0e-6 * AF
1217 <<
" " << std::setw(15) << 1.0e-6 * AFT
1220 os <<
"================================================================================" << endl;
1228 template<
typename T>
1241 #endif // IFPACK_ADDITIVESCHWARZ_H
virtual int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Applies the preconditioner to X, returns the result in Y.
int Setup()
Sets up the localized matrix and the singleton filter.
int OverlapLevel_
Level of overlap among the processors.
std::string ReorderingType_
Type of reordering of the local matrix.
Teuchos::RefCountPtr< const Epetra_RowMatrix > Matrix_
Pointers to the matrix to be preconditioned.
virtual const Epetra_Map & RowMatrixRowMap() const =0
double ComputeFlops_
Contains the number of flops for Compute().
virtual int SetUseTranspose(bool UseTranspose_in)
If set true, transpose of this operator will be applied (not implemented).
virtual const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
Ifpack_AMDReordering: approximate minimum degree reordering.
T & get(ParameterList &l, const std::string &name)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
bool UseTranspose_
If true, solve with the transpose (not supported by all solvers).
#define IFPACK_CHK_GLOBAL_ERR(ifpack_err)
virtual double ComputeFlops() const
Returns the number of flops in the computation phase.
bool ComputeCondest_
If true, compute the condition number estimate each time Compute() is called.
double Condest_
Contains the estimated condition number.
virtual bool IsInitialized() const
Returns true if the preconditioner has been successfully initialized.
Ifpack_METISReordering: A class to reorder a graph using METIS.
double ApplyInverseFlops_
Contain sthe number of flops for ApplyInverse().
virtual bool IsComputed() const
Returns true if the preconditioner has been successfully computed.
Teuchos::RefCountPtr< Ifpack_OverlappingRowMatrix > OverlappingMatrix_
Pointers to the overlapping matrix.
virtual int OverlapLevel() const
Returns the level of overlap.
virtual double ApplyInverseTime() const
Returns the time spent in ApplyInverse().
virtual bool IsOverlapping() const
Returns true is an overlapping matrix is present.
virtual double ComputeTime() const
Returns the time spent in Compute().
Teuchos::RefCountPtr< Ifpack_Reordering > Reordering_
Pointer to a reordering object.
virtual double ApplyInverseFlops() const
Returns the number of flops in the application of the preconditioner.
Teuchos::RefCountPtr< Ifpack_ReorderFilter > ReorderedLocalizedMatrix_
Pointer to the reorderd matrix.
ParameterEntry * getEntryPtr(const std::string &name)
Ifpack_CondestType
Ifpack_CondestType: enum to define the type of condition number estimate.
std::string Label_
Contains the label of this object.
virtual int NumInitialize() const
Returns the number of calls to Initialize().
bool IsOverlapping_
If true, overlapping is used.
virtual std::ostream & Print(std::ostream &) const
Prints major information about this preconditioner.
virtual int Initialize()
Initialized the preconditioner.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Ifpack_AdditiveSchwarz(const Ifpack_AdditiveSchwarz &RHS)
Copy constructor (should never be used)
double ComputeTime_
Contains the time for all successful calls to Compute().
Ifpack_RCMReordering: reverse Cuthill-McKee reordering.
Ifpack_AdditiveSchwarz: a class to define Additive Schwarz preconditioners of Epetra_RowMatrix's.
virtual int NumApplyInverse() const
Returns the number of calls to ApplyInverse().
virtual const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
Ifpack_Preconditioner: basic class for preconditioning in Ifpack.
virtual const Epetra_BlockMap & Map() const =0
virtual bool UseTranspose() const
Returns the current UseTranspose setting.
Ifpack_AdditiveSchwarz(Epetra_RowMatrix *Matrix_in, int OverlapLevel_in=0)
Ifpack_AdditiveSchwarz constructor with given Epetra_RowMatrix.
double InitializeFlops_
Contains the number of flops for Initialize().
virtual bool HasNormInf() const
Returns true if the this object can provide an approximate Inf-norm, false otherwise.
Ifpack_SingletonFilter: Filter based on matrix entries.
virtual int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Applies the matrix to X, returns the result in Y.
virtual ~Ifpack_AdditiveSchwarz()
Destructor.
Ifpack_OverlappingRowMatrix: matrix with ghost rows, based on Epetra_RowMatrix.
Ifpack_LocalFilter a class for light-weight extraction of the submatrix corresponding to local rows a...
bool IsInitialized_
If true, the preconditioner has been successfully initialized.
std::string Ifpack_toString(const int &x)
Converts an integer to std::string.
int NumInitialize_
Contains the number of successful calls to Initialize().
virtual int SetParameters(Teuchos::ParameterList &List)
Sets the parameters.
bool UseReordering_
If true, reorder the local matrix.
double Ifpack_Condest(const Ifpack_Preconditioner &IFP, const Ifpack_CondestType CT, const int MaxIters, const double Tol, Epetra_RowMatrix *Matrix)
virtual int Compute()
Computes the preconditioner.
int NumApplyInverse_
Contains the number of successful call to ApplyInverse().
Teuchos::RefCountPtr< Ifpack_SingletonFilter > SingletonFilter_
filtering object.
double InitializeTime_
Contains the time for all successful calls to Initialize().
virtual double InitializeFlops() const
Returns the number of flops in the initialization phase.
bool IsComputed_
If true, the preconditioner has been successfully computed.
virtual const Teuchos::ParameterList & List() const
Returns a reference to the internally stored list.
virtual const char * Label() const
Returns a character string describing the operator.
Teuchos::RefCountPtr< T > Inverse_
Pointer to the local solver.
int NumCompute_
Contains the number of successful call to Compute().
virtual const T * Inverse() const
virtual double InitializeTime() const
Returns the time spent in Initialize().
virtual const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this operator.
double ApplyInverseTime_
Contains the time for all successful calls to ApplyInverse().
#define IFPACK_CHK_ERR(ifpack_err)
Ifpack_ReorderFilter: a class for light-weight reorder of local rows and columns of an Epetra_RowMatr...
virtual double NormInf() const
Returns the infinity norm of the global matrix (not implemented)
Teuchos::ParameterList List_
Stores a copy of the list given in SetParameters()
Epetra_CombineMode CombineMode_
Combine mode for off-process elements (only if overlap is used)
virtual const Epetra_RowMatrix & Matrix() const
Returns a refernence to the internally stored matrix.
bool FilterSingletons_
Filter for singletons.
virtual double Condest() const
Returns the estimated condition number, or -1.0 if not computed.
virtual int NumCompute() const
Returns the number of calls to Compute().
Teuchos::RefCountPtr< Epetra_Time > Time_
Object used for timing purposes.
Teuchos::RefCountPtr< Ifpack_LocalFilter > LocalizedMatrix_
Localized version of Matrix_ or OverlappingMatrix_.