ExaDG
Loading...
Searching...
No Matches
multigrid_parameters.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef INCLUDE_SOLVERS_AND_PRECONDITIONERS_MULTIGRIDINPUTPARAMETERS_H_
23#define INCLUDE_SOLVERS_AND_PRECONDITIONERS_MULTIGRIDINPUTPARAMETERS_H_
24
25// C/C++
26#include <string>
27#include <vector>
28
29// deal.II
30#include <deal.II/lac/petsc_precondition.h>
31#include <deal.II/lac/trilinos_precondition.h>
32
33// ExaDG
34#include <exadg/solvers_and_preconditioners/solvers/solver_data.h>
35#include <exadg/utilities/print_functions.h>
36
37namespace ExaDG
38{
39enum class MultigridType
40{
41 Undefined,
42 hMG,
43 chMG,
44 hcMG,
45 cMG,
46 pMG,
47 cpMG,
48 pcMG,
49 hpMG,
50 chpMG,
51 hcpMG,
52 hpcMG,
53 phMG,
54 cphMG,
55 pchMG,
56 phcMG
57};
58
59enum class PSequenceType
60{
61 GoToOne,
62 DecreaseByOne,
63 Bisect,
64 Manual
65};
66
67enum class MultigridSmoother
68{
69 Chebyshev,
70 GMRES,
71 CG,
72 Jacobi
73};
74
75enum class AMGType
76{
77 ML,
78 BoomerAMG
79};
80
81enum class MultigridCoarseGridSolver
82{
83 Chebyshev,
84 CG,
85 GMRES,
86 AMG
87};
88
89enum class MultigridCoarseGridPreconditioner
90{
91 None,
92 PointJacobi,
93 BlockJacobi,
94 AMG
95};
96
97struct AMGData
98{
99 AMGData()
100 {
101 amg_type = AMGType::ML;
102
103#ifdef DEAL_II_WITH_TRILINOS
104 ml_data.smoother_sweeps = 1;
105 ml_data.n_cycles = 1;
106 ml_data.smoother_type = "ILU";
107#endif
108
109#ifdef DEAL_II_WITH_PETSC
110 boomer_data.n_sweeps_coarse = 1;
111 boomer_data.max_iter = 1;
112 boomer_data.relaxation_type_down =
113 dealii::PETScWrappers::PreconditionBoomerAMG::AdditionalData::RelaxationType::Chebyshev;
114 boomer_data.relaxation_type_up =
115 dealii::PETScWrappers::PreconditionBoomerAMG::AdditionalData::RelaxationType::Chebyshev;
116 boomer_data.relaxation_type_coarse =
117 dealii::PETScWrappers::PreconditionBoomerAMG::AdditionalData::RelaxationType::Chebyshev;
118#endif
119 };
120
121 void
122 print(dealii::ConditionalOStream const & pcout) const
123 {
124 print_parameter(pcout, " AMG type", amg_type);
125
126 if(amg_type == AMGType::ML)
127 {
128#ifdef DEAL_II_WITH_TRILINOS
129 print_parameter(pcout, " Smoother sweeps", ml_data.smoother_sweeps);
130 print_parameter(pcout, " Number of cycles", ml_data.n_cycles);
131 print_parameter(pcout, " Smoother type", ml_data.smoother_type);
132#endif
133 }
134 else if(amg_type == AMGType::BoomerAMG)
135 {
136#ifdef DEAL_II_WITH_PETSC
137 print_parameter(pcout, " Smoother sweeps", boomer_data.n_sweeps_coarse);
138 print_parameter(pcout, " Number of cycles", boomer_data.max_iter);
139 print_parameter(pcout, " Smoother type down", boomer_data.relaxation_type_down);
140 print_parameter(pcout, " Smoother type up", boomer_data.relaxation_type_up);
141 print_parameter(pcout, " Smoother type coarse", boomer_data.relaxation_type_coarse);
142#endif
143 }
144 else
145 {
146 AssertThrow(false, dealii::ExcNotImplemented());
147 }
148 }
149
150 AMGType amg_type;
151
152#ifdef DEAL_II_WITH_TRILINOS
153 dealii::TrilinosWrappers::PreconditionAMG::AdditionalData ml_data;
154#endif
155
156#ifdef DEAL_II_WITH_PETSC
157 dealii::PETScWrappers::PreconditionBoomerAMG::AdditionalData boomer_data;
158#endif
159};
160
161enum class PreconditionerSmoother
162{
163 None,
164 PointJacobi,
165 BlockJacobi,
166 AdditiveSchwarz
167};
168
170{
172 : smoother(MultigridSmoother::Chebyshev),
173 preconditioner(PreconditionerSmoother::PointJacobi),
174 iterations(5),
175 relaxation_factor(0.8),
176 smoothing_range(20),
177 iterations_eigenvalue_estimation(20)
178 {
179 }
180
181 void
182 print(dealii::ConditionalOStream const & pcout) const
183 {
184 print_parameter(pcout, "Smoother", smoother);
185 print_parameter(pcout, "Preconditioner smoother", preconditioner);
186 print_parameter(pcout, "Iterations smoother", iterations);
187
188 if(smoother == MultigridSmoother::Jacobi)
189 {
190 print_parameter(pcout, "Relaxation factor", relaxation_factor);
191 }
192
193 if(smoother == MultigridSmoother::Chebyshev)
194 {
195 print_parameter(pcout, "Smoothing range", smoothing_range);
196 print_parameter(pcout, "Iterations eigenvalue estimation", iterations_eigenvalue_estimation);
197 }
198 }
199
200 // Type of smoother
201 MultigridSmoother smoother;
202
203 // Preconditioner used for smoother
204 PreconditionerSmoother preconditioner;
205
206 // Number of iterations
207 unsigned int iterations;
208
209 // damping/relaxation factor for Jacobi smoother
210 double relaxation_factor;
211
212 // Chebyshev smmother: sets the smoothing range (range of eigenvalues to be smoothed)
213 double smoothing_range;
214
215 // Chebyshev smmother: number of CG iterations for estimation of eigenvalues
216 unsigned int iterations_eigenvalue_estimation;
217};
218
220{
222 : solver(MultigridCoarseGridSolver::Chebyshev),
223 preconditioner(MultigridCoarseGridPreconditioner::PointJacobi),
224 solver_data(SolverData(1e4, 1.e-12, 1.e-3)),
225 amg_data(AMGData())
226 {
227 }
228
229 void
230 print(dealii::ConditionalOStream const & pcout) const
231 {
232 print_parameter(pcout, "Coarse grid solver", solver);
233 print_parameter(pcout, "Coarse grid preconditioner", preconditioner);
234
235 solver_data.print(pcout);
236
237 if(solver == MultigridCoarseGridSolver::AMG or
238 preconditioner == MultigridCoarseGridPreconditioner::AMG)
239 {
240 amg_data.print(pcout);
241 }
242 }
243
244 // Coarse grid solver
245 MultigridCoarseGridSolver solver;
246
247 // Coarse grid preconditioner
248 MultigridCoarseGridPreconditioner preconditioner;
249
250 // Solver data for coarse grid solver
251 SolverData solver_data;
252
253 // Configuration of AMG settings
254 AMGData amg_data;
255};
256
257
259{
261 : type(MultigridType::hMG),
262 p_sequence(PSequenceType::Bisect),
263 smoother_data(SmootherData()),
264 coarse_problem(CoarseGridData())
265 {
266 }
267
268 void
269 print(dealii::ConditionalOStream const & pcout) const
270 {
271 print_parameter(pcout, "Multigrid type", type);
272
273 if(involves_p_transfer())
274 {
275 print_parameter(pcout, "p-sequence", p_sequence);
276 }
277
278 smoother_data.print(pcout);
279
280 coarse_problem.print(pcout);
281 }
282
283 bool
284 involves_h_transfer() const;
285
286 bool
287 involves_c_transfer() const;
288
289 bool
290 involves_p_transfer() const;
291
292 // Multigrid type: p-MG vs. h-MG
293 MultigridType type;
294
295 // Sequence of polynomial degrees during p-multigrid
296 PSequenceType p_sequence;
297
298 // Smoother data
299 SmootherData smoother_data;
300
301 // Coarse grid problem
302 CoarseGridData coarse_problem;
303};
304
305} // namespace ExaDG
306
307
308#endif /* INCLUDE_SOLVERS_AND_PRECONDITIONERS_MULTIGRIDINPUTPARAMETERS_H_ */
Definition driver.cpp:33
Definition multigrid_parameters.h:98
Definition multigrid_parameters.h:220
Definition multigrid_parameters.h:259
Definition multigrid_parameters.h:170
Definition solver_data.h:34