FlowLinearSolverParameters.hpp
Go to the documentation of this file.
1/*
2 Copyright 2015, 2020 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2015 IRIS AS
4 Copyright 2015 Dr. Blatt - HPC-Simulation-Software & Services
5 Copyright 2015 NTNU
6 Copyright 2015 Statoil AS
7
8 This file is part of the Open Porous Media project (OPM).
9
10 OPM is free software: you can redistribute it and/or modify
11 it under the terms of the GNU General Public License as published by
12 the Free Software Foundation, either version 3 of the License, or
13 (at your option) any later version.
14
15 OPM is distributed in the hope that it will be useful,
16 but WITHOUT ANY WARRANTY; without even the implied warranty of
17 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 GNU General Public License for more details.
19
20 You should have received a copy of the GNU General Public License
21 along with OPM. If not, see <http://www.gnu.org/licenses/>.
22*/
23
24#ifndef OPM_FLOWLINEARSOLVERPARAMETERS_HEADER_INCLUDED
25#define OPM_FLOWLINEARSOLVERPARAMETERS_HEADER_INCLUDED
26
28
29#include <opm/simulators/linalg/linalgproperties.hh>
30#include <opm/models/utils/parametersystem.hh>
31
32namespace Opm {
33template <class TypeTag>
34class ISTLSolverBda;
35template <class TypeTag>
36class ISTLSolver;
37}
38
39
40
41namespace Opm::Properties {
42
43namespace TTag {
45}
46
47template<class TypeTag, class MyTypeTag>
49 using type = UndefinedProperty;
50};
51
52template<class TypeTag, class MyTypeTag>
54 using type = UndefinedProperty;
55};
56
57template<class TypeTag, class MyTypeTag>
59 using type = UndefinedProperty;
60};
61template<class TypeTag, class MyTypeTag>
63 using type = UndefinedProperty;
64};
65//
66// LinearSolverVerbosity defined in opm-models
67//
68template<class TypeTag, class MyTypeTag>
70 using type = UndefinedProperty;
71};
72template<class TypeTag, class MyTypeTag>
74 using type = UndefinedProperty;
75};
76template<class TypeTag, class MyTypeTag>
78 using type = UndefinedProperty;
79};
80template<class TypeTag, class MyTypeTag>
82 using type = UndefinedProperty;
83};
84template<class TypeTag, class MyTypeTag>
86 using type = UndefinedProperty;
87};
88template<class TypeTag, class MyTypeTag>
89struct UseGmres {
90 using type = UndefinedProperty;
91};
92template<class TypeTag, class MyTypeTag>
94 using type = UndefinedProperty;
95};
96template<class TypeTag, class MyTypeTag>
98 using type = UndefinedProperty;
99};
100template<class TypeTag, class MyTypeTag>
102 using type = UndefinedProperty;
103};
104template<class TypeTag, class MyTypeTag>
106 using type = UndefinedProperty;
107};
108template<class TypeTag, class MyTypeTag>
110 using type = UndefinedProperty;
111};
112template<class TypeTag, class MyTypeTag>
114 using type = UndefinedProperty;
115};
116template<class TypeTag, class MyTypeTag>
118 using type = UndefinedProperty;
119};
120template<class TypeTag, class MyTypeTag>
122 using type = UndefinedProperty;
123};
124template<class TypeTag, class MyTypeTag>
126 using type = UndefinedProperty;
127};
128template<class TypeTag, class MyTypeTag>
130 using type = UndefinedProperty;
131};
132template<class TypeTag>
133struct LinearSolverReduction<TypeTag, TTag::FlowIstlSolverParams> {
134 using type = GetPropType<TypeTag, Scalar>;
135 static constexpr type value = 1e-2;
136};
137template<class TypeTag>
138struct RelaxedLinearSolverReduction<TypeTag, TTag::FlowIstlSolverParams> {
139 using type = GetPropType<TypeTag, Scalar>;
140 static constexpr type value = 1e-2;
141};
142template<class TypeTag>
143struct LinearSolverMaxIter<TypeTag, TTag::FlowIstlSolverParams> {
144 static constexpr int value = 200;
145};
146template<class TypeTag>
147struct LinearSolverRestart<TypeTag, TTag::FlowIstlSolverParams> {
148 static constexpr int value = 40;
149};
150template<class TypeTag>
151struct LinearSolverVerbosity<TypeTag, TTag::FlowIstlSolverParams> {
152 static constexpr int value = 0;
153};
154template<class TypeTag>
155struct IluRelaxation<TypeTag, TTag::FlowIstlSolverParams> {
156 using type = GetPropType<TypeTag, Scalar>;
157 static constexpr type value = 0.9;
158};
159template<class TypeTag>
160struct IluFillinLevel<TypeTag, TTag::FlowIstlSolverParams> {
161 static constexpr int value = 0;
162};
163template<class TypeTag>
164struct MiluVariant<TypeTag, TTag::FlowIstlSolverParams> {
165 static constexpr auto value = "ILU";
166};
167template<class TypeTag>
168struct IluRedblack<TypeTag, TTag::FlowIstlSolverParams> {
169 static constexpr bool value = false;
170};
171template<class TypeTag>
172struct IluReorderSpheres<TypeTag, TTag::FlowIstlSolverParams> {
173 static constexpr bool value = false;
174};
175template<class TypeTag>
176struct UseGmres<TypeTag, TTag::FlowIstlSolverParams> {
177 static constexpr bool value = false;
178};
179template<class TypeTag>
180struct LinearSolverIgnoreConvergenceFailure<TypeTag, TTag::FlowIstlSolverParams> {
181 static constexpr bool value = false;
182};
183template<class TypeTag>
184struct ScaleLinearSystem<TypeTag, TTag::FlowIstlSolverParams> {
185 static constexpr bool value = false;
186};
187template<class TypeTag>
188struct LinearSolver<TypeTag, TTag::FlowIstlSolverParams> {
189 static constexpr auto value = "ilu0";
190};
191template<class TypeTag>
192struct LinearSolverPrintJsonDefinition<TypeTag, TTag::FlowIstlSolverParams> {
193 static constexpr auto value = true;
194};
195template<class TypeTag>
196struct CprReuseSetup<TypeTag, TTag::FlowIstlSolverParams> {
197 static constexpr int value = 4;
198};
199template<class TypeTag>
200struct CprReuseInterval<TypeTag, TTag::FlowIstlSolverParams> {
201 static constexpr int value = 30;
202};
203template<class TypeTag>
204struct AcceleratorMode<TypeTag, TTag::FlowIstlSolverParams> {
205 static constexpr auto value = "none";
206};
207template<class TypeTag>
208struct BdaDeviceId<TypeTag, TTag::FlowIstlSolverParams> {
209 static constexpr int value = 0;
210};
211template<class TypeTag>
212struct OpenclPlatformId<TypeTag, TTag::FlowIstlSolverParams> {
213 static constexpr int value = 0;
214};
215template<class TypeTag>
216struct OpenclIluParallel<TypeTag, TTag::FlowIstlSolverParams> {
217 static constexpr bool value = true; // note: false should only be used in debug
218};
219
220// Set the backend to be used.
221template<class TypeTag>
222struct LinearSolverBackend<TypeTag, TTag::FlowIstlSolverParams> {
223#if COMPILE_BDA_BRIDGE
225#else
227#endif
228};
229} // namespace Opm::Properties
230
231namespace Opm
232{
233
234
237 {
251 std::string linsolver_;
255 std::string accelerator_mode_;
259
260 template <class TypeTag>
261 void init(bool cprRequestedInDataFile)
262 {
263 // TODO: these parameters have undocumented non-trivial dependencies
264 linear_solver_reduction_ = Parameters::get<TypeTag, Properties::LinearSolverReduction>();
265 relaxed_linear_solver_reduction_ = Parameters::get<TypeTag, Properties::RelaxedLinearSolverReduction>();
266 linear_solver_maxiter_ = Parameters::get<TypeTag, Properties::LinearSolverMaxIter>();
267 linear_solver_restart_ = Parameters::get<TypeTag, Properties::LinearSolverRestart>();
268 linear_solver_verbosity_ = Parameters::get<TypeTag, Properties::LinearSolverVerbosity>();
269 ilu_relaxation_ = Parameters::get<TypeTag, Properties::IluRelaxation>();
270 ilu_fillin_level_ = Parameters::get<TypeTag, Properties::IluFillinLevel>();
271 ilu_milu_ = convertString2Milu(Parameters::get<TypeTag, Properties::MiluVariant>());
272 ilu_redblack_ = Parameters::get<TypeTag, Properties::IluRedblack>();
273 ilu_reorder_sphere_ = Parameters::get<TypeTag, Properties::IluReorderSpheres>();
274 newton_use_gmres_ = Parameters::get<TypeTag, Properties::UseGmres>();
275 ignoreConvergenceFailure_ = Parameters::get<TypeTag, Properties::LinearSolverIgnoreConvergenceFailure>();
276 scale_linear_system_ = Parameters::get<TypeTag, Properties::ScaleLinearSystem>();
277 linsolver_ = Parameters::get<TypeTag, Properties::LinearSolver>();
278 linear_solver_print_json_definition_ = Parameters::get<TypeTag, Properties::LinearSolverPrintJsonDefinition>();
279 cpr_reuse_setup_ = Parameters::get<TypeTag, Properties::CprReuseSetup>();
280 cpr_reuse_interval_ = Parameters::get<TypeTag, Properties::CprReuseInterval>();
281
282 if (!Parameters::isSet<TypeTag, Properties::LinearSolver>() && cprRequestedInDataFile) {
283 linsolver_ = "cpr";
284 } else {
285 linsolver_ = Parameters::get<TypeTag, Properties::LinearSolver>();
286 }
287
288 accelerator_mode_ = Parameters::get<TypeTag, Properties::AcceleratorMode>();
289 bda_device_id_ = Parameters::get<TypeTag, Properties::BdaDeviceId>();
290 opencl_platform_id_ = Parameters::get<TypeTag, Properties::OpenclPlatformId>();
291 opencl_ilu_parallel_ = Parameters::get<TypeTag, Properties::OpenclIluParallel>();
292 }
293
294 template <class TypeTag>
295 static void registerParameters()
296 {
297 Parameters::registerParam<TypeTag, Properties::LinearSolverReduction>
298 ("The minimum reduction of the residual which the linear solver must achieve");
299 Parameters::registerParam<TypeTag, Properties::RelaxedLinearSolverReduction>
300 ("The minimum reduction of the residual which the linear solver need to "
301 "achieve for the solution to be accepted");
302 Parameters::registerParam<TypeTag, Properties::LinearSolverMaxIter>
303 ("The maximum number of iterations of the linear solver");
304 Parameters::registerParam<TypeTag, Properties::LinearSolverRestart>
305 ("The number of iterations after which GMRES is restarted");
306 Parameters::registerParam<TypeTag, Properties::LinearSolverVerbosity>
307 ("The verbosity level of the linear solver (0: off, 2: all)");
308 Parameters::registerParam<TypeTag, Properties::IluRelaxation>
309 ("The relaxation factor of the linear solver's ILU preconditioner");
310 Parameters::registerParam<TypeTag, Properties::IluFillinLevel>
311 ("The fill-in level of the linear solver's ILU preconditioner");
312 Parameters::registerParam<TypeTag, Properties::MiluVariant>
313 ("Specify which variant of the modified-ILU preconditioner ought to be used. "
314 "Possible variants are: ILU (default, plain ILU), "
315 "MILU_1 (lump diagonal with dropped row entries), "
316 "MILU_2 (lump diagonal with the sum of the absolute values of the dropped row entries), "
317 "MILU_3 (if diagonal is positive add sum of dropped row entries, otherwise subtract them), "
318 "MILU_4 (if diagonal is positive add sum of dropped row entries, otherwise do nothing");
319 Parameters::registerParam<TypeTag, Properties::IluRedblack>
320 ("Use red-black partitioning for the ILU preconditioner");
321 Parameters::registerParam<TypeTag, Properties::IluReorderSpheres>
322 ("Whether to reorder the entries of the matrix in the red-black "
323 "ILU preconditioner in spheres starting at an edge. "
324 "If false the original ordering is preserved in each color. "
325 "Otherwise why try to ensure D4 ordering (in a 2D structured grid, "
326 "the diagonal elements are consecutive).");
327 Parameters::registerParam<TypeTag, Properties::UseGmres>
328 ("Use GMRES as the linear solver");
329 Parameters::registerParam<TypeTag, Properties::LinearSolverIgnoreConvergenceFailure>
330 ("Continue with the simulation like nothing happened "
331 "after the linear solver did not converge");
332 Parameters::registerParam<TypeTag, Properties::ScaleLinearSystem>
333 ("Scale linear system according to equation scale and primary variable types");
334 Parameters::registerParam<TypeTag, Properties::LinearSolver>
335 ("Configuration of solver. Valid options are: ilu0 (default), "
336 "dilu, cprw, cpr (an alias for cprw), cpr_quasiimpes, "
337 "cpr_trueimpes, cpr_trueimpesanalytic, amg or hybrid (experimental). "
338 "Alternatively, you can request a configuration to be read from a "
339 "JSON file by giving the filename here, ending with '.json.'");
340 Parameters::registerParam<TypeTag, Properties::LinearSolverPrintJsonDefinition>
341 ("Write the JSON definition of the linear solver setup to the DBG file.");
342 Parameters::registerParam<TypeTag, Properties::CprReuseSetup>
343 ("Reuse preconditioner setup. Valid options are "
344 "0: recreate the preconditioner for every linear solve, "
345 "1: recreate once every timestep, "
346 "2: recreate if last linear solve took more than 10 iterations, "
347 "3: never recreate, "
348 "4: recreated every CprReuseInterval");
349 Parameters::registerParam<TypeTag, Properties::CprReuseInterval>
350 ("Reuse preconditioner interval. Used when CprReuseSetup is set to 4, "
351 "then the preconditioner will be fully recreated instead of reused "
352 "every N linear solve, where N is this parameter.");
353 Parameters::registerParam<TypeTag, Properties::AcceleratorMode>
354 ("Choose a linear solver, usage: "
355 "'--accelerator-mode=[none|cusparse|opencl|amgcl|rocalution|rocsparse]'");
356 Parameters::registerParam<TypeTag, Properties::BdaDeviceId>
357 ("Choose device ID for cusparseSolver or openclSolver, "
358 "use 'nvidia-smi' or 'clinfo' to determine valid IDs");
359 Parameters::registerParam<TypeTag, Properties::OpenclPlatformId>
360 ("Choose platform ID for openclSolver, use 'clinfo' "
361 "to determine valid platform IDs");
362 Parameters::registerParam<TypeTag, Properties::OpenclIluParallel>
363 ("Parallelize ILU decomposition and application on GPU");
364 }
365
367
368 // set default values
369 void reset()
370 {
376 ilu_relaxation_ = 0.9;
379 ilu_redblack_ = false;
380 ilu_reorder_sphere_ = false;
381 newton_use_gmres_ = false;
383 scale_linear_system_ = false;
384 linsolver_ = "ilu0";
388 accelerator_mode_ = "none";
389 bda_device_id_ = 0;
392 }
393 };
394
395
396} // namespace Opm
397
398
399
400
401#endif // OPM_FLOWLINEARSOLVERPARAMETERS_HEADER_INCLUDED
Definition: ISTLSolverBda.hpp:102
Definition: ISTLSolver.hpp:143
Definition: AluGridVanguard.hpp:57
Definition: BlackoilPhases.hpp:27
MILU_VARIANT
Definition: MILU.hpp:34
@ ILU
Do not perform modified ILU.
MILU_VARIANT convertString2Milu(const std::string &milu)
This class carries all parameters for the NewtonIterationBlackoilInterleaved class.
Definition: FlowLinearSolverParameters.hpp:237
bool ilu_redblack_
Definition: FlowLinearSolverParameters.hpp:246
int cpr_reuse_interval_
Definition: FlowLinearSolverParameters.hpp:254
int ilu_fillin_level_
Definition: FlowLinearSolverParameters.hpp:244
bool scale_linear_system_
Definition: FlowLinearSolverParameters.hpp:250
void init(bool cprRequestedInDataFile)
Definition: FlowLinearSolverParameters.hpp:261
double linear_solver_reduction_
Definition: FlowLinearSolverParameters.hpp:238
std::string accelerator_mode_
Definition: FlowLinearSolverParameters.hpp:255
bool newton_use_gmres_
Definition: FlowLinearSolverParameters.hpp:248
bool ignoreConvergenceFailure_
Definition: FlowLinearSolverParameters.hpp:249
int cpr_reuse_setup_
Definition: FlowLinearSolverParameters.hpp:253
int linear_solver_restart_
Definition: FlowLinearSolverParameters.hpp:241
double relaxed_linear_solver_reduction_
Definition: FlowLinearSolverParameters.hpp:239
int linear_solver_maxiter_
Definition: FlowLinearSolverParameters.hpp:240
int linear_solver_verbosity_
Definition: FlowLinearSolverParameters.hpp:242
int opencl_platform_id_
Definition: FlowLinearSolverParameters.hpp:257
FlowLinearSolverParameters()
Definition: FlowLinearSolverParameters.hpp:366
double ilu_relaxation_
Definition: FlowLinearSolverParameters.hpp:243
bool opencl_ilu_parallel_
Definition: FlowLinearSolverParameters.hpp:258
bool linear_solver_print_json_definition_
Definition: FlowLinearSolverParameters.hpp:252
bool ilu_reorder_sphere_
Definition: FlowLinearSolverParameters.hpp:247
static void registerParameters()
Definition: FlowLinearSolverParameters.hpp:295
void reset()
Definition: FlowLinearSolverParameters.hpp:369
int bda_device_id_
Definition: FlowLinearSolverParameters.hpp:256
MILU_VARIANT ilu_milu_
Definition: FlowLinearSolverParameters.hpp:245
std::string linsolver_
Definition: FlowLinearSolverParameters.hpp:251
Definition: FlowLinearSolverParameters.hpp:117
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:118
Definition: FlowLinearSolverParameters.hpp:121
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:122
Definition: FlowLinearSolverParameters.hpp:113
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:114
Definition: FlowLinearSolverParameters.hpp:109
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:110
Definition: FlowLinearSolverParameters.hpp:73
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:74
Definition: FlowLinearSolverParameters.hpp:81
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:82
GetPropType< TypeTag, Scalar > type
Definition: FlowLinearSolverParameters.hpp:156
Definition: FlowLinearSolverParameters.hpp:69
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:70
Definition: FlowLinearSolverParameters.hpp:85
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:86
Definition: FlowLinearSolverParameters.hpp:93
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:94
Definition: FlowLinearSolverParameters.hpp:58
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:59
Definition: FlowLinearSolverParameters.hpp:101
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:102
Definition: FlowLinearSolverParameters.hpp:105
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:106
GetPropType< TypeTag, Scalar > type
Definition: FlowLinearSolverParameters.hpp:134
Definition: FlowLinearSolverParameters.hpp:48
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:49
Definition: FlowLinearSolverParameters.hpp:62
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:63
Definition: FlowLinearSolverParameters.hpp:77
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:78
Definition: FlowLinearSolverParameters.hpp:129
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:130
Definition: FlowLinearSolverParameters.hpp:125
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:126
GetPropType< TypeTag, Scalar > type
Definition: FlowLinearSolverParameters.hpp:139
Definition: FlowLinearSolverParameters.hpp:53
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:54
Definition: FlowLinearSolverParameters.hpp:97
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:98
Definition: FlowLinearSolverParameters.hpp:44
Definition: FlowLinearSolverParameters.hpp:89
UndefinedProperty type
Definition: FlowLinearSolverParameters.hpp:90