3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_LINEAR_SOLVER_PARAMETERS_HH
26#define DUMUX_LINEAR_SOLVER_PARAMETERS_HH
27
28#include <string>
29#include <array>
30#include <vector>
31
32#include <dune/common/parametertree.hh>
34
35namespace Dumux {
40template<class LinearSolverTraits>
42{
43public:
44
46 static const std::vector<std::array<std::string, 2>> dumuxToIstlSolverParams;
47
49 static Dune::ParameterTree createParameterTree(const std::string& paramGroup = "")
50 {
51 Dune::ParameterTree params;
52 setDefaultParameters(params, paramGroup);
53 fillValuesForIstlKeys(params, paramGroup);
54 return params;
55 }
56
58 static void setDefaultParameters(Dune::ParameterTree& params, const std::string& paramGroup = "")
59 {
60 params["restart"] = "10";
61 params["maxit"] = "250";
62 params["reduction"] = "1e-13";
63 params["verbose"] = "0";
64 params["preconditioner.iterations"] = "1";
65 params["preconditioner.relaxation"] = "1.0";
66 params["preconditioner.verbosity"] = "0";
67 params["preconditioner.defaultAggregationSizeMode"] = "isotropic";
68 params["preconditioner.defaultAggregationDimension"] = std::to_string(LinearSolverTraits::GridView::dimension);
69 params["preconditioner.maxLevel"] = "100";
70 params["ParameterGroup"] = paramGroup;
71 params["preconditioner.ParameterGroup"] = paramGroup;
72 }
73
75 static void fillValuesForIstlKeys(Dune::ParameterTree& params, const std::string& paramGroup = "")
76 {
77 const auto linearSolverGroups = getParamSubGroups("LinearSolver", paramGroup);
78 if (linearSolverGroups.empty()) // no linear solver parameters were specified
79 return;
80
81 for (const auto& [dumuxKey, istlKey] : dumuxToIstlSolverParams)
82 {
83 for (const auto& group : linearSolverGroups)
84 {
85 const auto fullDumuxKey = group + "." + dumuxKey;
86 const auto value = getParam<std::string>(fullDumuxKey, "");
87 if (!value.empty())
88 {
89 params[istlKey] = value;
90 break; // skip groups with smaller depth in the tree
91 }
92 }
93 }
94 }
95};
96
99template<class LinearSolverTraits>
100const std::vector<std::array<std::string, 2>>
102{
103 // solver params
104 {"Verbosity", "verbose"},
105 {"MaxIterations", "maxit"},
106 {"ResidualReduction", "reduction"},
107 {"Type", "type"},
108 {"GMResRestart", "restart"}, // cycles before restarting
109 {"Restart", "restart"}, // cycles before restarting
110 {"MaxOrthogonalizationVectors", "mmax"},
111
112 // preconditioner params
113 {"Preconditioner.Verbosity", "preconditioner.verbosity"},
114 {"Preconditioner.Type", "preconditioner.type"},
115 {"Preconditioner.Iterations", "preconditioner.iterations"},
116 {"Preconditioner.Relaxation", "preconditioner.relaxation"},
117 {"Preconditioner.ILUOrder", "preconditioner.n"},
118 {"Preconditioner.ILUResort", "preconditioner.resort"},
119 {"Preconditioner.AmgSmootherRelaxation", "preconditioner.smootherRelaxation"},
120 {"Preconditioner.AmgSmootherIterations", "preconditioner.smootherIterations"},
121 {"Preconditioner.AmgMaxLevel", "preconditioner.maxLevel"},
122 {"Preconditioner.AmgCoarsenTarget", "preconditioner.coarsenTarget"},
123 {"Preconditioner.AmgMinCoarseningRate", "preconditioner.minCoarseningRate"},
124 {"Preconditioner.AmgAccumulationMode", "preconditioner.accumulationMode"},
125 {"Preconditioner.AmgProlongationDampingFactor", "preconditioner.prolongationDampingFactor"},
126 {"Preconditioner.AmgAlpha", "preconditioner.alpha"},
127 {"Preconditioner.AmgBeta", "preconditioner.beta"},
128 {"Preconditioner.AmgAdditive", "preconditioner.additive"},
129 {"Preconditioner.AmgGamma", "preconditioner.gamma"},
130 {"Preconditioner.AmgPreSmoothingSteps", "preconditioner.preSteps"},
131 {"Preconditioner.AmgPostSmoothingSteps", "preconditioner.postSteps"},
132 {"Preconditioner.AmgCriterionSymmetric", "preconditioner.criterionSymmetric"},
133 {"Preconditioner.AmgStrengthMeasure", "preconditioner.strengthMeasure"},
134 {"Preconditioner.AmgDiagonalRowIndex", "preconditioner.diagonalRowIndex"},
135 {"Preconditioner.AmgDefaultAggregationSizeMode", "preconditioner.defaultAggregationSizeMode"},
136 {"Preconditioner.AmgDefaultAggregationDimension", "preconditioner.defaultAggregationDimension"},
137 {"Preconditioner.AmgMaxAggregateDistance", "preconditioner.maxAggregateDistance"},
138 {"Preconditioner.AmgMinAggregateSize", "preconditioner.minAggregateSize"},
139 {"Preconditioner.AmgMaxAggregateSize", "preconditioner.maxAggregateSize"}
140};
141
142} // end namespace Dumux
143
144#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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:185
Definition: adapt.hh:29
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Definition: linearsolverparameters.hh:42
static const std::vector< std::array< std::string, 2 > > dumuxToIstlSolverParams
Translation table for solver parameters.
Definition: linearsolverparameters.hh:46
static void setDefaultParameters(Dune::ParameterTree &params, const std::string &paramGroup="")
Set some defaults for the solver parameters.
Definition: linearsolverparameters.hh:58
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:75
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:49