version 3.8
linearsolverparameters.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_LINEAR_SOLVER_PARAMETERS_HH
14#define DUMUX_LINEAR_SOLVER_PARAMETERS_HH
15
16#include <string>
17#include <array>
18#include <vector>
19
20#include <dune/common/parametertree.hh>
21#include <dune/common/std/type_traits.hh>
23
25template <typename T>
26using GVDetector = typename T::GridView;
27
28template <typename T>
29constexpr bool hasGridView = Dune::Std::is_detected<GVDetector, T>::value;
30} // end namespace Dumux::Detail::LinearSolverParameters
31
32namespace Dumux {
33
38template<class LinearSolverTraits>
40{
41public:
42
44 static const std::vector<std::array<std::string, 2>> dumuxToIstlSolverParams;
45
47 static Dune::ParameterTree createParameterTree(const std::string& paramGroup = "")
48 {
49 Dune::ParameterTree params;
50 setDefaultParameters(params, paramGroup);
51 fillValuesForIstlKeys(params, paramGroup);
52 return params;
53 }
54
56 static void setDefaultParameters(Dune::ParameterTree& params, const std::string& paramGroup = "")
57 {
58 params["restart"] = "10";
59 params["maxit"] = "250";
60 params["reduction"] = "1e-13";
61 params["verbose"] = "0";
62 params["preconditioner.iterations"] = "1";
63 params["preconditioner.relaxation"] = "1.0";
64 params["preconditioner.verbosity"] = "0";
65 params["preconditioner.defaultAggregationSizeMode"] = "isotropic";
66 if constexpr (Detail::LinearSolverParameters::hasGridView<LinearSolverTraits>)
67 params["preconditioner.defaultAggregationDimension"] = std::to_string(LinearSolverTraits::GridView::dimension);
68 else
69 params["preconditioner.defaultAggregationDimension"] = "3";
70 params["preconditioner.maxLevel"] = "100";
71 params["ParameterGroup"] = paramGroup;
72 params["preconditioner.ParameterGroup"] = paramGroup;
73 }
74
76 static void fillValuesForIstlKeys(Dune::ParameterTree& params, const std::string& paramGroup = "")
77 {
78 const auto linearSolverGroups = getParamSubGroups("LinearSolver", paramGroup);
79 if (linearSolverGroups.empty()) // no linear solver parameters were specified
80 return;
81
82 for (const auto& [dumuxKey, istlKey] : dumuxToIstlSolverParams)
83 {
84 for (const auto& group : linearSolverGroups)
85 {
86 const auto fullDumuxKey = group + "." + dumuxKey;
87 const auto value = getParam<std::string>(fullDumuxKey, "");
88 if (!value.empty())
89 {
90 params[istlKey] = value;
91 break; // skip groups with smaller depth in the tree
92 }
93 }
94 }
95 }
96};
97
100template<class LinearSolverTraits>
101const std::vector<std::array<std::string, 2>>
103{
104 // solver params
105 {"Verbosity", "verbose"},
106 {"MaxIterations", "maxit"},
107 {"ResidualReduction", "reduction"},
108 {"Type", "type"},
109 {"GMResRestart", "restart"}, // cycles before restarting
110 {"Restart", "restart"}, // cycles before restarting
111 {"MaxOrthogonalizationVectors", "mmax"},
112
113 // preconditioner params
114 {"Preconditioner.Verbosity", "preconditioner.verbosity"},
115 {"Preconditioner.Type", "preconditioner.type"},
116 {"Preconditioner.Iterations", "preconditioner.iterations"},
117 {"Preconditioner.Relaxation", "preconditioner.relaxation"},
118 {"Preconditioner.ILUOrder", "preconditioner.n"},
119 {"Preconditioner.ILUResort", "preconditioner.resort"},
120 {"Preconditioner.AmgSmootherRelaxation", "preconditioner.smootherRelaxation"},
121 {"Preconditioner.AmgSmootherIterations", "preconditioner.smootherIterations"},
122 {"Preconditioner.AmgMaxLevel", "preconditioner.maxLevel"},
123 {"Preconditioner.AmgCoarsenTarget", "preconditioner.coarsenTarget"},
124 {"Preconditioner.AmgMinCoarseningRate", "preconditioner.minCoarseningRate"},
125 {"Preconditioner.AmgAccumulationMode", "preconditioner.accumulationMode"},
126 {"Preconditioner.AmgProlongationDampingFactor", "preconditioner.prolongationDampingFactor"},
127 {"Preconditioner.AmgAlpha", "preconditioner.alpha"},
128 {"Preconditioner.AmgBeta", "preconditioner.beta"},
129 {"Preconditioner.AmgAdditive", "preconditioner.additive"},
130 {"Preconditioner.AmgGamma", "preconditioner.gamma"},
131 {"Preconditioner.AmgPreSmoothingSteps", "preconditioner.preSteps"},
132 {"Preconditioner.AmgPostSmoothingSteps", "preconditioner.postSteps"},
133 {"Preconditioner.AmgCriterionSymmetric", "preconditioner.criterionSymmetric"},
134 {"Preconditioner.AmgStrengthMeasure", "preconditioner.strengthMeasure"},
135 {"Preconditioner.AmgDiagonalRowIndex", "preconditioner.diagonalRowIndex"},
136 {"Preconditioner.AmgDefaultAggregationSizeMode", "preconditioner.defaultAggregationSizeMode"},
137 {"Preconditioner.AmgDefaultAggregationDimension", "preconditioner.defaultAggregationDimension"},
138 {"Preconditioner.AmgMaxAggregateDistance", "preconditioner.maxAggregateDistance"},
139 {"Preconditioner.AmgMinAggregateSize", "preconditioner.minAggregateSize"},
140 {"Preconditioner.AmgMaxAggregateSize", "preconditioner.maxAggregateSize"}
141};
142
143} // end namespace Dumux
144
145#endif
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Definition: linearsolverparameters.hh:40
static const std::vector< std::array< std::string, 2 > > dumuxToIstlSolverParams
Translation table for solver parameters.
Definition: linearsolverparameters.hh:44
static void setDefaultParameters(Dune::ParameterTree &params, const std::string &paramGroup="")
Set some defaults for the solver parameters.
Definition: linearsolverparameters.hh:56
static void fillValuesForIstlKeys(Dune::ParameterTree &params, const std::string &paramGroup="")
Iterate over all keys required by the ISTL, translate them to Dumux syntax and add values to tree.
Definition: linearsolverparameters.hh:76
static Dune::ParameterTree createParameterTree(const std::string &paramGroup="")
Create a tree containing parameters required for the linear solvers and precondioners of the Dune IST...
Definition: linearsolverparameters.hh:47
std::vector< std::string > getParamSubGroups(const std::string &subGroupName, const std::string &paramGroup)
Get a list of sub groups from the parameter tree sorted by relevance.
Definition: parameters.hh:173
Definition: linearsolverparameters.hh:24
typename T::GridView GVDetector
Definition: linearsolverparameters.hh:26
constexpr bool hasGridView
Definition: linearsolverparameters.hh:29
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.