ExaDG
Loading...
Searching...
No Matches
multigrid_preconditioner_base.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_MULTIGRID_PRECONDITIONER_ADAPTER_BASE_H_
23#define INCLUDE_SOLVERS_AND_PRECONDITIONERS_MULTIGRID_PRECONDITIONER_ADAPTER_BASE_H_
24
25// deal.II
26#include <deal.II/base/mg_level_object.h>
27#include <deal.II/distributed/tria.h>
28#include <deal.II/fe/mapping.h>
29#include <deal.II/multigrid/mg_constrained_dofs.h>
30
31// ExaDG
32#include <exadg/grid/grid.h>
33#include <exadg/matrix_free/matrix_free_data.h>
34#include <exadg/operators/multigrid_operator_base.h>
35#include <exadg/solvers_and_preconditioners/multigrid/coarse_grid_solvers.h>
36#include <exadg/solvers_and_preconditioners/multigrid/levels_hybrid_multigrid.h>
37#include <exadg/solvers_and_preconditioners/multigrid/multigrid_parameters.h>
38#include <exadg/solvers_and_preconditioners/multigrid/smoothers/smoother_base.h>
39#include <exadg/solvers_and_preconditioners/multigrid/transfer.h>
40#include <exadg/solvers_and_preconditioners/preconditioners/preconditioner_base.h>
41
42// forward declarations
43namespace ExaDG
44{
45template<typename VectorType, typename Operator, typename Smoother>
46class MultigridAlgorithm;
47
48template<int dim, typename Number>
49class MappingDoFVector;
50} // namespace ExaDG
51
52namespace dealii
53{
54template<typename VectorType>
56}
57
58namespace ExaDG
59{
60template<int dim, typename Number, typename MultigridNumber_ = float>
62{
63public:
64 typedef MultigridNumber_ MultigridNumber;
65
66protected:
67 typedef std::map<dealii::types::boundary_id, std::shared_ptr<dealii::Function<dim>>> Map_DBC;
68 typedef std::map<dealii::types::boundary_id, dealii::ComponentMask> Map_DBC_ComponentMask;
69
70 typedef std::vector<
71 dealii::GridTools::PeriodicFacePair<typename dealii::Triangulation<dim>::cell_iterator>>
72 PeriodicFacePairs;
73
74 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
75 typedef dealii::LinearAlgebra::distributed::Vector<MultigridNumber> VectorTypeMG;
76
77private:
79
80 typedef std::vector<std::pair<unsigned int, unsigned int>> Levels;
81
83
84public:
85 /*
86 * Constructor.
87 */
88 MultigridPreconditionerBase(MPI_Comm const & comm);
89
90 /*
91 * Destructor.
92 */
94 {
95 }
96
97 /*
98 * Initialization function.
99 */
100 void
101 initialize(MultigridData const & data,
102 std::shared_ptr<Grid<dim> const> grid,
103 std::shared_ptr<MultigridMappings<dim, Number>> const multigrid_mappings,
104 dealii::FiniteElement<dim> const & fe,
105 bool const operator_is_singular,
106 Map_DBC const & dirichlet_bc,
107 Map_DBC_ComponentMask const & dirichlet_bc_component_mask,
108 bool const initialize_preconditioners);
109
110 /*
111 * Update of multigrid preconditioner including operators, smoothers, etc. (e.g. for problems
112 * with time-dependent coefficients).
113 */
114 void
115 update() override;
116
117 /*
118 * This function applies the multigrid preconditioner dst = P^{-1} src.
119 */
120 void
121 vmult(VectorType & dst, VectorType const & src) const override;
122
123 /*
124 * Use multigrid as a solver.
125 */
126 unsigned int
127 solve(VectorType & dst, VectorType const & src) const;
128
129 /*
130 * This function applies the smoother on the fine level as a means to test the
131 * multigrid ingredients.
132 */
133 virtual void
134 apply_smoother_on_fine_level(VectorTypeMG & dst, VectorTypeMG const & src) const;
135
136 std::shared_ptr<TimerTree>
137 get_timings() const override;
138
139protected:
140 /*
141 * Initialization of mapping depending on multigrid transfer type. Note that the mapping needs to
142 * be re-initialized if the domain changes over time.
143 */
144 void
145 initialize_mapping();
146
147 /*
148 * This function initializes the matrix-free objects for all multigrid levels.
149 */
150 virtual void
151 initialize_matrix_free_objects();
152
153 /*
154 * This function updates the matrix-free objects for all multigrid levels, which
155 * is necessary if the domain changes over time.
156 */
157 void
158 update_matrix_free_objects();
159
164 void
166
171 void
173
174 /*
175 * Dof-handlers and constraints.
176 */
177 virtual void
178 initialize_dof_handler_and_constraints(bool is_singular,
179 unsigned int const n_components,
180 Map_DBC const & dirichlet_bc,
181 Map_DBC_ComponentMask const & dirichlet_bc_component_mask);
182
183 void
184 do_initialize_dof_handler_and_constraints(
185 bool is_singular,
186 unsigned int const n_components,
187 Map_DBC const & dirichlet_bc,
188 Map_DBC_ComponentMask const & dirichlet_bc_component_mask,
189 dealii::MGLevelObject<std::shared_ptr<dealii::DoFHandler<dim> const>> & dofhandlers,
190 dealii::MGLevelObject<std::shared_ptr<dealii::AffineConstraints<MultigridNumber>>> &
191 constraints);
192
193 /*
194 * Transfer operators.
195 */
196 virtual void
197 initialize_transfer_operators();
198
199 void
200 do_initialize_transfer_operators(
202 unsigned int const dof_index);
203
210 unsigned int
211 get_number_of_levels() const;
212
217 void
218 for_all_levels(std::function<void(unsigned int const)> const & function_on_level)
219 {
220 for(unsigned int level = 0; level < this->get_number_of_levels(); ++level)
221 function_on_level(level);
222 }
223
228 void
229 for_all_smoothing_levels(std::function<void(unsigned int const)> const & function_on_level)
230 {
231 // level l = 0 is the coarse problem where we do not have a smoother,
232 // so we skip the coarsest level
233 for(unsigned int level = 1; level < this->get_number_of_levels(); ++level)
234 function_on_level(level);
235 }
236
243 void
245 std::function<void(unsigned int const, unsigned int const)> const & levelwise_transfer)
246 {
247 for(unsigned int fine_level = this->get_number_of_levels() - 1; fine_level > 0; --fine_level)
248 levelwise_transfer(fine_level, fine_level - 1);
249 }
250
251 // Pointer to grid class.
252 std::shared_ptr<Grid<dim> const> grid;
253
254 std::shared_ptr<MultigridMappings<dim, Number>> multigrid_mappings;
255
256 dealii::MGLevelObject<std::shared_ptr<dealii::DoFHandler<dim> const>> dof_handlers;
257 dealii::MGLevelObject<std::shared_ptr<dealii::AffineConstraints<MultigridNumber>>> constraints;
258
259 dealii::MGLevelObject<std::shared_ptr<MatrixFreeData<dim, MultigridNumber>>>
260 matrix_free_data_objects;
261 dealii::MGLevelObject<std::shared_ptr<dealii::MatrixFree<dim, MultigridNumber>>>
262 matrix_free_objects;
263
264 dealii::MGLevelObject<std::shared_ptr<Operator>> operators;
265
266 std::shared_ptr<MultigridTransfer<dim, MultigridNumber, VectorTypeMG>> transfers;
267
268 std::vector<MGLevelInfo> level_info;
269
270private:
274 void
275 initialize_levels(unsigned int const degree, bool const is_dg);
276
277 /*
278 * Returns the correct mapping depending on the multigrid transfer type and the current h-level.
279 */
280 dealii::Mapping<dim> const &
281 get_mapping(unsigned int const h_level) const;
282
283 /*
284 * Data structures needed for matrix-free operator evaluation.
285 */
286 virtual void
287 fill_matrix_free_data(MatrixFreeData<dim, MultigridNumber> & matrix_free_data,
288 unsigned int const level,
289 unsigned int const h_level) = 0;
290
291 /*
292 * Initializes the multigrid operators for all multigrid levels.
293 */
294 void
295 initialize_operators();
296
297 /*
298 * This function initializes an operator for a specified level. It needs to be implemented by
299 * derived classes.
300 */
301 virtual std::shared_ptr<Operator>
302 initialize_operator(unsigned int const level);
303
304 /*
305 * Smoother.
306 */
307 void
308 initialize_smoothers(bool const initialize_preconditioner);
309
310 void
311 initialize_smoother(Operator & matrix, unsigned int level, bool const initialize_preconditioner);
312
313 /*
314 * Coarse grid solver.
315 */
316 void
317 initialize_coarse_solver(bool const operator_is_singular, bool const initialize_preconditioners);
318
319 /*
320 * Initialization of actual multigrid algorithm.
321 */
322 virtual void
323 initialize_multigrid_algorithm();
324
325 MPI_Comm const mpi_comm;
326
327 MultigridData data;
328
329 // TODO try to avoid this private member variable by extracting this information from level_info
330 // when needed.
331 std::vector<MGDoFHandlerIdentifier> p_levels;
332
333 dealii::MGLevelObject<std::shared_ptr<Smoother>> smoothers;
334
335 std::shared_ptr<CoarseGridSolverBase<Operator>> coarse_grid_solver;
336
337 std::shared_ptr<MultigridAlgorithm<VectorTypeMG, Operator, Smoother>> multigrid_algorithm;
338};
339} // namespace ExaDG
340
341#endif /* INCLUDE_SOLVERS_AND_PRECONDITIONERS_MULTIGRID_PRECONDITIONER_ADAPTER_BASE_H_ \
342 */
Definition grid.h:40
Definition grid.h:84
Definition multigrid_operator_base.h:35
Definition multigrid_preconditioner_base.h:62
void for_all_smoothing_levels(std::function< void(unsigned int const)> const &function_on_level)
Definition multigrid_preconditioner_base.h:229
void transfer_from_fine_to_coarse_levels(std::function< void(unsigned int const, unsigned int const)> const &levelwise_transfer)
Definition multigrid_preconditioner_base.h:244
void for_all_levels(std::function< void(unsigned int const)> const &function_on_level)
Definition multigrid_preconditioner_base.h:218
unsigned int get_number_of_levels() const
Definition multigrid_preconditioner_base.cpp:346
void update_smoothers()
Definition multigrid_preconditioner_base.cpp:769
void update_coarse_solver()
Definition multigrid_preconditioner_base.cpp:776
Definition transfer.h:39
Definition preconditioner_base.h:35
Definition smoother_base.h:29
Definition multigrid_preconditioner_base.h:55
Definition driver.cpp:33
Definition matrix_free_data.h:40
Definition multigrid_parameters.h:259