13#ifndef DUMUX_LINEAR_SOLVER_PARAMETERS_HH
14#define DUMUX_LINEAR_SOLVER_PARAMETERS_HH
20#include <dune/common/parametertree.hh>
21#include <dune/common/std/type_traits.hh>
29constexpr bool hasGridView = Dune::Std::is_detected<GVDetector, T>::value;
38template<
class LinearSolverTraits>
49 Dune::ParameterTree params;
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);
69 params[
"preconditioner.defaultAggregationDimension"] =
"3";
70 params[
"preconditioner.maxLevel"] =
"100";
71 params[
"ParameterGroup"] = paramGroup;
72 params[
"preconditioner.ParameterGroup"] = paramGroup;
79 if (linearSolverGroups.empty())
84 for (
const auto& group : linearSolverGroups)
86 const auto fullDumuxKey = group +
"." + dumuxKey;
87 const auto value = getParam<std::string>(fullDumuxKey,
"");
90 params[istlKey] = value;
100template<
class LinearSolverTraits>
101const std::vector<std::array<std::string, 2>>
105 {
"Verbosity",
"verbose"},
106 {
"MaxIterations",
"maxit"},
107 {
"ResidualReduction",
"reduction"},
109 {
"GMResRestart",
"restart"},
110 {
"Restart",
"restart"},
111 {
"MaxOrthogonalizationVectors",
"mmax"},
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"}
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 ¶ms, const std::string ¶mGroup="")
Set some defaults for the solver parameters.
Definition: linearsolverparameters.hh:56
static void fillValuesForIstlKeys(Dune::ParameterTree ¶ms, const std::string ¶mGroup="")
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 ¶mGroup="")
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 ¶mGroup)
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
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.