46 #ifndef MUELU_MLPARAMETERLISTINTERPRETER_DEF_HPP
47 #define MUELU_MLPARAMETERLISTINTERPRETER_DEF_HPP
52 #if defined(HAVE_MUELU_ML) && defined(HAVE_MUELU_EPETRA)
53 #include <ml_ValidateParameters.h>
64 #include "MueLu_Hierarchy.hpp"
65 #include "MueLu_FactoryManager.hpp"
67 #include "MueLu_TentativePFactory.hpp"
68 #include "MueLu_SaPFactory.hpp"
69 #include "MueLu_PgPFactory.hpp"
70 #include "MueLu_AmalgamationFactory.hpp"
71 #include "MueLu_TransPFactory.hpp"
72 #include "MueLu_GenericRFactory.hpp"
73 #include "MueLu_SmootherPrototype.hpp"
74 #include "MueLu_SmootherFactory.hpp"
75 #include "MueLu_TrilinosSmoother.hpp"
77 #include "MueLu_DirectSolver.hpp"
78 #include "MueLu_HierarchyUtils.hpp"
79 #include "MueLu_RAPFactory.hpp"
80 #include "MueLu_CoalesceDropFactory.hpp"
81 #include "MueLu_CoupledAggregationFactory.hpp"
82 #include "MueLu_UncoupledAggregationFactory.hpp"
83 #include "MueLu_HybridAggregationFactory.hpp"
84 #include "MueLu_NullspaceFactory.hpp"
87 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
88 #include "MueLu_CoalesceDropFactory_kokkos.hpp"
92 #include "MueLu_SaPFactory_kokkos.hpp"
93 #include "MueLu_TentativePFactory_kokkos.hpp"
94 #include "MueLu_UncoupledAggregationFactory_kokkos.hpp"
97 #if defined(HAVE_MUELU_ISORROPIA) && defined(HAVE_MPI)
98 #include "MueLu_IsorropiaInterface.hpp"
99 #include "MueLu_RepartitionHeuristicFactory.hpp"
100 #include "MueLu_RepartitionFactory.hpp"
101 #include "MueLu_RebalanceTransferFactory.hpp"
102 #include "MueLu_RepartitionInterface.hpp"
103 #include "MueLu_RebalanceAcFactory.hpp"
113 #define MUELU_READ_PARAM(paramList, paramStr, varType, defaultValue, varName) \
114 varType varName = defaultValue; if (paramList.isParameter(paramStr)) varName = paramList.get<varType>(paramStr);
117 #define MUELU_COPY_PARAM(paramList, paramStr, varType, defaultValue, outParamList, outParamStr) \
118 if (paramList.isParameter(paramStr)) \
119 outParamList.set(outParamStr, paramList.get<varType>(paramStr)); \
120 else outParamList.set(outParamStr, static_cast<varType>(defaultValue)); \
124 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
125 MLParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MLParameterListInterpreter(
Teuchos::ParameterList & paramList,
Teuchos::RCP<
const Teuchos::Comm<int> > comm, std::vector<
RCP<FactoryBase> > factoryList) : nullspace_(NULL), xcoord_(NULL), ycoord_(NULL), zcoord_(NULL),TransferFacts_(factoryList), blksize_(1) {
128 std::string filename = paramList.
get(
"xml parameter file",
"");
129 if (filename.length() != 0) {
133 paramList2.remove(
"xml parameter file");
143 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
149 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
164 MUELU_READ_PARAM(paramList,
"aggregation: type", std::string,
"Uncoupled", agg_type);
166 MUELU_READ_PARAM(paramList,
"aggregation: damping factor",
double, (
double)4/(
double)3, agg_damping);
168 MUELU_READ_PARAM(paramList,
"aggregation: nodes per aggregate",
int, 1, minPerAgg);
169 MUELU_READ_PARAM(paramList,
"aggregation: keep Dirichlet bcs",
bool,
false, bKeepDirichletBcs);
170 MUELU_READ_PARAM(paramList,
"aggregation: max neighbours already aggregated",
int, 0, maxNbrAlreadySelected);
171 MUELU_READ_PARAM(paramList,
"aggregation: aux: enable",
bool,
false, agg_use_aux);
172 MUELU_READ_PARAM(paramList,
"aggregation: aux: threshold",
double,
false, agg_aux_thresh);
174 MUELU_READ_PARAM(paramList,
"null space: type", std::string,
"default vectors", nullspaceType);
175 MUELU_READ_PARAM(paramList,
"null space: dimension",
int, -1, nullspaceDim);
176 MUELU_READ_PARAM(paramList,
"null space: vectors",
double*, NULL, nullspaceVec);
178 MUELU_READ_PARAM(paramList,
"energy minimization: enable",
bool,
false, bEnergyMinimization);
180 MUELU_READ_PARAM(paramList,
"RAP: fix diagonal",
bool,
false, bFixDiagonal);
195 paramList = paramListWithSubList;
198 bool setKokkosRefactor =
false;
199 #if ( defined(HAVE_MUELU_KOKKOS_REFACTOR) && defined(HAVE_MUELU_KOKKOS_REFACTOR_USE_BY_DEFAULT) )
200 bool useKokkosRefactor =
true;
202 bool useKokkosRefactor =
false;
204 if (paramList.
isType<
bool>(
"use kokkos refactor")) {
205 useKokkosRefactor = paramList.
get<
bool>(
"use kokkos refactor");
206 setKokkosRefactor =
true;
207 paramList.
remove(
"use kokkos refactor");
215 bool validate = paramList.
get(
"ML validate parameter list",
true);
218 #if defined(HAVE_MUELU_ML) && defined(HAVE_MUELU_EPETRA)
220 int depth = paramList.
get(
"ML validate depth", 5);
222 "ERROR: ML's Teuchos::ParameterList contains incorrect parameter!");
225 this->GetOStream(
Warnings0) <<
"Warning: MueLu_ENABLE_ML=OFF. The parameter list cannot be validated." << std::endl;
226 paramList.
set(
"ML validate parameter list",
false);
228 #endif // HAVE_MUELU_ML
234 blksize_ = nDofsPerNode;
240 if (verbosityLevel == 0) eVerbLevel =
None;
241 if (verbosityLevel >= 1) eVerbLevel =
Low;
242 if (verbosityLevel >= 5) eVerbLevel =
Medium;
243 if (verbosityLevel >= 10) eVerbLevel =
High;
244 if (verbosityLevel >= 11) eVerbLevel =
Extreme;
245 if (verbosityLevel >= 42) eVerbLevel =
Test;
246 this->verbosity_ = eVerbLevel;
250 "MueLu::MLParameterListInterpreter::SetParameterList(): parameter \"aggregation: type\": only 'Uncoupled' or 'Coupled' aggregation is supported.");
254 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
255 if(useKokkosRefactor)
256 dropFact =
rcp(
new CoalesceDropFactory_kokkos() );
267 if (agg_type ==
"Uncoupled") {
270 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
271 if(useKokkosRefactor) {
272 MyUncoupledAggFact =
rcp(
new UncoupledAggregationFactory_kokkos() );
278 MyUncoupledAggFact->
SetFactory(
"Graph", dropFact);
279 MyUncoupledAggFact->
SetFactory(
"DofsPerNode", dropFact);
285 AggFact = MyUncoupledAggFact;
288 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
289 if(useKokkosRefactor) {
290 AggFact =
rcp(
new UncoupledAggregationFactory_kokkos() );
297 CoupledAggFact2->
SetFactory(
"Graph", dropFact);
298 CoupledAggFact2->
SetFactory(
"DofsPerNode", dropFact);
299 AggFact = CoupledAggFact2;
308 CoupledAggFact2->
SetFactory(
"Graph", dropFact);
309 CoupledAggFact2->
SetFactory(
"DofsPerNode", dropFact);
310 AggFact = CoupledAggFact2;
313 if (verbosityLevel > 3) {
314 std::ostringstream oss;
315 oss <<
"========================= Aggregate option summary  =========================" << std::endl;
316 oss <<
"min Nodes per aggregate : Â Â Â Â Â Â Â " << minPerAgg << std::endl;
317 oss <<
"min # of root nbrs already aggregated : " << maxNbrAlreadySelected << std::endl;
318 oss <<
"aggregate ordering : Â Â Â Â Â Â Â Â Â Â natural" << std::endl;
319 oss <<
"=============================================================================" << std::endl;
320 this->GetOStream(
Runtime1) << oss.str();
326 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
327 if(useKokkosRefactor)
328 PtentFact =
rcp(
new TentativePFactory_kokkos() );
332 if (agg_damping == 0.0 && bEnergyMinimization ==
false) {
336 }
else if (agg_damping != 0.0 && bEnergyMinimization ==
false) {
339 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
340 if(useKokkosRefactor)
341 SaPFact =
rcp(
new SaPFactory_kokkos() );
348 }
else if (bEnergyMinimization ==
true) {
356 for (
size_t i = 0; i<TransferFacts_.size(); i++) {
363 #if defined(HAVE_MUELU_ISORROPIA) && defined(HAVE_MPI)
370 if (bDoRepartition == 1) {
380 rebAmalgFact->SetFactory(
"A", AcFact);
382 MUELU_READ_PARAM(paramList,
"repartition: max min ratio",
double, 1.3, maxminratio);
383 MUELU_READ_PARAM(paramList,
"repartition: min per proc",
int, 512, minperproc);
389 paramListRepFact.
set(
"repartition: min rows per proc", minperproc);
390 paramListRepFact.
set(
"repartition: max imbalance", maxminratio);
393 RepartitionHeuristicFact->
SetFactory(
"A", AcFact);
397 isoInterface->SetFactory(
"A", AcFact);
398 isoInterface->SetFactory(
"number of partitions", RepartitionHeuristicFact);
399 isoInterface->SetFactory(
"UnAmalgamationInfo", rebAmalgFact);
403 repInterface->SetFactory(
"A", AcFact);
404 repInterface->SetFactory(
"number of partitions", RepartitionHeuristicFact);
405 repInterface->SetFactory(
"AmalgamatedPartition", isoInterface);
411 RepartitionFact->
SetFactory(
"number of partitions", RepartitionHeuristicFact);
412 RepartitionFact->
SetFactory(
"Partition", repInterface);
418 RebalancedPFact->
SetFactory(
"Nullspace", PtentFact);
419 RebalancedPFact->
SetFactory(
"Importer", RepartitionFact);
424 RebalancedRFact->
SetFactory(
"Importer", RepartitionFact);
430 #else // #ifdef HAVE_MUELU_ISORROPIA
444 if (nullspaceType !=
"default vectors") {
449 nullspaceDim_ = nullspaceDim;
450 nullspace_ = nullspaceVec;
469 this->numDesiredLevel_ = maxLevels;
470 this->maxCoarseSize_ = maxCoarseSize;
478 coarseList.
set(
"smoother: type",
"Amesos-KLU");
490 for (
int levelID=0; levelID < maxLevels; levelID++) {
497 if (setKokkosRefactor)
515 manager->
SetFactory(
"Smoother", smootherFact);
522 manager->
SetFactory(
"CoarseSolver", coarseFact);
528 #if defined(HAVE_MUELU_ISORROPIA) && defined(HAVE_MPI)
529 if (bDoRepartition == 1) {
533 manager->
SetFactory(
"Nullspace", RebalancedPFact);
534 manager->
SetFactory(
"Importer", RepartitionFact);
536 #endif // #ifdef HAVE_MUELU_ISORROPIA
541 #if defined(HAVE_MUELU_ISORROPIA) && defined(HAVE_MPI)
545 this->AddFactoryManager(levelID, 1, manager);
550 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
554 if (nullspace_ != NULL) {
560 RCP<MultiVector> nullspace = MultiVectorFactory::Build(rowMap, nullspaceDim_,
true);
562 for (
size_t i=0; i < Teuchos::as<size_t>(nullspaceDim_); i++) {
564 const size_t myLength = nullspace->getLocalLength();
566 for (
size_t j = 0; j < myLength; j++) {
567 nullspacei[j] = nullspace_[i*myLength + j];
571 fineLevel->
Set(
"Nullspace", nullspace);
576 size_t num_coords = 0;
577 double * coordPTR[3];
579 coordPTR[0] = xcoord_;
582 coordPTR[1] = ycoord_;
585 coordPTR[2] = zcoord_;
598 for (
size_t i=0; i < num_coords; i++) {
600 const size_t myLength = coordinates->getLocalLength();
601 for (
size_t j = 0; j < myLength; j++) {
602 coordsi[j] = coordPTR[0][j];
605 fineLevel->
Set(
"Coordinates",coordinates);
613 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
622 std::string type =
"symmetric Gauss-Seidel";
644 if (paramList.
isParameter(
"smoother: type")) type = paramList.
get<std::string>(
"smoother: type");
652 std::string ifpackType;
655 if (type ==
"Jacobi" || type ==
"Gauss-Seidel" || type ==
"symmetric Gauss-Seidel") {
656 if (type ==
"symmetric Gauss-Seidel") type =
"Symmetric Gauss-Seidel";
658 ifpackType =
"RELAXATION";
659 smootherParamList.
set(
"relaxation: type", type);
661 MUELU_COPY_PARAM(paramList,
"smoother: sweeps",
int, 2, smootherParamList,
"relaxation: sweeps");
662 MUELU_COPY_PARAM(paramList,
"smoother: damping factor", Scalar, one, smootherParamList,
"relaxation: damping factor");
665 smooProto->SetFactory(
"A", AFact);
667 }
else if (type ==
"Chebyshev" || type ==
"MLS") {
669 ifpackType =
"CHEBYSHEV";
671 MUELU_COPY_PARAM(paramList,
"smoother: sweeps",
int, 2, smootherParamList,
"chebyshev: degree");
672 if (paramList.
isParameter(
"smoother: MLS alpha")) {
673 MUELU_COPY_PARAM(paramList,
"smoother: MLS alpha",
double, 20, smootherParamList,
"chebyshev: ratio eigenvalue");
675 MUELU_COPY_PARAM(paramList,
"smoother: Chebyshev alpha",
double, 20, smootherParamList,
"chebyshev: ratio eigenvalue");
680 smooProto->SetFactory(
"A", AFact);
682 }
else if (type ==
"IFPACK") {
684 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_IFPACK)
685 ifpackType = paramList.
get<std::string>(
"smoother: ifpack type");
687 if (ifpackType ==
"ILU") {
690 if (paramList.
isParameter(
"smoother: ifpack level-of-fill"))
691 smootherParamList.
set(
"fact: level-of-fill", Teuchos::as<int>(paramList.
get<
double>(
"smoother: ifpack level-of-fill")));
692 else smootherParamList.
set(
"fact: level-of-fill", as<int>(0));
694 MUELU_COPY_PARAM(paramList,
"smoother: ifpack overlap",
int, 2, smootherParamList,
"partitioner: overlap");
698 MueLu::GetIfpackSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node> (ifpackType,
700 paramList.
get<
int> (
"smoother: ifpack overlap"));
701 smooProto->SetFactory(
"A", AFact);
709 }
else if (type.length() > strlen(
"Amesos") && type.substr(0, strlen(
"Amesos")) ==
"Amesos") {
710 std::string solverType = type.substr(strlen(
"Amesos")+1);
714 const int validatorSize = 5;
715 std::string validator[validatorSize] = {
"Superlu",
"Superludist",
"KLU",
"UMFPACK"};
716 for (
int i=0; i < validatorSize; i++) {
if (validator[i] == solverType) valid =
true; }
720 std::transform(solverType.begin()+1, solverType.end(), solverType.begin()+1,
::tolower);
723 smooProto->SetFactory(
"A", AFact);
739 MUELU_READ_PARAM(paramList,
"smoother: pre or post", std::string,
"both", preOrPost);
740 if (preOrPost ==
"both") {
742 }
else if (preOrPost ==
"pre") {
744 }
else if (preOrPost ==
"post") {
751 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
754 TEUCHOS_TEST_FOR_EXCEPTION(Teuchos::rcp_dynamic_cast<TwoLevelFactoryBase>(factory) == Teuchos::null,
Exceptions::BadCast,
"Transfer factory is not derived from TwoLevelFactoryBase. Since transfer factories will be handled by the RAPFactory they have to be derived from TwoLevelFactoryBase!");
755 TransferFacts_.push_back(factory);
758 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
760 return TransferFacts_.size();
763 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
766 Matrix& A = dynamic_cast<Matrix&>(Op);
767 if (A.GetFixedBlockSize() != blksize_)
768 this->GetOStream(
Warnings0) <<
"Setting matrix block size to " << blksize_ <<
" (value of the parameter in the list) "
769 <<
"instead of " << A.GetFixedBlockSize() <<
" (provided matrix)." << std::endl;
771 A.SetFixedBlockSize(blksize_);
773 }
catch (std::bad_cast& e) {
774 this->GetOStream(
Warnings0) <<
"Skipping setting block size as the operator is not a matrix" << std::endl;
780 #define MUELU_MLPARAMETERLISTINTERPRETER_SHORT