MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_MLParameterListInterpreter_def.hpp
Go to the documentation of this file.
1// @HEADER
2//
3// ***********************************************************************
4//
5// MueLu: A package for multigrid based preconditioning
6// Copyright 2012 Sandia Corporation
7//
8// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9// the U.S. Government retains certain rights in this software.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact
39// Jonathan Hu (jhu@sandia.gov)
40// Andrey Prokopenko (aprokop@sandia.gov)
41// Ray Tuminaro (rstumin@sandia.gov)
42//
43// ***********************************************************************
44//
45// @HEADER
46#ifndef MUELU_MLPARAMETERLISTINTERPRETER_DEF_HPP
47#define MUELU_MLPARAMETERLISTINTERPRETER_DEF_HPP
48
49#include <Teuchos_XMLParameterListHelpers.hpp>
50
51#include "MueLu_ConfigDefs.hpp"
52#if defined(HAVE_MUELU_ML)
53#include <ml_ValidateParameters.h>
54#endif
55
56#include <Xpetra_Matrix.hpp>
57#include <Xpetra_MatrixUtils.hpp>
58#include <Xpetra_MultiVector.hpp>
59#include <Xpetra_MultiVectorFactory.hpp>
60#include <Xpetra_Operator.hpp>
61
63
64#include "MueLu_Level.hpp"
65#include "MueLu_Hierarchy.hpp"
66#include "MueLu_FactoryManager.hpp"
67
68#include "MueLu_TentativePFactory.hpp"
69#include "MueLu_SaPFactory.hpp"
70#include "MueLu_PgPFactory.hpp"
71#include "MueLu_AmalgamationFactory.hpp"
72#include "MueLu_TransPFactory.hpp"
73#include "MueLu_GenericRFactory.hpp"
74#include "MueLu_SmootherPrototype.hpp"
75#include "MueLu_SmootherFactory.hpp"
76#include "MueLu_TrilinosSmoother.hpp"
78#include "MueLu_DirectSolver.hpp"
79#include "MueLu_HierarchyUtils.hpp"
80#include "MueLu_RAPFactory.hpp"
81#include "MueLu_CoalesceDropFactory.hpp"
82#include "MueLu_UncoupledAggregationFactory.hpp"
83#include "MueLu_NullspaceFactory.hpp"
85
86#include "MueLu_CoalesceDropFactory_kokkos.hpp"
87// #include "MueLu_CoordinatesTransferFactory_kokkos.hpp"
88// #include "MueLu_NullspaceFactory_kokkos.hpp"
89#include "MueLu_SaPFactory_kokkos.hpp"
90#include "MueLu_TentativePFactory_kokkos.hpp"
91#include "MueLu_UncoupledAggregationFactory_kokkos.hpp"
92
93#if defined(HAVE_MUELU_ISORROPIA) && defined(HAVE_MPI)
94#include "MueLu_IsorropiaInterface.hpp"
95#include "MueLu_RepartitionHeuristicFactory.hpp"
96#include "MueLu_RepartitionFactory.hpp"
97#include "MueLu_RebalanceTransferFactory.hpp"
98#include "MueLu_RepartitionInterface.hpp"
99#include "MueLu_RebalanceAcFactory.hpp"
100//#include "MueLu_RebalanceMapFactory.hpp"
101#endif
102
103// Note: do not add options that are only recognized by MueLu.
104
105// TODO: this parameter list interpreter should force MueLu to use default ML parameters
106// - Ex: smoother sweep=2 by default for ML
107
108// Read a parameter value from a parameter list and store it into a variable named 'varName'
109#define MUELU_READ_PARAM(paramList, paramStr, varType, defaultValue, varName) \
110 varType varName = defaultValue; if (paramList.isParameter(paramStr)) varName = paramList.get<varType>(paramStr);
111
112// Read a parameter value from a paraeter list and copy it into a new parameter list (with another parameter name)
113#define MUELU_COPY_PARAM(paramList, paramStr, varType, defaultValue, outParamList, outParamStr) \
114 if (paramList.isParameter(paramStr)) \
115 outParamList.set(outParamStr, paramList.get<varType>(paramStr)); \
116 else outParamList.set(outParamStr, static_cast<varType>(defaultValue)); \
117
118namespace MueLu {
119
120 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
121 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) {
122
123 if (paramList.isParameter("xml parameter file")){
124 std::string filename = paramList.get("xml parameter file","");
125 if (filename.length() != 0) {
126 TEUCHOS_TEST_FOR_EXCEPTION(comm.is_null(), Exceptions::RuntimeError, "xml parameter file requires a valid comm");
127 Teuchos::ParameterList paramList2 = paramList;
128 Teuchos::updateParametersFromXmlFileAndBroadcast(filename, Teuchos::Ptr<Teuchos::ParameterList>(&paramList2),*comm);
129 paramList2.remove("xml parameter file");
130 SetParameterList(paramList2);
131 }
132 else
133 SetParameterList(paramList);
134 }
135 else
136 SetParameterList(paramList);
137 }
138
139 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
140 MLParameterListInterpreter<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MLParameterListInterpreter(const std::string & xmlFileName, std::vector<RCP<FactoryBase> > factoryList) : nullspace_(NULL), TransferFacts_(factoryList), blksize_(1) {
141 Teuchos::RCP<Teuchos::ParameterList> paramList = Teuchos::getParametersFromXmlFile(xmlFileName);
142 SetParameterList(*paramList);
143 }
144
145 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
147 Teuchos::ParameterList paramList = paramList_in;
148
149 //
150 // Read top-level of the parameter list
151 //
152
153 // hard-coded default values == ML defaults according to the manual
154 MUELU_READ_PARAM(paramList, "ML output", int, 0, verbosityLevel);
155 MUELU_READ_PARAM(paramList, "max levels", int, 10, maxLevels);
156 MUELU_READ_PARAM(paramList, "PDE equations", int, 1, nDofsPerNode);
157
158 MUELU_READ_PARAM(paramList, "coarse: max size", int, 128, maxCoarseSize);
159
160 MUELU_READ_PARAM(paramList, "aggregation: type", std::string, "Uncoupled", agg_type);
161 //MUELU_READ_PARAM(paramList, "aggregation: threshold", double, 0.0, agg_threshold);
162 MUELU_READ_PARAM(paramList, "aggregation: damping factor", double, (double)4/(double)3, agg_damping);
163 //MUELU_READ_PARAM(paramList, "aggregation: smoothing sweeps", int, 1, agg_smoothingsweeps);
164 MUELU_READ_PARAM(paramList, "aggregation: nodes per aggregate", int, 1, minPerAgg);
165 MUELU_READ_PARAM(paramList, "aggregation: keep Dirichlet bcs", bool, false, bKeepDirichletBcs); // This is a MueLu specific extension that does not exist in ML
166 MUELU_READ_PARAM(paramList, "aggregation: max neighbours already aggregated", int, 0, maxNbrAlreadySelected); // This is a MueLu specific extension that does not exist in M
167 MUELU_READ_PARAM(paramList, "aggregation: aux: enable", bool, false, agg_use_aux);
168 MUELU_READ_PARAM(paramList, "aggregation: aux: threshold", double, false, agg_aux_thresh);
169
170 MUELU_READ_PARAM(paramList, "null space: type", std::string, "default vectors", nullspaceType);
171 MUELU_READ_PARAM(paramList, "null space: dimension", int, -1, nullspaceDim); // TODO: ML default not in documentation
172 MUELU_READ_PARAM(paramList, "null space: vectors", double*, NULL, nullspaceVec); // TODO: ML default not in documentation
173
174 MUELU_READ_PARAM(paramList, "energy minimization: enable", bool, false, bEnergyMinimization);
175
176 MUELU_READ_PARAM(paramList, "RAP: fix diagonal", bool, false, bFixDiagonal); // This is a MueLu specific extension that does not exist in ML
177
178 MUELU_READ_PARAM(paramList, "x-coordinates", double*, NULL, xcoord);
179 MUELU_READ_PARAM(paramList, "y-coordinates", double*, NULL, ycoord);
180 MUELU_READ_PARAM(paramList, "z-coordinates", double*, NULL, zcoord);
181
182
183 //
184 // Move smoothers/aggregation/coarse parameters to sublists
185 //
186
187 // ML allows to have level-specific smoothers/aggregation/coarse parameters at the top level of the list or/and defined in sublists:
188 // See also: ML Guide section 6.4.1, MueLu::CreateSublists, ML_CreateSublists
189 ParameterList paramListWithSubList;
190 MueLu::CreateSublists(paramList, paramListWithSubList);
191 paramList = paramListWithSubList; // swap
192
193 // pull out "use kokkos refactor"
194 bool setKokkosRefactor = false;
195 bool useKokkosRefactor;
196# ifdef HAVE_MUELU_SERIAL
197 if (typeid(Node).name() == typeid(Tpetra::KokkosCompat::KokkosSerialWrapperNode).name())
198 useKokkosRefactor = false;
199# endif
200# ifdef HAVE_MUELU_OPENMP
201 if (typeid(Node).name() == typeid(Tpetra::KokkosCompat::KokkosOpenMPWrapperNode).name())
202 useKokkosRefactor = true;
203# endif
204# ifdef HAVE_MUELU_CUDA
205 if (typeid(Node).name() == typeid(Tpetra::KokkosCompat::KokkosCudaWrapperNode).name())
206 useKokkosRefactor = true;
207# endif
208# ifdef HAVE_MUELU_HIP
209 if (typeid(Node).name() == typeid(Tpetra::KokkosCompat::KokkosHIPWrapperNode).name())
210 useKokkosRefactor = true;
211# endif
212# ifdef HAVE_MUELU_SYCL
213 if (typeid(Node).name() == typeid(Tpetra::KokkosCompat::KokkosSYCLWrapperNode).name())
214 useKokkosRefactor = true;
215#endif
216 if (paramList.isType<bool>("use kokkos refactor")) {
217 useKokkosRefactor = paramList.get<bool>("use kokkos refactor");
218 setKokkosRefactor = true;
219 paramList.remove("use kokkos refactor");
220 }
221
222 //
223 // Validate parameter list
224 //
225
226 {
227 bool validate = paramList.get("ML validate parameter list", true); /* true = default in ML */
228 if (validate) {
229
230#if defined(HAVE_MUELU_ML) && defined(HAVE_MUELU_EPETRA)
231 // Validate parameter list using ML validator
232 int depth = paramList.get("ML validate depth", 5); /* 5 = default in ML */
233 TEUCHOS_TEST_FOR_EXCEPTION(! ML_Epetra::ValidateMLPParameters(paramList, depth), Exceptions::RuntimeError,
234 "ERROR: ML's Teuchos::ParameterList contains incorrect parameter!");
235#else
236 // If no validator available: issue a warning and set parameter value to false in the output list
237 this->GetOStream(Warnings0) << "Warning: MueLu_ENABLE_ML=OFF. The parameter list cannot be validated." << std::endl;
238 paramList.set("ML validate parameter list", false);
239
240#endif // HAVE_MUELU_ML
241 } // if(validate)
242 } // scope
243
244
245 // Matrix option
246 blksize_ = nDofsPerNode;
247
248 // Translate verbosity parameter
249
250 // Translate verbosity parameter
251 MsgType eVerbLevel = None;
252 if (verbosityLevel == 0) eVerbLevel = None;
253 if (verbosityLevel >= 1) eVerbLevel = Low;
254 if (verbosityLevel >= 5) eVerbLevel = Medium;
255 if (verbosityLevel >= 10) eVerbLevel = High;
256 if (verbosityLevel >= 11) eVerbLevel = Extreme;
257 if (verbosityLevel >= 42) eVerbLevel = Test;
258 if (verbosityLevel >= 43) eVerbLevel = InterfaceTest;
259 this->verbosity_ = eVerbLevel;
260
261
262 TEUCHOS_TEST_FOR_EXCEPTION(agg_type != "Uncoupled", Exceptions::RuntimeError,
263 "MueLu::MLParameterListInterpreter::SetParameterList(): parameter \"aggregation: type\": only 'Uncoupled' aggregation is supported.");
264
265 // Create MueLu factories
266 RCP<Factory> dropFact;
267 if(useKokkosRefactor)
268 dropFact = rcp( new CoalesceDropFactory_kokkos() );
269 else
270 dropFact = rcp( new CoalesceDropFactory() );
271
272 if (agg_use_aux) {
273 dropFact->SetParameter("aggregation: drop scheme",Teuchos::ParameterEntry(std::string("distance laplacian")));
274 dropFact->SetParameter("aggregation: drop tol",Teuchos::ParameterEntry(agg_aux_thresh));
275 }
276
277 // Uncoupled aggregation
278 RCP<Factory> AggFact = Teuchos::null;
279 if(useKokkosRefactor) {
280 AggFact = rcp( new UncoupledAggregationFactory_kokkos() );
281 }
282 else
283 AggFact = rcp( new UncoupledAggregationFactory() );
284
285 AggFact->SetFactory("Graph", dropFact);
286 AggFact->SetFactory("DofsPerNode", dropFact);
287 AggFact->SetParameter("aggregation: preserve Dirichlet points", Teuchos::ParameterEntry(bKeepDirichletBcs));
288 AggFact->SetParameter("aggregation: ordering", Teuchos::ParameterEntry(std::string("natural")));
289 AggFact->SetParameter("aggregation: max selected neighbors", Teuchos::ParameterEntry(maxNbrAlreadySelected));
290 AggFact->SetParameter("aggregation: min agg size", Teuchos::ParameterEntry(minPerAgg));
291
292
293 if (verbosityLevel > 3) {
294 std::ostringstream oss;
295 oss << "========================= Aggregate option summary  =========================" << std::endl;
296 oss << "min Nodes per aggregate :              " << minPerAgg << std::endl;
297 oss << "min # of root nbrs already aggregated : " << maxNbrAlreadySelected << std::endl;
298 oss << "aggregate ordering :                    natural" << std::endl;
299 oss << "=============================================================================" << std::endl;
300 this->GetOStream(Runtime1) << oss.str();
301 }
302
303 RCP<Factory> PFact;
304 RCP<Factory> RFact;
305 RCP<Factory> PtentFact;
306 if(useKokkosRefactor)
307 PtentFact = rcp( new TentativePFactory_kokkos() );
308 else
309 PtentFact = rcp( new TentativePFactory() );
310 if (agg_damping == 0.0 && bEnergyMinimization == false) {
311 // tentative prolongation operator (PA-AMG)
312 PFact = PtentFact;
313 RFact = rcp( new TransPFactory() );
314 } else if (agg_damping != 0.0 && bEnergyMinimization == false) {
315 // smoothed aggregation (SA-AMG)
316 RCP<Factory> SaPFact;
317 if(useKokkosRefactor)
318 SaPFact = rcp( new SaPFactory_kokkos() );
319 else
320 SaPFact = rcp( new SaPFactory() );
321 SaPFact->SetParameter("sa: damping factor", ParameterEntry(agg_damping));
322 PFact = SaPFact;
323 RFact = rcp( new TransPFactory() );
324 } else if (bEnergyMinimization == true) {
325 // Petrov Galerkin PG-AMG smoothed aggregation (energy minimization in ML)
326 PFact = rcp( new PgPFactory() );
327 RFact = rcp( new GenericRFactory() );
328 }
329
330 RCP<RAPFactory> AcFact = rcp( new RAPFactory() );
331 AcFact->SetParameter("RepairMainDiagonal", Teuchos::ParameterEntry(bFixDiagonal));
332 for (size_t i = 0; i<TransferFacts_.size(); i++) {
333 AcFact->AddTransferFactory(TransferFacts_[i]);
334 }
335
336 //
337 // introduce rebalancing
338 //
339#if defined(HAVE_MUELU_ISORROPIA) && defined(HAVE_MPI)
340 Teuchos::RCP<Factory> RebalancedPFact = Teuchos::null;
341 Teuchos::RCP<Factory> RebalancedRFact = Teuchos::null;
342 Teuchos::RCP<Factory> RepartitionFact = Teuchos::null;
343 Teuchos::RCP<RebalanceAcFactory> RebalancedAFact = Teuchos::null;
344
345 MUELU_READ_PARAM(paramList, "repartition: enable", int, 0, bDoRepartition);
346 if (bDoRepartition == 1) {
347 // The Factory Manager will be configured to return the rebalanced versions of P, R, A by default.
348 // Everytime we want to use the non-rebalanced versions, we need to explicitly define the generating factory.
349 RFact->SetFactory("P", PFact);
350 //
351 AcFact->SetFactory("P", PFact);
352 AcFact->SetFactory("R", RFact);
353
354 // define rebalancing factory for coarse matrix
355 Teuchos::RCP<MueLu::AmalgamationFactory<SC, LO, GO, NO> > rebAmalgFact = Teuchos::rcp(new MueLu::AmalgamationFactory<SC, LO, GO, NO>());
356 rebAmalgFact->SetFactory("A", AcFact);
357
358 MUELU_READ_PARAM(paramList, "repartition: max min ratio", double, 1.3, maxminratio);
359 MUELU_READ_PARAM(paramList, "repartition: min per proc", int, 512, minperproc);
360
361 // Repartitioning heuristic
362 RCP<RepartitionHeuristicFactory> RepartitionHeuristicFact = Teuchos::rcp(new RepartitionHeuristicFactory());
363 {
364 Teuchos::ParameterList paramListRepFact;
365 paramListRepFact.set("repartition: min rows per proc", minperproc);
366 paramListRepFact.set("repartition: max imbalance", maxminratio);
367 RepartitionHeuristicFact->SetParameterList(paramListRepFact);
368 }
369 RepartitionHeuristicFact->SetFactory("A", AcFact);
370
371 // create "Partition"
372 Teuchos::RCP<MueLu::IsorropiaInterface<LO, GO, NO> > isoInterface = Teuchos::rcp(new MueLu::IsorropiaInterface<LO, GO, NO>());
373 isoInterface->SetFactory("A", AcFact);
374 isoInterface->SetFactory("number of partitions", RepartitionHeuristicFact);
375 isoInterface->SetFactory("UnAmalgamationInfo", rebAmalgFact);
376
377 // create "Partition" by unamalgamtion
378 Teuchos::RCP<MueLu::RepartitionInterface<LO, GO, NO> > repInterface = Teuchos::rcp(new MueLu::RepartitionInterface<LO, GO, NO>());
379 repInterface->SetFactory("A", AcFact);
380 repInterface->SetFactory("number of partitions", RepartitionHeuristicFact);
381 repInterface->SetFactory("AmalgamatedPartition", isoInterface);
382 //repInterface->SetFactory("UnAmalgamationInfo", rebAmalgFact); // not necessary?
383
384 // Repartitioning (creates "Importer" from "Partition")
385 RepartitionFact = Teuchos::rcp(new RepartitionFactory());
386 RepartitionFact->SetFactory("A", AcFact);
387 RepartitionFact->SetFactory("number of partitions", RepartitionHeuristicFact);
388 RepartitionFact->SetFactory("Partition", repInterface);
389
390 // Reordering of the transfer operators
391 RebalancedPFact = Teuchos::rcp(new RebalanceTransferFactory());
392 RebalancedPFact->SetParameter("type", Teuchos::ParameterEntry(std::string("Interpolation")));
393 RebalancedPFact->SetFactory("P", PFact);
394 RebalancedPFact->SetFactory("Nullspace", PtentFact);
395 RebalancedPFact->SetFactory("Importer", RepartitionFact);
396
397 RebalancedRFact = Teuchos::rcp(new RebalanceTransferFactory());
398 RebalancedRFact->SetParameter("type", Teuchos::ParameterEntry(std::string("Restriction")));
399 RebalancedRFact->SetFactory("R", RFact);
400 RebalancedRFact->SetFactory("Importer", RepartitionFact);
401
402 // Compute Ac from rebalanced P and R
403 RebalancedAFact = Teuchos::rcp(new RebalanceAcFactory());
404 RebalancedAFact->SetFactory("A", AcFact);
405 }
406#else // #ifdef HAVE_MUELU_ISORROPIA
407 // Get rid of [-Wunused] warnings
408 //(void)
409 //
410 // ^^^ FIXME (mfh 17 Nov 2013) That definitely doesn't compile.
411#endif
412
413 //
414 // Nullspace factory
415 //
416
417 // Set fine level nullspace
418 // extract pre-computed nullspace from ML parameter list
419 // store it in nullspace_ and nullspaceDim_
420 if (nullspaceType != "default vectors") {
421 TEUCHOS_TEST_FOR_EXCEPTION(nullspaceType != "pre-computed", Exceptions::RuntimeError, "MueLu::MLParameterListInterpreter: no valid nullspace (no pre-computed null space). error.");
422 TEUCHOS_TEST_FOR_EXCEPTION(nullspaceDim == -1, Exceptions::RuntimeError, "MueLu::MLParameterListInterpreter: no valid nullspace (nullspace dim == -1). error.");
423 TEUCHOS_TEST_FOR_EXCEPTION(nullspaceVec == NULL, Exceptions::RuntimeError, "MueLu::MLParameterListInterpreter: no valid nullspace (nullspace == NULL). You have to provide a valid fine-level nullspace in \'null space: vectors\'");
424
425 nullspaceDim_ = nullspaceDim;
426 nullspace_ = nullspaceVec;
427 }
428
429 Teuchos::RCP<NullspaceFactory> nspFact = Teuchos::rcp(new NullspaceFactory("Nullspace"));
430 nspFact->SetFactory("Nullspace", PtentFact);
431
432
433 // Stash coordinates
434 xcoord_ = xcoord;
435 ycoord_ = ycoord;
436 zcoord_ = zcoord;
437
438
439
440 //
441 // Hierarchy + FactoryManager
442 //
443
444 // Hierarchy options
445 this->numDesiredLevel_ = maxLevels;
446 this->maxCoarseSize_ = maxCoarseSize;
447
448 //
449 // Coarse Smoother
450 //
451 ParameterList& coarseList = paramList.sublist("coarse: list");
452 // check whether coarse solver is set properly. If not, set default coarse solver.
453 if (!coarseList.isParameter("smoother: type"))
454 coarseList.set("smoother: type", "Amesos-KLU"); // set default coarse solver according to ML 5.0 guide
455 RCP<SmootherFactory> coarseFact = GetSmootherFactory(coarseList, Teuchos::null);
456
457 // Smoothers Top Level Parameters
458
459 RCP<ParameterList> topLevelSmootherParam = ExtractSetOfParameters(paramList, "smoother");
460
461 //
462
463 // Prepare factory managers
464 // TODO: smootherFact can be reuse accross level if same parameters/no specific parameterList
465
466 for (int levelID=0; levelID < maxLevels; levelID++) {
467
468 //
469 // Level FactoryManager
470 //
471
472 RCP<FactoryManager> manager = rcp(new FactoryManager());
473 if (setKokkosRefactor)
474 manager->SetKokkosRefactor(useKokkosRefactor);
475
476 //
477 // Smoothers
478 //
479
480 {
481 // Merge level-specific parameters with global parameters. level-specific parameters takes precedence.
482 // TODO: unit-test this part alone
483
484 ParameterList levelSmootherParam = GetMLSubList(paramList, "smoother", levelID); // copy
485 MergeParameterList(*topLevelSmootherParam, levelSmootherParam, false); /* false = do no overwrite levelSmootherParam parameters by topLevelSmootherParam parameters */
486 // std::cout << std::endl << "Merged List for level " << levelID << std::endl;
487 // std::cout << levelSmootherParam << std::endl;
488
489 RCP<SmootherFactory> smootherFact = GetSmootherFactory(levelSmootherParam, Teuchos::null); // TODO: missing AFact input arg.
490
491 manager->SetFactory("Smoother", smootherFact);
492 }
493
494 //
495 // Misc
496 //
497
498 manager->SetFactory("CoarseSolver", coarseFact); // TODO: should not be done in the loop
499 manager->SetFactory("Graph", dropFact);
500 manager->SetFactory("Aggregates", AggFact);
501 manager->SetFactory("DofsPerNode", dropFact);
502 manager->SetFactory("Ptent", PtentFact);
503
504#if defined(HAVE_MUELU_ISORROPIA) && defined(HAVE_MPI)
505 if (bDoRepartition == 1) {
506 manager->SetFactory("A", RebalancedAFact);
507 manager->SetFactory("P", RebalancedPFact);
508 manager->SetFactory("R", RebalancedRFact);
509 manager->SetFactory("Nullspace", RebalancedPFact);
510 manager->SetFactory("Importer", RepartitionFact);
511 } else {
512#endif // #ifdef HAVE_MUELU_ISORROPIA
513 manager->SetFactory("Nullspace", nspFact); // use same nullspace factory throughout all multigrid levels
514 manager->SetFactory("A", AcFact); // same RAP factory for all levels
515 manager->SetFactory("P", PFact); // same prolongator and restrictor factories for all levels
516 manager->SetFactory("R", RFact); // same prolongator and restrictor factories for all levels
517#if defined(HAVE_MUELU_ISORROPIA) && defined(HAVE_MPI)
518 }
519#endif
520
521 this->AddFactoryManager(levelID, 1, manager);
522 } // for (level loop)
523
524 }
525
526 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
528 // if nullspace_ has already been extracted from ML parameter list
529 // make nullspace available for MueLu
530 if (nullspace_ != NULL) {
531 RCP<Level> fineLevel = H.GetLevel(0);
532 RCP<Operator> Op = fineLevel->Get<RCP<Operator> >("A");
533 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Op);
534 if (!A.is_null()) {
535 const RCP<const Map> rowMap = fineLevel->Get< RCP<Matrix> >("A")->getRowMap();
536 RCP<MultiVector> nullspace = MultiVectorFactory::Build(rowMap, nullspaceDim_, true);
537
538 for ( size_t i=0; i < Teuchos::as<size_t>(nullspaceDim_); i++) {
539 Teuchos::ArrayRCP<Scalar> nullspacei = nullspace->getDataNonConst(i);
540 const size_t myLength = nullspace->getLocalLength();
541
542 for (size_t j = 0; j < myLength; j++) {
543 nullspacei[j] = nullspace_[i*myLength + j];
544 }
545 }
546
547 fineLevel->Set("Nullspace", nullspace);
548 }
549 }
550
551 // Do the same for coordinates
552 size_t num_coords = 0;
553 double * coordPTR[3];
554 if (xcoord_) {
555 coordPTR[0] = xcoord_;
556 num_coords++;
557 if (ycoord_) {
558 coordPTR[1] = ycoord_;
559 num_coords++;
560 if (zcoord_) {
561 coordPTR[2] = zcoord_;
562 num_coords++;
563 }
564 }
565 }
566 if (num_coords){
567 Teuchos::RCP<Level> fineLevel = H.GetLevel(0);
568 Teuchos::RCP<Operator> Op = fineLevel->Get<RCP<Operator> >("A");
569 Teuchos::RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Op);
570 if (!A.is_null()) {
571 const Teuchos::RCP<const Map> rowMap = fineLevel->Get< RCP<Matrix> >("A")->getRowMap();
572 Teuchos::RCP<MultiVector> coordinates = MultiVectorFactory::Build(rowMap, num_coords, true);
573
574 for ( size_t i=0; i < num_coords; i++) {
575 Teuchos::ArrayRCP<Scalar> coordsi = coordinates->getDataNonConst(i);
576 const size_t myLength = coordinates->getLocalLength();
577 for (size_t j = 0; j < myLength; j++) {
578 coordsi[j] = coordPTR[i][j];
579 }
580 }
581 fineLevel->Set("Coordinates",coordinates);
582 }
583 }
584
586 }
587
588 // TODO: code factorization with MueLu_ParameterListInterpreter.
589 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
590 RCP<MueLu::SmootherFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
592 GetSmootherFactory (const Teuchos::ParameterList & paramList,
593 const RCP<FactoryBase> & AFact)
594 {
595 typedef Teuchos::ScalarTraits<Scalar> STS;
596 SC one = STS::one();
597
598 std::string type = "symmetric Gauss-Seidel"; // default
599
600 //
601 // Get 'type'
602 //
603
604// //TODO: fix defaults!!
605
606// // Default coarse grid smoother
607// std::string type;
608// if ("smoother" == "coarse") {
609// #if (defined(HAVE_MUELU_EPETRA) && defined( HAVE_MUELU_AMESOS)) || (defined(HAVE_MUELU_AMESOS2)) // FIXME: test is wrong (ex: compiled with Epetra&&Tpetra&&Amesos2 but without Amesos => error running Epetra problem)
610// type = ""; // use default defined by AmesosSmoother or Amesos2Smoother
611// #else
612// type = "symmetric Gauss-Seidel"; // use a sym Gauss-Seidel (with no damping) as fallback "coarse solver" (TODO: needs Ifpack(2))
613// #endif
614// } else {
615// // TODO: default smoother?
616// type = "";
617// }
618
619
620 if (paramList.isParameter("smoother: type")) type = paramList.get<std::string>("smoother: type");
621 TEUCHOS_TEST_FOR_EXCEPTION(type.empty(), Exceptions::RuntimeError, "MueLu::MLParameterListInterpreter: no \"smoother: type\" in the smoother parameter list" << std::endl << paramList);
622
623 //
624 // Create the smoother prototype
625 //
626
627 RCP<SmootherPrototype> smooProto;
628 std::string ifpackType;
629 Teuchos::ParameterList smootherParamList;
630
631 if (type == "Jacobi" || type == "Gauss-Seidel" || type == "symmetric Gauss-Seidel") {
632 if (type == "symmetric Gauss-Seidel") type = "Symmetric Gauss-Seidel"; // FIXME
633
634 ifpackType = "RELAXATION";
635 smootherParamList.set("relaxation: type", type);
636
637 MUELU_COPY_PARAM(paramList, "smoother: sweeps", int, 2, smootherParamList, "relaxation: sweeps");
638 MUELU_COPY_PARAM(paramList, "smoother: damping factor", Scalar, one, smootherParamList, "relaxation: damping factor");
639
640 smooProto = rcp( new TrilinosSmoother(ifpackType, smootherParamList, 0) );
641 smooProto->SetFactory("A", AFact);
642
643 } else if (type == "Chebyshev" || type == "MLS") {
644
645 ifpackType = "CHEBYSHEV";
646
647 MUELU_COPY_PARAM(paramList, "smoother: sweeps", int, 2, smootherParamList, "chebyshev: degree");
648 if (paramList.isParameter("smoother: MLS alpha")) {
649 MUELU_COPY_PARAM(paramList, "smoother: MLS alpha", double, 20, smootherParamList, "chebyshev: ratio eigenvalue");
650 } else {
651 MUELU_COPY_PARAM(paramList, "smoother: Chebyshev alpha", double, 20, smootherParamList, "chebyshev: ratio eigenvalue");
652 }
653
654
655 smooProto = rcp( new TrilinosSmoother(ifpackType, smootherParamList, 0) );
656 smooProto->SetFactory("A", AFact);
657
658 } else if (type == "Hiptmair") {
659 ifpackType = "HIPTMAIR";
660 std::string subSmootherType = "Chebyshev";
661 if (paramList.isParameter("subsmoother: type"))
662 subSmootherType = paramList.get<std::string>("subsmoother: type");
663 std::string subSmootherIfpackType;
664 if (subSmootherType == "Chebyshev")
665 subSmootherIfpackType = "CHEBYSHEV";
666 else if (subSmootherType == "Jacobi" || subSmootherType == "Gauss-Seidel" || subSmootherType == "symmetric Gauss-Seidel") {
667 if (subSmootherType == "symmetric Gauss-Seidel") subSmootherType = "Symmetric Gauss-Seidel"; // FIXME
668 subSmootherIfpackType = "RELAXATION";
669 } else
670 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "MueLu::MLParameterListInterpreter: unknown smoother type. '" << subSmootherType << "' not supported by MueLu.");
671
672 smootherParamList.set("hiptmair: smoother type 1", subSmootherIfpackType);
673 smootherParamList.set("hiptmair: smoother type 2", subSmootherIfpackType);
674
675 auto smoother1ParamList = smootherParamList.sublist("hiptmair: smoother list 1");
676 auto smoother2ParamList = smootherParamList.sublist("hiptmair: smoother list 2");
677
678 if (subSmootherType == "Chebyshev") {
679 MUELU_COPY_PARAM(paramList, "subsmoother: edge sweeps", int, 2, smoother1ParamList, "chebyshev: degree");
680 MUELU_COPY_PARAM(paramList, "subsmoother: node sweeps", int, 2, smoother2ParamList, "chebyshev: degree");
681
682 MUELU_COPY_PARAM(paramList, "subsmoother: Chebyshev", double, 20, smoother1ParamList, "chebyshev: ratio eigenvalue");
683 MUELU_COPY_PARAM(paramList, "subsmoother: Chebyshev", double, 20, smoother2ParamList, "chebyshev: ratio eigenvalue");
684 } else {
685 MUELU_COPY_PARAM(paramList, "subsmoother: edge sweeps", int, 2, smoother1ParamList, "relaxation: sweeps");
686 MUELU_COPY_PARAM(paramList, "subsmoother: node sweeps", int, 2, smoother2ParamList, "relaxation: sweeps");
687
688 MUELU_COPY_PARAM(paramList, "subsmoother: SGS damping factor", double, 0.8, smoother2ParamList, "relaxation: damping factor");
689 }
690
691
692 smooProto = rcp( new TrilinosSmoother(ifpackType, smootherParamList, 0) );
693 smooProto->SetFactory("A", AFact);
694
695 } else if (type == "IFPACK") { // TODO: this option is not described in the ML Guide v5.0
696
697#if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_IFPACK)
698 ifpackType = paramList.get<std::string>("smoother: ifpack type");
699
700 if (ifpackType == "ILU") {
701 // TODO fix this (type mismatch double vs. int)
702 //MUELU_COPY_PARAM(paramList, "smoother: ifpack level-of-fill", double /*int*/, 0.0 /*2*/, smootherParamList, "fact: level-of-fill");
703 if (paramList.isParameter("smoother: ifpack level-of-fill"))
704 smootherParamList.set("fact: level-of-fill", Teuchos::as<int>(paramList.get<double>("smoother: ifpack level-of-fill")));
705 else smootherParamList.set("fact: level-of-fill", as<int>(0));
706
707 MUELU_COPY_PARAM(paramList, "smoother: ifpack overlap", int, 2, smootherParamList, "partitioner: overlap");
708
709 // TODO change to TrilinosSmoother as soon as Ifpack2 supports all preconditioners from Ifpack
710 smooProto =
712 smootherParamList,
713 paramList.get<int> ("smoother: ifpack overlap"));
714 smooProto->SetFactory("A", AFact);
715 } else {
716 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "MueLu::MLParameterListInterpreter: unknown ML smoother type " + type + " (IFPACK) not supported by MueLu. Only ILU is supported.");
717 }
718#else
719 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "MueLu::MLParameterListInterpreter: MueLu compiled without Ifpack support");
720#endif
721
722 } else if (type.length() > strlen("Amesos") && type.substr(0, strlen("Amesos")) == "Amesos") { /* catch Amesos-* */
723 std::string solverType = type.substr(strlen("Amesos")+1); /* ("Amesos-KLU" -> "KLU") */
724
725 // Validator: following upper/lower case is what is allowed by ML
726 bool valid = false;
727 const int validatorSize = 5;
728 std::string validator[validatorSize] = {"Superlu", "Superludist", "KLU", "UMFPACK", "MUMPS"}; /* TODO: should "" be allowed? */
729 for (int i=0; i < validatorSize; i++) { if (validator[i] == solverType) valid = true; }
730 TEUCHOS_TEST_FOR_EXCEPTION(!valid, Exceptions::RuntimeError, "MueLu::MLParameterListInterpreter: unknown smoother type. '" << type << "' not supported.");
731
732 // FIXME: MueLu should accept any Upper/Lower case. Not the case for the moment
733 std::transform(solverType.begin()+1, solverType.end(), solverType.begin()+1, ::tolower);
734
735 smooProto = Teuchos::rcp( new DirectSolver(solverType, Teuchos::ParameterList()) );
736 smooProto->SetFactory("A", AFact);
737
738 } else {
739
740 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "MueLu::MLParameterListInterpreter: unknown smoother type. '" << type << "' not supported by MueLu.");
741
742 }
743 TEUCHOS_TEST_FOR_EXCEPTION(smooProto == Teuchos::null, Exceptions::RuntimeError, "MueLu::MLParameterListInterpreter: no smoother prototype. fatal error.");
744
745 //
746 // Create the smoother factory
747 //
748
749 RCP<SmootherFactory> SmooFact = rcp( new SmootherFactory() );
750
751 // Set parameters of the smoother factory
752 MUELU_READ_PARAM(paramList, "smoother: pre or post", std::string, "both", preOrPost);
753 if (preOrPost == "both") {
754 SmooFact->SetSmootherPrototypes(smooProto, smooProto);
755 } else if (preOrPost == "pre") {
756 SmooFact->SetSmootherPrototypes(smooProto, Teuchos::null);
757 } else if (preOrPost == "post") {
758 SmooFact->SetSmootherPrototypes(Teuchos::null, smooProto);
759 }
760
761 return SmooFact;
762 }
763
764 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
766 // check if it's a TwoLevelFactoryBase based transfer factory
767 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!");
768 TransferFacts_.push_back(factory);
769 }
770
771 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
775
776 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
778 try {
779 Matrix& A = dynamic_cast<Matrix&>(Op);
780 if (A.IsFixedBlockSizeSet() && (A.GetFixedBlockSize() != blksize_))
781 this->GetOStream(Warnings0) << "Setting matrix block size to " << blksize_ << " (value of the parameter in the list) "
782 << "instead of " << A.GetFixedBlockSize() << " (provided matrix)." << std::endl;
783
784 A.SetFixedBlockSize(blksize_);
785
786#ifdef HAVE_MUELU_DEBUG
787 MatrixUtils::checkLocalRowMapMatchesColMap(A);
788#endif // HAVE_MUELU_DEBUG
789
790 } catch (std::bad_cast&) {
791 this->GetOStream(Warnings0) << "Skipping setting block size as the operator is not a matrix" << std::endl;
792 }
793 }
794
795} // namespace MueLu
796
797#define MUELU_MLPARAMETERLISTINTERPRETER_SHORT
798#endif /* MUELU_MLPARAMETERLISTINTERPRETER_DEF_HPP */
799
800//TODO: see if it can be factorized with ML interpreter (ex: generation of Ifpack param list)
#define MUELU_COPY_PARAM(paramList, paramStr, varType, defaultValue, outParamList, outParamStr)
#define MUELU_READ_PARAM(paramList, paramStr, varType, defaultValue, varName)
MueLu::DefaultScalar Scalar
MueLu::DefaultNode Node
AmalgamationFactory for subblocks of strided map based amalgamation data.
Factory for creating a graph based on a given matrix.
Factory for creating a graph based on a given matrix.
Exception indicating invalid cast attempted.
Exception throws to report errors in the internal logical of the program.
This class specifies the default factory that should generate some data on a Level if the data does n...
Factory for building restriction operators using a prolongator factory.
virtual void SetupHierarchy(Hierarchy &H) const
Setup Hierarchy object.
void AddFactoryManager(int startLevel, int numDesiredLevel, RCP< FactoryManagerBase > manager)
RCP< Level > & GetLevel(const int levelID=0)
Retrieve a certain level from hierarchy.
void AddTransferFactory(const RCP< FactoryBase > &factory)
Add transfer factory in the end of list of transfer factories for RAPFactory.
double * xcoord_
coordinates can be embedded in the ML parameter list
void SetParameterList(const Teuchos::ParameterList &paramList)
static RCP< SmootherFactory > GetSmootherFactory(const Teuchos::ParameterList &paramList, const RCP< FactoryBase > &AFact=Teuchos::null)
Read smoother options and build the corresponding smoother factory.
size_t NumTransferFactories() const
Returns number of transfer factories.
MueLu::DirectSolver< Scalar, LocalOrdinal, GlobalOrdinal, Node > DirectSolver
virtual void SetupHierarchy(Hierarchy &H) const
Setup Hierarchy object.
MueLu::RepartitionHeuristicFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node > RepartitionHeuristicFactory
MueLu::RebalanceTransferFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node > RebalanceTransferFactory
int nullspaceDim_
nullspace can be embedded in the ML parameter list
MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > Hierarchy
MueLu::RebalanceAcFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node > RebalanceAcFactory
MueLu::RepartitionFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node > RepartitionFactory
virtual void SetupOperator(Operator &Op) const
Setup Operator object.
Factory for generating nullspace.
Factory for building Petrov-Galerkin Smoothed Aggregation prolongators.
Factory for building coarse matrices.
Helper class which transforms an "AmalgamatedPartition" array to an unamalgamated "Partition".
Factory for building Smoothed Aggregation prolongators.
Factory for building Smoothed Aggregation prolongators.
Generic Smoother Factory for generating the smoothers of the MG hierarchy.
Factory for building tentative prolongator.
Factory for building restriction operators.
Class that encapsulates external library smoothers.
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
Namespace for MueLu classes and methods.
Teuchos::RCP< Teuchos::ParameterList > ExtractSetOfParameters(const Teuchos::ParameterList &paramList, const std::string &str)
@ Warnings0
Important warning messages (one line).
@ Runtime1
Description of what is happening (more verbose).
void CreateSublists(const ParameterList &List, ParameterList &newList)
RCP< MueLu::SmootherPrototype< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetIfpackSmoother(const std::string &="", const Teuchos::ParameterList &=Teuchos::ParameterList(), const LocalOrdinal &=0)
void MergeParameterList(const Teuchos::ParameterList &source, Teuchos::ParameterList &dest, bool overWrite)
: merge two parameter lists
const Teuchos::ParameterList & GetMLSubList(const Teuchos::ParameterList &paramList, const std::string &type, int levelID)