PRISMS-PF Manual
Loading...
Searching...
No Matches
linear_solve_parameters.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: © 2025 PRISMS Center at the University of Michigan
2// SPDX-License-Identifier: GNU Lesser General Public Version 2.1
3
4#pragma once
5
6#include <deal.II/base/exceptions.h>
7#include <deal.II/base/parameter_handler.h>
8#include <deal.II/lac/precondition.h>
9#include <deal.II/lac/solver_selector.h>
10
13#include <prismspf/core/types.h>
14
15#include <prismspf/config.h>
16
17#include <map>
18
20
25{
26 // Solver type. richardson|cg|bicgstab|gmres|fgmres|minres
27 std::string solver_type = "cg";
28
29 // Solver tolerance
31
32 // Solver tolerance type
34
35 // Max number of iterations for the linear solve
37
38 // Preconditioner
40
41 dealii::PreconditionChebyshev<>::AdditionalData chebyshev_parameters;
42
43 dealii::SolverRichardson<>::AdditionalData richardson_parameters;
44 dealii::SolverBicgstab<>::AdditionalData bicgstab_parameters;
45 dealii::SolverGMRES<>::AdditionalData gmres_parameters;
46
47 // MinRes and CG do not have additional parameters, so we do not need to store them here
48 // dealii::SolverMinRes<>::AdditionalData minres_parameters;
49 // dealii::SolverCG<>::AdditionalData cg_parameters;
50 // FGMRES additional parameters are a subset of GMRES parameters, so we can just use
51 // gmres_parameters for both solvers
52 // dealii::SolverFGMRES<>::AdditionalData fgmres_parameters;
53
54 // The minimum multigrid level
55 unsigned int min_mg_level = 0;
56};
57
62{
66 void
67 validate();
68
72 void
74
78 void
79 declare_parameters(dealii::ParameterHandler &parameter_handler,
80 unsigned int max_criteria = 5) const;
81
85 void
86 assign_parameters(dealii::ParameterHandler &parameter_handler,
87 unsigned int max_criteria = 5);
88
89 // Map of linear solve parameters for fields that require them
90 std::map<unsigned int, LinearSolverParameters> linear_solvers;
91};
92
93inline void
95{
96 // Nothing to do here for now
97}
98
99inline void
101{
102 if (!linear_solvers.empty())
103 {
105 << "================================================\n"
106 << " Linear Solve Parameters\n"
107 << "================================================\n";
108
109 for (const auto &[index, linear_solver_parameters] : linear_solvers)
110 {
112 << "Index: " << index << "\n"
113 << " Tolerance: " << linear_solver_parameters.tolerance << "\n"
114 << " Type: " << to_string(linear_solver_parameters.tolerance_type) << "\n"
115 << " Max iterations: " << linear_solver_parameters.max_iterations << "\n"
116 << " Preconditioner: " << to_string(linear_solver_parameters.preconditioner)
117 << "\n";
118
119 if (linear_solver_parameters.preconditioner == PreconditionerType::GMG)
120 {
121 // ConditionalOStreams::pout_summary();
122 }
123 }
124
125 ConditionalOStreams::pout_summary() << "\n" << std::flush;
126 }
127}
128
129inline void
130LinearSolveParameters::declare_parameters(dealii::ParameterHandler &parameter_handler,
131 unsigned int max_criteria) const
132{
133 for (unsigned int criterion_id = 0; criterion_id < max_criteria; criterion_id++)
134 {
135 // For linear solves
136 std::string subsection_text =
137 "linear solver parameters: " + std::to_string(criterion_id);
138 parameter_handler.enter_subsection(subsection_text);
139 {
140 parameter_handler.declare_entry(
141 "solver_ids",
142 "",
143 dealii::Patterns::Anything(),
144 "The ids of the solvers that will use these settings.");
145 parameter_handler.declare_entry(
146 "solver type",
147 "cg",
148 dealii::Patterns::Selection(dealii::SolverSelector<>::get_solver_names()),
149 "The type of iterative solver to use for linear solves.");
150 parameter_handler.declare_entry(
151 "tolerance type",
152 "AbsoluteResidual",
153 dealii::Patterns::Selection("AbsoluteResidual|RMSEPerField|IntegratedPerField|"
154 "RMSETotal|IntegratedTotal"),
155 "The tolerance type for the linear solver.");
156 parameter_handler.declare_entry("tolerance value",
157 "1.0e-10",
158 dealii::Patterns::Double(DBL_MIN, DBL_MAX),
159 "The value of for the linear solver tolerance.");
160 parameter_handler.declare_entry(
161 "max iterations",
162 "100",
163 dealii::Patterns::Integer(1, INT_MAX),
164 "The maximum number of linear solver iterations before the loop "
165 "is stopped.");
166 parameter_handler.declare_entry("preconditioner type",
167 "None",
168 dealii::Patterns::Selection("None|Chebyshev"),
169 "The preconditioner type for the linear solver.");
170 parameter_handler.enter_subsection("Chebyshev");
171 {
172 parameter_handler.declare_entry("smoothing range",
173 "15.0",
174 dealii::Patterns::Double(DBL_MIN, DBL_MAX),
175 "The smoothing range for eigenvalues.");
176 parameter_handler.declare_entry("smoother degree",
177 "5",
178 dealii::Patterns::Integer(1, INT_MAX),
179 "The smoother polynomial degree.");
180 parameter_handler.declare_entry(
181 "eigenvalue cg iterations",
182 "10",
183 dealii::Patterns::Integer(1, INT_MAX),
184 "The maximum number of CG iterations used to find the maximum eigenvalue.");
185 }
186 parameter_handler.leave_subsection();
187
188 parameter_handler.enter_subsection("Richardson");
189 {
190 parameter_handler.declare_entry("omega",
191 "1.0",
192 dealii::Patterns::Double(),
193 "Damping factor.");
194 parameter_handler.declare_entry("use preconditioned residual",
195 "false",
196 dealii::Patterns::Bool(),
197 "Whether to use the preconditioned residual l2 "
198 "norm in the stopping criterion.");
199 }
200 parameter_handler.leave_subsection();
201 parameter_handler.enter_subsection("BiCGStab");
202 {
203 parameter_handler.declare_entry("exact residual",
204 "true",
205 dealii::Patterns::Bool(),
206 "Flag for exact computation of residual.");
207 parameter_handler.declare_entry("breakdown",
208 std::to_string(
209 std::numeric_limits<double>::min()),
210 dealii::Patterns::Double(),
211 "Breakdown threshold.");
212 }
213 parameter_handler.leave_subsection();
214 parameter_handler.enter_subsection("GMRES");
215 {
216 parameter_handler.declare_entry(
217 "max basis size",
218 "30",
219 dealii::Patterns::Integer(1, INT_MAX),
220 "The maximum size of the Krylov basis used in GMRES before restarting.");
221 parameter_handler.declare_entry(
222 "orthogonalization strategy",
223 "delayed_classical_gram_schmidt",
224 dealii::Patterns::Selection("classical_gram_schmidt|modified_gram_schmidt|"
225 "delayed_classical_gram_schmidt"),
226 "The orthogonalization strategy to use in GMRES.");
227 parameter_handler.declare_entry("right preconditioning",
228 "false",
229 dealii::Patterns::Bool(),
230 "Whether to use right preconditioning.");
231 parameter_handler.declare_entry(
232 "use default residual",
233 "true",
234 dealii::Patterns::Bool(),
235 "Whether to use the default residual computation in GMRES.");
236 parameter_handler.declare_entry(
237 "force re-orthogonalization",
238 "false",
239 dealii::Patterns::Bool(),
240 "Whether to force re-orthogonalization of the Krylov basis in GMRES.");
241 parameter_handler.declare_entry("batched mode",
242 "false",
243 dealii::Patterns::Bool(),
244 "Whether to use batched mode in GMRES.");
245 }
246 parameter_handler.leave_subsection();
247
248 parameter_handler.declare_entry("min mg level",
249 "0",
250 dealii::Patterns::Integer(0, INT_MAX),
251 "The minimum multigrid level.");
252 parameter_handler.declare_alias("tolerance value", "tolerance");
253 parameter_handler.declare_alias("solver_ids", "solve blocks");
254 parameter_handler.declare_alias("solver_ids", "solve_blocks");
255 parameter_handler.declare_alias("solver_ids", "solve block ids");
256 parameter_handler.declare_alias("solver_ids", "solve_block_ids");
257 parameter_handler.declare_alias("solver type", "solver_type");
258 parameter_handler.declare_alias("solver type", "linear solver type");
259 parameter_handler.declare_alias("solver type", "linear_solver_type");
260 parameter_handler.declare_alias("solver type", "linear solver");
261 parameter_handler.declare_alias("solver type", "linear_solver");
262 parameter_handler.declare_alias("solver type", "type");
263 parameter_handler.declare_alias("preconditioner type", "preconditioner_type");
264 parameter_handler.declare_alias("preconditioner type", "preconditioner");
265 }
266 parameter_handler.leave_subsection();
267 }
268}
269
270inline void
271LinearSolveParameters::assign_parameters(dealii::ParameterHandler &parameter_handler,
272 unsigned int max_criteria)
273{
274 for (unsigned int criterion_id = 0; criterion_id < max_criteria; criterion_id++)
275 {
276 // For linear solves
277 std::string subsection_text =
278 "linear solver parameters: " + std::to_string(criterion_id);
279 parameter_handler.enter_subsection(subsection_text);
280 {
281 std::vector<int> solver_ids = dealii::Utilities::string_to_int(
282 dealii::Utilities::split_string_list(parameter_handler.get("solver_ids")));
283
284 LinearSolverParameters linear_solver_parameters;
285
286 // Set the linear solver type
287 linear_solver_parameters.solver_type = parameter_handler.get("solver type");
288
289 // Set the tolerance type
290 static const std::map<std::string, SolverToleranceType> tolerance_types = {
291 {"AbsoluteResidual", AbsoluteResidual },
292 {"RMSEPerField", RMSEPerField },
293 {"IntegratedPerField", IntegratedPerField},
294 {"RMSETotal", RMSETotal },
295 {"IntegratedTotal", IntegratedTotal },
296 };
297 const std::string type_string = parameter_handler.get("tolerance type");
298 const auto &iter = tolerance_types.find(type_string);
299 if (iter != tolerance_types.end())
300 {
301 linear_solver_parameters.tolerance_type = iter->second;
302 }
303 else
304 {
305 AssertThrow(false,
306 dealii::ExcMessage("Invalid tolerance type: " + type_string));
307 }
308
309 // Set the tolerance value
310 linear_solver_parameters.tolerance =
311 parameter_handler.get_double("tolerance value");
312
313 // Set the maximum number of iterations
314 linear_solver_parameters.max_iterations =
315 static_cast<unsigned int>(parameter_handler.get_integer("max iterations"));
316
317 // Set preconditioner type and related parameters
318 static const std::map<std::string, PreconditionerType> preconditioner_map = {
319 {"None", None },
320 {"none", None },
321 {"", None },
322 {"Chebyshev", Chebyshev},
323 {"chebyshev", Chebyshev},
324 {"GMG", GMG },
325 {"gmg", GMG }
326 };
327 linear_solver_parameters.preconditioner =
328 preconditioner_map.at(parameter_handler.get("preconditioner type"));
329
330 parameter_handler.enter_subsection("Chebyshev");
331 {
332 linear_solver_parameters.chebyshev_parameters.degree =
333 static_cast<unsigned int>(parameter_handler.get_integer("smoother degree"));
334
335 linear_solver_parameters.chebyshev_parameters.smoothing_range =
336 parameter_handler.get_double("smoothing range");
337
338 linear_solver_parameters.chebyshev_parameters.eig_cg_n_iterations =
339 static_cast<unsigned int>(
340 parameter_handler.get_integer("eigenvalue cg iterations"));
341 }
342 parameter_handler.leave_subsection();
343 parameter_handler.enter_subsection("Richardson");
344 {
345 linear_solver_parameters.richardson_parameters.omega =
346 parameter_handler.get_double("omega");
347 linear_solver_parameters.richardson_parameters.use_preconditioned_residual =
348 parameter_handler.get_bool("use preconditioned residual");
349 }
350 parameter_handler.leave_subsection();
351 parameter_handler.enter_subsection("BiCGStab");
352 {
353 linear_solver_parameters.bicgstab_parameters.exact_residual =
354 parameter_handler.get_bool("exact residual");
355 linear_solver_parameters.bicgstab_parameters.breakdown =
356 parameter_handler.get_double("breakdown");
357 }
358 parameter_handler.leave_subsection();
359 parameter_handler.enter_subsection("GMRES");
360 {
361 using OStrat = dealii::LinearAlgebra::OrthogonalizationStrategy;
362 static const std::map<std::string, OStrat> orthogonalization_strategy_names = {
363 {"classical_gram_schmidt", OStrat::classical_gram_schmidt },
364 {"modified_gram_schmidt", OStrat::modified_gram_schmidt },
365 {"delayed_classical_gram_schmidt", OStrat::delayed_classical_gram_schmidt}
366 };
367 linear_solver_parameters.gmres_parameters.orthogonalization_strategy =
368 orthogonalization_strategy_names.at(
369 parameter_handler.get("orthogonalization strategy"));
370 linear_solver_parameters.gmres_parameters.right_preconditioning =
371 parameter_handler.get_bool("right preconditioning");
372 linear_solver_parameters.gmres_parameters.use_default_residual =
373 parameter_handler.get_bool("use default residual");
374 linear_solver_parameters.gmres_parameters.force_re_orthogonalization =
375 parameter_handler.get_bool("force re-orthogonalization");
376 linear_solver_parameters.gmres_parameters.batched_mode =
377 parameter_handler.get_bool("batched mode");
378 }
379 parameter_handler.leave_subsection();
380
381 linear_solver_parameters.min_mg_level =
382 static_cast<unsigned int>(parameter_handler.get_integer("min mg level"));
383 for (auto solver_id : solver_ids)
384 {
385 linear_solvers[static_cast<unsigned int>(solver_id)] =
386 linear_solver_parameters;
387 }
388 }
389 parameter_handler.leave_subsection();
390 }
391}
392
393PRISMS_PF_END_NAMESPACE
static dealii::ConditionalOStream & pout_summary()
Log output stream for writing a summary.log file.
Definition conditional_ostreams.cc:35
static const unsigned int iterations
Default iterations.
Definition types.h:67
static const double tolerance
Default tolerance.
Definition types.h:57
Definition conditional_ostreams.cc:20
Struct that holds linear solver parameters.
Definition linear_solve_parameters.h:62
void validate()
Postprocess and validate parameters.
Definition linear_solve_parameters.h:94
std::map< unsigned int, LinearSolverParameters > linear_solvers
Definition linear_solve_parameters.h:90
void assign_parameters(dealii::ParameterHandler &parameter_handler, unsigned int max_criteria=5)
Assign the parameters read from an input file to this object.
Definition linear_solve_parameters.h:271
void print_parameter_summary() const
Print parameters to summary.log.
Definition linear_solve_parameters.h:100
void declare_parameters(dealii::ParameterHandler &parameter_handler, unsigned int max_criteria=5) const
Declare the parameters to be read from an input file.
Definition linear_solve_parameters.h:130
Struct that stores relevant linear solve information of a certain solve block.
Definition linear_solve_parameters.h:25
unsigned int max_iterations
Definition linear_solve_parameters.h:36
unsigned int min_mg_level
Definition linear_solve_parameters.h:55
double tolerance
Definition linear_solve_parameters.h:30
dealii::SolverGMRES ::AdditionalData gmres_parameters
Definition linear_solve_parameters.h:45
SolverToleranceType tolerance_type
Definition linear_solve_parameters.h:33
dealii::PreconditionChebyshev ::AdditionalData chebyshev_parameters
Definition linear_solve_parameters.h:41
PreconditionerType preconditioner
Definition linear_solve_parameters.h:39
std::string solver_type
Definition linear_solve_parameters.h:27
dealii::SolverBicgstab ::AdditionalData bicgstab_parameters
Definition linear_solve_parameters.h:44
dealii::SolverRichardson ::AdditionalData richardson_parameters
Definition linear_solve_parameters.h:43
SolverToleranceType
Solver tolerance type.
Definition type_enums.h:112
@ RMSEPerField
The mean local error averaged over each field is lower than the tolerance.
Definition type_enums.h:124
@ IntegratedTotal
The sum of the integrated errors of each field is lower than the tolerance.
Definition type_enums.h:136
@ AbsoluteResidual
Legacy.
Definition type_enums.h:116
@ RMSETotal
The sum of the average local errors of each field is lower than the tolerance.
Definition type_enums.h:132
@ IntegratedPerField
The integrated error averaged over each field is lower than the tolerance.
Definition type_enums.h:128
std::string to_string(ElasticityModel type)
Enum to string for ElasticityModel.
Definition type_enums.h:164
PreconditionerType
Preconditioner type.
Definition type_enums.h:143
@ Chebyshev
Definition type_enums.h:145
@ None
Definition type_enums.h:144
@ GMG
Definition type_enums.h:146