44 #ifndef IFPACK_ADDITIVESCHWARZ_H
45 #define IFPACK_ADDITIVESCHWARZ_H
47 #include "Ifpack_ConfigDefs.h"
48 #include "Ifpack_Preconditioner.h"
49 #include "Ifpack_ConfigDefs.h"
50 #include "Ifpack_Preconditioner.h"
51 #include "Ifpack_Reordering.h"
52 #include "Ifpack_RCMReordering.h"
53 #include "Ifpack_METISReordering.h"
54 #include "Ifpack_LocalFilter.h"
55 #include "Ifpack_NodeFilter.h"
56 #include "Ifpack_SingletonFilter.h"
57 #include "Ifpack_ReorderFilter.h"
59 #include "Ifpack_OverlappingRowMatrix.h"
61 #include "Epetra_MultiVector.h"
62 #include "Epetra_Map.h"
63 #include "Epetra_Comm.h"
64 #include "Epetra_Time.h"
65 #include "Epetra_LinearProblem.h"
66 #include "Epetra_RowMatrix.h"
67 #include "Epetra_CrsMatrix.h"
68 #include "Teuchos_ParameterList.hpp"
69 #include "Teuchos_RefCountPtr.hpp"
71 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
72 #include "Ifpack_SubdomainFilter.h"
73 #include "EpetraExt_Reindex_CrsMatrix.h"
74 #include "EpetraExt_Reindex_MultiVector.h"
76 #ifdef IFPACK_NODE_AWARE_CODE
77 #include "EpetraExt_OperatorOut.h"
78 #include "EpetraExt_RowMatrixOut.h"
79 #include "EpetraExt_BlockMapOut.h"
82 #ifdef HAVE_IFPACK_AMESOS
83 #include "Ifpack_AMDReordering.h"
160 int OverlapLevel_in = 0);
211 virtual double NormInf()
const;
217 virtual const char *
Label()
const;
280 virtual double Condest(
const Ifpack_CondestType CT = Ifpack_Cheap,
281 const int MaxIters = 1550,
282 const double Tol = 1e-9,
304 virtual std::ostream&
Print(std::ostream&)
const;
306 virtual const T* Inverse()
const
370 virtual const Teuchos::ParameterList&
List()
const
393 Teuchos::RefCountPtr<const Epetra_RowMatrix>
Matrix_;
402 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
406 Teuchos::RCP<Epetra_CrsMatrix> SubdomainCrsMatrix_;
408 Teuchos::RCP<Epetra_CrsMatrix> ReindexedCrsMatrix_;
412 Teuchos::RCP<Epetra_Map> tempMap_;
413 Teuchos::RCP<Epetra_Map> tempDomainMap_;
414 Teuchos::RCP<Epetra_Map> tempRangeMap_;
415 Teuchos::RCP<EpetraExt::CrsMatrix_Reindex> SubdomainMatrixReindexer_;
416 Teuchos::RCP<EpetraExt::MultiVector_Reindex> DomainVectorReindexer_;
417 Teuchos::RCP<EpetraExt::MultiVector_Reindex> RangeVectorReindexer_;
418 mutable Teuchos::RCP<Epetra_MultiVector> tempX_;
419 mutable Teuchos::RCP<Epetra_MultiVector> tempY_;
421 # ifdef IFPACK_NODE_AWARE_CODE
427 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
434 int NumMpiProcsPerSubdomain_;
491 Teuchos::RefCountPtr<Epetra_Time>
Time_;
495 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
496 mutable Teuchos::RefCountPtr<Epetra_MultiVector> OverlappingX;
497 mutable Teuchos::RefCountPtr<Epetra_MultiVector> OverlappingY;
499 # ifdef IFPACK_NODE_AWARE_CODE
500 mutable Teuchos::RefCountPtr<Epetra_MultiVector> OverlappingX;
501 mutable Teuchos::RefCountPtr<Epetra_MultiVector> OverlappingY;
510 int OverlapLevel_in) :
511 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
514 NumMpiProcsPerSubdomain_(1),
518 IsInitialized_(false),
520 UseTranspose_(false),
521 IsOverlapping_(false),
522 OverlapLevel_(OverlapLevel_in),
525 ComputeCondest_(true),
526 UseReordering_(false),
527 ReorderingType_(
"none"),
528 FilterSingletons_(false),
532 InitializeTime_(0.0),
534 ApplyInverseTime_(0.0),
535 InitializeFlops_(0.0),
537 ApplyInverseFlops_(0.0)
540 Matrix_ = Teuchos::rcp( Matrix_in,
false );
542 if (
Matrix_->Comm().NumProc() == 1)
548 Teuchos::ParameterList List_in;
552 #ifdef IFPACK_NODE_AWARE_CODE
553 extern int ML_NODE_ID;
565 #ifndef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
566 # ifdef IFPACK_NODE_AWARE_CODE
584 try{ nodeID = List_.get(
"ML node id",0);}
585 catch(...){fprintf(stderr,
"%s",
"Ifpack_AdditiveSchwarz<T>::Setup(): no parameter \"ML node id\"\n\n");
586 std::cout << List_ << std::endl;}
591 if (OverlappingMatrix_ != Teuchos::null)
593 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
594 if (NumMpiProcsPerSubdomain_ > 1) {
595 LocalizedMatrix_ = Teuchos::rcp(
new Ifpack_SubdomainFilter(OverlappingMatrix_, SubdomainId_));
600 # ifdef IFPACK_NODE_AWARE_CODE
601 Ifpack_NodeFilter *tt =
new Ifpack_NodeFilter(OverlappingMatrix_,nodeID);
602 LocalizedMatrix_ = Teuchos::rcp(tt);
611 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
613 if (NumMpiProcsPerSubdomain_ > 1) {
614 LocalizedMatrix_ = Teuchos::rcp(
new Ifpack_SubdomainFilter(Matrix_, SubdomainId_));
619 # ifdef IFPACK_NODE_AWARE_CODE
620 Ifpack_NodeFilter *tt =
new Ifpack_NodeFilter(Matrix_,nodeID);
621 LocalizedMatrix_ = Teuchos::rcp(tt);
630 fprintf(stderr,
"%s",
"AdditiveSchwarz Setup: problem creating local filter matrix.\n");
633 if (LocalizedMatrix_ == Teuchos::null)
637 if (FilterSingletons_) {
639 MatrixPtr = &*SingletonFilter_;
642 MatrixPtr = &*LocalizedMatrix_;
644 if (UseReordering_) {
647 if (ReorderingType_ ==
"rcm")
649 else if (ReorderingType_ ==
"metis")
651 #ifdef HAVE_IFPACK_AMESOS
652 else if (ReorderingType_ ==
"amd" )
656 cerr <<
"reordering type not correct (" << ReorderingType_ <<
")" << endl;
659 if (Reordering_ == Teuchos::null) IFPACK_CHK_ERR(-5);
661 IFPACK_CHK_ERR(Reordering_->SetParameters(List_));
662 IFPACK_CHK_ERR(Reordering_->Compute(*MatrixPtr));
665 ReorderedLocalizedMatrix_ =
668 if (ReorderedLocalizedMatrix_ == Teuchos::null) IFPACK_CHK_ERR(-5);
670 MatrixPtr = &*ReorderedLocalizedMatrix_;
673 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
680 Teuchos::RCP<Epetra_Import> tempImporter = Teuchos::rcp(
new Epetra_Import(SubdomainCrsMatrix_->Map(), MatrixPtr->
Map()));
681 SubdomainCrsMatrix_->Import(*MatrixPtr, *tempImporter,
Insert);
682 SubdomainCrsMatrix_->FillComplete();
684 if (NumMpiProcsPerSubdomain_ > 1) {
685 tempMap_.reset(
new Epetra_Map(SubdomainCrsMatrix_->RowMap().NumGlobalElements(),
686 SubdomainCrsMatrix_->RowMap().NumMyElements(),
687 0, SubdomainCrsMatrix_->Comm()));
688 tempRangeMap_.reset(
new Epetra_Map(SubdomainCrsMatrix_->OperatorRangeMap().NumGlobalElements(),
689 SubdomainCrsMatrix_->OperatorRangeMap().NumMyElements(),
690 0, SubdomainCrsMatrix_->Comm()));
691 tempDomainMap_.reset(
new Epetra_Map(SubdomainCrsMatrix_->OperatorDomainMap().NumGlobalElements(),
692 SubdomainCrsMatrix_->OperatorDomainMap().NumMyElements(),
693 0, SubdomainCrsMatrix_->Comm()));
695 SubdomainMatrixReindexer_.reset(
new EpetraExt::CrsMatrix_Reindex(*tempMap_));
696 DomainVectorReindexer_.reset(
new EpetraExt::MultiVector_Reindex(*tempDomainMap_));
697 RangeVectorReindexer_.reset(
new EpetraExt::MultiVector_Reindex(*tempRangeMap_));
699 ReindexedCrsMatrix_.reset(&((*SubdomainMatrixReindexer_)(*SubdomainCrsMatrix_)),
false);
701 MatrixPtr = &*ReindexedCrsMatrix_;
703 Inverse_ = Teuchos::rcp(
new T(&*ReindexedCrsMatrix_) );
705 Inverse_ = Teuchos::rcp(
new T(&*SubdomainCrsMatrix_) );
708 Inverse_ = Teuchos::rcp(
new T(MatrixPtr) );
711 if (Inverse_ == Teuchos::null)
721 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
722 MpiRank_ = Matrix_->Comm().MyPID();
723 NumMpiProcs_ = Matrix_->Comm().NumProc();
724 NumMpiProcsPerSubdomain_ = List_in.get(
"subdomain: number-of-processors", 1);
725 NumSubdomains_ = NumMpiProcs_ / NumMpiProcsPerSubdomain_;
726 SubdomainId_ = MpiRank_ / NumMpiProcsPerSubdomain_;
728 if (NumSubdomains_ == 1) {
729 IsOverlapping_ =
false;
734 ComputeCondest_ = List_in.get(
"schwarz: compute condest", ComputeCondest_);
736 if( Teuchos::ParameterEntry *combineModeEntry = List_in.getEntryPtr(
"schwarz: combine mode") )
738 if(
typeid(std::string) == combineModeEntry->getAny().type() )
740 std::string mode = List_in.get(
"schwarz: combine mode",
"Add");
743 else if (mode ==
"Zero")
745 else if (mode ==
"Insert")
747 else if (mode ==
"InsertAdd")
749 else if (mode ==
"Average")
751 else if (mode ==
"AbsMax")
755 TEUCHOS_TEST_FOR_EXCEPTION(
756 true,std::logic_error
757 ,
"Error, The (Epetra) combine mode of \""<<mode<<
"\" is not valid! Only the values"
758 " \"Add\", \"Zero\", \"Insert\", \"InsertAdd\", \"Average\", and \"AbsMax\" are accepted!"
764 CombineMode_ = Teuchos::any_cast<Epetra_CombineMode>(combineModeEntry->getAny());
769 Teuchos::getParameter<std::string>(List_in,
"schwarz: combine mode");
775 List_in.get(
"schwarz: combine mode",
"Zero");
778 ReorderingType_ = List_in.get(
"schwarz: reordering type", ReorderingType_);
779 if (ReorderingType_ ==
"none")
780 UseReordering_ =
false;
782 UseReordering_ =
true;
787 FilterSingletons_ = List_in.get(
"schwarz: filter singletons", FilterSingletons_);
799 IsInitialized_ =
false;
803 if (Time_ == Teuchos::null)
806 Time_->ResetStartTime();
809 if (IsOverlapping_) {
810 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
811 if (NumMpiProcsPerSubdomain_ > 1) {
817 # ifdef IFPACK_NODE_AWARE_CODE
819 try{ myNodeID = List_.get(
"ML node id",-1);}
820 catch(...){fprintf(stderr,
"pid %d: no such entry (returned %d)\n",Comm().MyPID(),myNodeID);}
832 if (OverlappingMatrix_ == Teuchos::null) {
837 # ifdef IFPACK_NODE_AWARE_CODE
844 IFPACK_CHK_ERR(Setup());
846 # ifdef IFPACK_NODE_AWARE_CODE
853 if (Inverse_ == Teuchos::null)
856 if (LocalizedMatrix_ == Teuchos::null)
859 IFPACK_CHK_ERR(Inverse_->SetUseTranspose(UseTranspose()));
860 IFPACK_CHK_ERR(Inverse_->SetParameters(List_));
861 IFPACK_CHK_ERR(Inverse_->Initialize());
864 Label_ =
"Ifpack_AdditiveSchwarz, ";
866 Label_ +=
", transp";
868 +
", local solver = \n\t\t***** `" + std::string(Inverse_->Label()) +
"'";
870 IsInitialized_ =
true;
872 InitializeTime_ += Time_->ElapsedTime();
874 #ifdef IFPACK_FLOPCOUNTERS
878 double partial = Inverse_->InitializeFlops();
880 Comm().SumAll(&partial, &total, 1);
881 InitializeFlops_ += total;
891 if (IsInitialized() ==
false)
892 IFPACK_CHK_ERR(Initialize());
894 Time_->ResetStartTime();
898 IFPACK_CHK_ERR(Inverse_->Compute());
902 ComputeTime_ += Time_->ElapsedTime();
904 #ifdef IFPACK_FLOPCOUNTERS
906 double partial = Inverse_->ComputeFlops();
908 Comm().SumAll(&partial, &total, 1);
909 ComputeFlops_ += total;
915 R = ReorderingType_ +
" reord, ";
918 Condest(Ifpack_Cheap);
921 Label_ =
"Ifpack_AdditiveSchwarz, ov = " +
Ifpack_toString(OverlapLevel_)
922 +
", local solver = \n\t\t***** `" + std::string(Inverse_->Label()) +
"'"
923 +
"\n\t\t***** " + R +
"Condition number estimate = "
935 UseTranspose_ = UseTranspose_in;
938 if (Inverse_!=Teuchos::null)
939 IFPACK_CHK_ERR(Inverse_->SetUseTranspose(UseTranspose_in));
948 IFPACK_CHK_ERR(Matrix_->Apply(X,Y));
963 return(Label_.c_str());
970 return(UseTranspose_);
984 return(Matrix_->Comm());
991 return(Matrix_->OperatorDomainMap());
998 return(Matrix_->OperatorRangeMap());
1002 template<
typename T>
1015 Time_->ResetStartTime();
1017 Teuchos::RefCountPtr<Epetra_MultiVector> OverlappingX;
1018 Teuchos::RefCountPtr<Epetra_MultiVector> OverlappingY;
1019 Teuchos::RefCountPtr<Epetra_MultiVector> Xtmp;
1022 #ifdef IFPACK_FLOPCOUNTERS
1023 double pre_partial = Inverse_->ApplyInverseFlops();
1025 Comm().SumAll(&pre_partial, &pre_total, 1);
1029 if (IsOverlapping()) {
1030 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
1031 if (OverlappingX == Teuchos::null) {
1033 if (OverlappingX == Teuchos::null) IFPACK_CHK_ERR(-5);
1034 }
else assert(OverlappingX->NumVectors() == X.
NumVectors());
1035 if (OverlappingY == Teuchos::null) {
1037 if (OverlappingY == Teuchos::null) IFPACK_CHK_ERR(-5);
1038 }
else assert(OverlappingY->NumVectors() == Y.
NumVectors());
1040 # ifdef IFPACK_NODE_AWARE_CODE
1041 if (OverlappingX == Teuchos::null) {
1042 OverlappingX = Teuchos::rcp(
new Epetra_MultiVector(OverlappingMatrix_->RowMatrixRowMap(),
1044 if (OverlappingX == Teuchos::null) IFPACK_CHK_ERR(-5);
1045 }
else assert(OverlappingX->NumVectors() == X.
NumVectors());
1046 if (OverlappingY == Teuchos::null) {
1047 OverlappingY = Teuchos::rcp(
new Epetra_MultiVector(OverlappingMatrix_->RowMatrixRowMap(),
1049 if (OverlappingY == Teuchos::null) IFPACK_CHK_ERR(-5);
1050 }
else assert(OverlappingY->NumVectors() == Y.
NumVectors());
1052 OverlappingX = Teuchos::rcp(
new Epetra_MultiVector(OverlappingMatrix_->RowMatrixRowMap(),
1054 OverlappingY = Teuchos::rcp(
new Epetra_MultiVector(OverlappingMatrix_->RowMatrixRowMap(),
1056 if (OverlappingY == Teuchos::null) IFPACK_CHK_ERR(-5);
1059 OverlappingY->PutScalar(0.0);
1060 OverlappingX->PutScalar(0.0);
1061 IFPACK_CHK_ERR(OverlappingMatrix_->ImportMultiVector(X,*OverlappingX,
Insert));
1068 OverlappingX = Xtmp;
1069 OverlappingY = Teuchos::rcp( &Y,
false );
1072 if (FilterSingletons_) {
1076 IFPACK_CHK_ERR(SingletonFilter_->SolveSingletons(*OverlappingX,*OverlappingY));
1077 IFPACK_CHK_ERR(SingletonFilter_->CreateReducedRHS(*OverlappingY,*OverlappingX,ReducedX));
1080 if (!UseReordering_) {
1081 IFPACK_CHK_ERR(Inverse_->ApplyInverse(ReducedX,ReducedY));
1086 IFPACK_CHK_ERR(Reordering_->P(ReducedX,ReorderedX));
1087 IFPACK_CHK_ERR(Inverse_->ApplyInverse(ReorderedX,ReorderedY));
1088 IFPACK_CHK_ERR(Reordering_->Pinv(ReorderedY,ReducedY));
1092 IFPACK_CHK_ERR(SingletonFilter_->UpdateLHS(ReducedY,*OverlappingY));
1096 if (!UseReordering_) {
1097 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
1098 if (NumMpiProcsPerSubdomain_ > 1) {
1099 tempX_.reset(&((*RangeVectorReindexer_)(*OverlappingX)),
false);
1100 tempY_.reset(&((*DomainVectorReindexer_)(*OverlappingY)),
false);
1101 IFPACK_CHK_ERR(Inverse_->ApplyInverse(*tempX_,*tempY_));
1103 IFPACK_CHK_ERR(Inverse_->ApplyInverse(*OverlappingX, *OverlappingY));
1106 IFPACK_CHK_ERR(Inverse_->ApplyInverse(*OverlappingX,*OverlappingY));
1112 IFPACK_CHK_ERR(Reordering_->P(*OverlappingX,ReorderedX));
1113 IFPACK_CHK_ERR(Inverse_->ApplyInverse(ReorderedX,ReorderedY));
1114 IFPACK_CHK_ERR(Reordering_->Pinv(ReorderedY,*OverlappingY));
1118 if (IsOverlapping()) {
1119 IFPACK_CHK_ERR(OverlappingMatrix_->ExportMultiVector(*OverlappingY,Y,
1123 #ifdef IFPACK_FLOPCOUNTERS
1126 double partial = Inverse_->ApplyInverseFlops();
1128 Comm().SumAll(&partial, &total, 1);
1129 ApplyInverseFlops_ += total - pre_total;
1134 ApplyInverseTime_ += Time_->ElapsedTime();
1141 template<
typename T>
1147 #ifdef IFPACK_FLOPCOUNTERS
1148 double IF = InitializeFlops();
1149 double CF = ComputeFlops();
1150 double AF = ApplyInverseFlops();
1152 double IFT = 0.0, CFT = 0.0, AFT = 0.0;
1153 if (InitializeTime() != 0.0)
1154 IFT = IF / InitializeTime();
1155 if (ComputeTime() != 0.0)
1156 CFT = CF / ComputeTime();
1157 if (ApplyInverseTime() != 0.0)
1158 AFT = AF / ApplyInverseTime();
1161 if (Matrix().Comm().MyPID())
1165 os <<
"================================================================================" << endl;
1166 os <<
"Ifpack_AdditiveSchwarz, overlap level = " << OverlapLevel_ << endl;
1167 if (CombineMode_ ==
Insert)
1168 os <<
"Combine mode = Insert" << endl;
1169 else if (CombineMode_ ==
Add)
1170 os <<
"Combine mode = Add" << endl;
1171 else if (CombineMode_ ==
Zero)
1172 os <<
"Combine mode = Zero" << endl;
1173 else if (CombineMode_ ==
Average)
1174 os <<
"Combine mode = Average" << endl;
1175 else if (CombineMode_ ==
AbsMax)
1176 os <<
"Combine mode = AbsMax" << endl;
1178 os <<
"Condition number estimate = " << Condest_ << endl;
1179 os <<
"Global number of rows = " << Matrix_->NumGlobalRows64() << endl;
1181 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
1183 os <<
"================================================================================" << endl;
1184 os <<
"Subcommunicator stats" << endl;
1185 os <<
"Number of MPI processes in simulation: " << NumMpiProcs_ << endl;
1186 os <<
"Number of subdomains: " << NumSubdomains_ << endl;
1187 os <<
"Number of MPI processes per subdomain: " << NumMpiProcsPerSubdomain_ << endl;
1191 os <<
"Phase # calls Total Time (s) Total MFlops MFlops/s" << endl;
1192 os <<
"----- ------- -------------- ------------ --------" << endl;
1193 os <<
"Initialize() " << std::setw(5) << NumInitialize()
1194 <<
" " << std::setw(15) << InitializeTime()
1195 #ifdef IFPACK_FLOPCOUNTERS
1196 <<
" " << std::setw(15) << 1.0e-6 * IF
1197 <<
" " << std::setw(15) << 1.0e-6 * IFT
1200 os <<
"Compute() " << std::setw(5) << NumCompute()
1201 <<
" " << std::setw(15) << ComputeTime()
1202 #ifdef IFPACK_FLOPCOUNTERS
1203 <<
" " << std::setw(15) << 1.0e-6 * CF
1204 <<
" " << std::setw(15) << 1.0e-6 * CFT
1207 os <<
"ApplyInverse() " << std::setw(5) << NumApplyInverse()
1208 <<
" " << std::setw(15) << ApplyInverseTime()
1209 #ifdef IFPACK_FLOPCOUNTERS
1210 <<
" " << std::setw(15) << 1.0e-6 * AF
1211 <<
" " << std::setw(15) << 1.0e-6 * AFT
1214 os <<
"================================================================================" << endl;
1220 #include "Ifpack_Condest.h"
1222 template<
typename T>
1224 Condest(
const Ifpack_CondestType CT,
const int MaxIters,
1230 Condest_ = Ifpack_Condest(*
this, CT, MaxIters, Tol, Matrix_in);
1235 #endif // IFPACK_ADDITIVESCHWARZ_H