ExaDG
Loading...
Searching...
No Matches
solution_projection_between_triangulations.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2025 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 EXADG_OPERATORS_SOLUTION_PROJECTION_BETWEEN_TRIANGULATIONS_H_
23#define EXADG_OPERATORS_SOLUTION_PROJECTION_BETWEEN_TRIANGULATIONS_H_
24
25// C/C++
26#include <algorithm>
27
28// deal.II
29#include <deal.II/base/conditional_ostream.h>
30#include <deal.II/base/mpi_remote_point_evaluation.h>
31#include <deal.II/base/timer.h>
32#include <deal.II/matrix_free/fe_point_evaluation.h>
33#include <deal.II/matrix_free/matrix_free.h>
34#include <deal.II/numerics/vector_tools.h>
35
36// ExaDG
37#include <exadg/grid/grid_data.h>
38#include <exadg/matrix_free/integrators.h>
39#include <exadg/matrix_free/matrix_free_data.h>
40#include <exadg/operators/inverse_mass_operator.h>
41#include <exadg/operators/inverse_mass_parameters.h>
42#include <exadg/operators/mapping_flags.h>
43#include <exadg/operators/quadrature.h>
44#include <exadg/postprocessor/write_output.h>
45
46namespace ExaDG
47{
48namespace GridToGridProjection
49{
50// Parameters controlling the grid-to-grid projection.
51template<int dim>
52struct GridToGridProjectionData
53{
54 GridToGridProjectionData()
55 : rpe_data(),
56 solver_data(SolverData(1e3, 1e-20, 1e-6, LinearSolver::CG)),
57 preconditioner(PreconditionerMass::PointJacobi),
58 amg_data(AMGData()),
59 additional_quadrature_points(1)
60 {
61 }
62
63 void
64 print(dealii::ConditionalOStream const & pcout) const
65 {
66 print_parameter(pcout, "RPE tolerance", rpe_data.tolerance);
67 print_parameter(pcout, "RPE enforce unique mapping", rpe_data.enforce_unique_mapping);
68 print_parameter(pcout, "RPE rtree level", rpe_data.rtree_level);
69
70 // These parameters play only a role if an iterative scheme is used for projection.
71 // That is for `InverseMassType != InverseMassType::MatrixfreeOperator` determined at runtime.
72 solver_data.print(pcout);
73 print_parameter(pcout, "Preconditioner", preconditioner);
74 amg_data.print(pcout);
75 }
76
77 typename dealii::Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData rpe_data;
78 SolverData solver_data;
79 PreconditionerMass preconditioner;
80 AMGData amg_data;
81
82 // Number of additional integration points used for sampling the source grid.
83 // The default `additional_quadrature_points = 1` considers `fe_degree + 1` quadrature points in
84 // 1D using the `fe_degree` of the target grid's finite element.
85 unsigned int additional_quadrature_points;
86};
87
92template<int dim, int n_components, typename Number>
93std::vector<dealii::Point<dim>>
94collect_integration_points(
95 dealii::MatrixFree<dim, Number, dealii::VectorizedArray<Number>> const & matrix_free,
96 unsigned int const dof_index,
97 unsigned int const quad_index)
98{
99 CellIntegrator<dim, n_components, Number> fe_eval(matrix_free, dof_index, quad_index);
100
101 // Conservative estimate for the number of points.
102 std::vector<dealii::Point<dim>> integration_points;
103 integration_points.reserve(
104 matrix_free.get_dof_handler(dof_index).get_triangulation().n_active_cells() *
105 fe_eval.n_q_points);
106
107 for(unsigned int cell_batch_idx = 0; cell_batch_idx < matrix_free.n_cell_batches();
108 ++cell_batch_idx)
109 {
110 fe_eval.reinit(cell_batch_idx);
111 for(const unsigned int q : fe_eval.quadrature_point_indices())
112 {
113 dealii::Point<dim, dealii::VectorizedArray<Number>> const cell_batch_points =
114 fe_eval.quadrature_point(q);
115 for(unsigned int i = 0; i < dealii::VectorizedArray<Number>::size(); ++i)
116 {
117 dealii::Point<dim> p;
118 for(unsigned int d = 0; d < dim; ++d)
119 {
120 p[d] = cell_batch_points[d][i];
121 }
122 integration_points.push_back(p);
123 }
124 }
125 }
126
127 return integration_points;
128}
129
134template<int dim, int n_components, typename Number, typename VectorType>
135VectorType
136assemble_projection_rhs(
137 dealii::MatrixFree<dim, Number, dealii::VectorizedArray<Number>> const & matrix_free,
138 std::vector<
139 typename dealii::FEPointEvaluation<n_components, dim, dim, Number>::value_type> const &
140 values_source_in_q_points_target,
141 unsigned int const dof_index,
142 unsigned int const quad_index)
143{
144 VectorType system_rhs;
145 matrix_free.initialize_dof_vector(system_rhs, dof_index);
146
147 CellIntegrator<dim, n_components, Number> fe_eval(matrix_free, dof_index, quad_index);
148
149 unsigned int idx_q_point = 0;
150
151 for(unsigned int cell_batch_idx = 0; cell_batch_idx < matrix_free.n_cell_batches();
152 ++cell_batch_idx)
153 {
154 fe_eval.reinit(cell_batch_idx);
155 for(unsigned int const q : fe_eval.quadrature_point_indices())
156 {
157 dealii::Tensor<1, n_components, dealii::VectorizedArray<Number>> tmp;
158
159 for(unsigned int i = 0; i < dealii::VectorizedArray<Number>::size(); ++i)
160 {
161 typename dealii::FEPointEvaluation<n_components, dim, dim, Number>::value_type const
162 values = values_source_in_q_points_target[idx_q_point];
163
164 // Increment index into `values_source_in_q_points_target`, meaning that the sequence of
165 // function values in integration points need to match the particular sequence here.
166 ++idx_q_point;
167
168 if constexpr(n_components == 1)
169 {
170 tmp[0][i] = values;
171 }
172 else
173 {
174 for(unsigned int c = 0; c < n_components; ++c)
175 {
176 tmp[c][i] = values[c];
177 }
178 }
179 }
180
181 fe_eval.submit_value(tmp, q);
182 }
183 fe_eval.integrate(dealii::EvaluationFlags::values);
184 fe_eval.distribute_local_to_global(system_rhs);
185 }
186 system_rhs.compress(dealii::VectorOperation::add);
187
188 return system_rhs;
189}
190
195template<int dim, typename Number, int n_components, typename VectorType>
196void
197project_vectors(
198 std::shared_ptr<dealii::Mapping<dim> const> const & source_mapping,
199 dealii::DoFHandler<dim> const & source_dof_handler,
200 std::vector<VectorType *> const & source_vectors,
201 dealii::MatrixFree<dim, Number, dealii::VectorizedArray<Number>> const & target_matrix_free,
202 std::vector<VectorType *> const & target_vectors,
203 unsigned int const dof_index,
204 unsigned int const quad_index,
206{
207 // Setup inverse mass operator.
208 InverseMassOperatorData<Number> inverse_mass_operator_data;
209 inverse_mass_operator_data.dof_index = dof_index;
210 inverse_mass_operator_data.quad_index = quad_index;
211 inverse_mass_operator_data.parameters.preconditioner = data.preconditioner;
212 inverse_mass_operator_data.parameters.solver_data = data.solver_data;
213 inverse_mass_operator_data.parameters.amg_data = data.amg_data;
214
215 const bool cartesian_or_affine_mapping =
216 std::all_of(target_matrix_free.get_mapping_info().cell_type.begin(),
217 target_matrix_free.get_mapping_info().cell_type.end(),
218 [](auto g) {
219 return g <= dealii::internal::MatrixFreeFunctions::GeometryType::affine;
220 });
221 inverse_mass_operator_data.parameters.implementation_type =
222 InverseMassOperatorData<Number>::get_optimal_inverse_mass_type(
223 target_matrix_free.get_dof_handler(dof_index).get_fe(), cartesian_or_affine_mapping);
224
225 InverseMassOperator<dim, n_components, Number> inverse_mass_operator;
226 inverse_mass_operator.initialize(target_matrix_free, inverse_mass_operator_data);
227
228 dealii::Utilities::MPI::RemotePointEvaluation<dim> rpe_source(data.rpe_data);
229
230 // The sequence of integration points follows from the sequence of points as encountered during
231 // cell batch loop.
232 std::vector<dealii::Point<dim>> integration_points_target =
233 collect_integration_points<dim, n_components, Number>(target_matrix_free,
234 dof_index,
235 quad_index);
236
237 rpe_source.reinit(integration_points_target,
238 source_dof_handler.get_triangulation(),
239 *source_mapping);
240
241 if(not rpe_source.all_points_found())
242 {
243 write_points_in_dummy_triangulation(
244 integration_points_target, "./", "points_all", 0, source_dof_handler.get_mpi_communicator());
245
246 std::vector<dealii::Point<dim>> points_not_found;
247 points_not_found.reserve(integration_points_target.size());
248 for(unsigned int i = 0; i < integration_points_target.size(); ++i)
249 {
250 if(not rpe_source.point_found(i))
251 {
252 points_not_found.push_back(integration_points_target[i]);
253 }
254 }
255
256 write_points_in_dummy_triangulation(
257 points_not_found, "./", "points_not_found", 0, source_dof_handler.get_mpi_communicator());
258
259 AssertThrow(
260 rpe_source.all_points_found(),
261 dealii::ExcMessage(
262 "Could not interpolate source grid vector in target grid. "
263 "Points exported to `./points_all_points.pvtu` and `./points_not_found_points.pvtu`"));
264 }
265
266 // Loop over vectors and project.
267 for(unsigned int i = 0; i < target_vectors.size(); ++i)
268 {
269 // Evaluate the source vector at the target integration points.
270 VectorType const & source_vector = *source_vectors.at(i);
271 source_vector.update_ghost_values();
272
273 std::vector<
274 typename dealii::FEPointEvaluation<n_components, dim, dim, Number>::value_type> const
275 values_source_in_q_points_target = dealii::VectorTools::point_values<n_components>(
276 rpe_source, source_dof_handler, source_vector, dealii::VectorTools::EvaluationFlags::avg);
277
278 // Assemble right-hand side vector for the projection.
279 VectorType system_rhs = assemble_projection_rhs<dim, n_components, Number, VectorType>(
280 target_matrix_free, values_source_in_q_points_target, dof_index, quad_index);
281
282 // Solve linear system starting from zero initial guess.
283 VectorType sol;
284 sol.reinit(system_rhs, false /* omit_zeroing_entries */);
285
286 inverse_mass_operator.apply(sol, system_rhs);
287
288 // Copy solution to target vector.
289 *target_vectors[i] = sol;
290 }
291}
292
298template<int dim, typename Number, typename VectorType>
299void
300do_grid_to_grid_projection(
301 std::shared_ptr<dealii::Mapping<dim> const> const & source_mapping,
302 std::vector<dealii::DoFHandler<dim> const *> const & source_dof_handlers,
303 std::vector<std::vector<VectorType *>> const & source_vectors_per_dof_handler,
304 std::vector<dealii::DoFHandler<dim> const *> const & target_dof_handlers,
305 dealii::MatrixFree<dim, Number, dealii::VectorizedArray<Number>> const & target_matrix_free,
306 std::vector<std::vector<VectorType *>> & target_vectors_per_dof_handler,
308{
309 // Check input dimensions.
310 AssertThrow(source_vectors_per_dof_handler.size() == source_dof_handlers.size(),
311 dealii::ExcMessage("First dimension of source vector of vectors "
312 "has to match source DoFHandler count."));
313 AssertThrow(target_vectors_per_dof_handler.size() == target_dof_handlers.size(),
314 dealii::ExcMessage("First dimension of target vector of vectors "
315 "has to match target DoFHandler count."));
316 AssertThrow(source_dof_handlers.size() == target_dof_handlers.size(),
317 dealii::ExcMessage("Target and source DoFHandler counts have to match"));
318 AssertThrow(source_vectors_per_dof_handler.size() > 0,
319 dealii::ExcMessage("Vector of source vectors empty."));
320 for(unsigned int i = 0; i < source_vectors_per_dof_handler.size(); ++i)
321 {
322 AssertThrow(source_vectors_per_dof_handler[i].size() ==
323 target_vectors_per_dof_handler.at(i).size(),
324 dealii::ExcMessage("Vectors of source and target vectors need to have same size."));
325 }
326
327 // Project vectors per `dealii::DoFHandler`.
328 for(unsigned int i = 0; i < target_dof_handlers.size(); ++i)
329 {
330 unsigned int const n_components = target_dof_handlers[i]->get_fe().n_components();
331 if(n_components == 1)
332 {
333 project_vectors<dim, Number, 1 /* n_components */, VectorType>(
334 source_mapping,
335 *source_dof_handlers.at(i),
336 source_vectors_per_dof_handler.at(i),
337 target_matrix_free,
338 target_vectors_per_dof_handler.at(i),
339 i /* dof_index */,
340 i /* quad_index */,
341 data);
342 }
343 else if(n_components == dim)
344 {
345 project_vectors<dim, Number, dim /* n_components */, VectorType>(
346 source_mapping,
347 *source_dof_handlers.at(i),
348 source_vectors_per_dof_handler.at(i),
349 target_matrix_free,
350 target_vectors_per_dof_handler.at(i),
351 i /* dof_index */,
352 i /* quad_index */,
353 data);
354 }
355 else if(n_components == dim + 2)
356 {
357 project_vectors<dim, Number, dim + 2 /* n_components */, VectorType>(
358 source_mapping,
359 *source_dof_handlers.at(i),
360 source_vectors_per_dof_handler.at(i),
361 target_matrix_free,
362 target_vectors_per_dof_handler.at(i),
363 i /* dof_index */,
364 i /* quad_index */,
365 data);
366 }
367 else
368 {
369 AssertThrow(n_components == 1 or n_components == dim,
370 dealii::ExcMessage("The requested number of components is not"
371 "supported in `grid_to_grid_projection()`."));
372 }
373 }
374}
375
379template<int dim, typename Number, typename VectorType>
380void
381grid_to_grid_projection(
382 std::shared_ptr<dealii::Mapping<dim> const> const & source_mapping,
383 std::vector<dealii::DoFHandler<dim> const *> const & source_dof_handlers,
384 std::vector<std::vector<VectorType *>> const & source_vectors_per_dof_handler,
385 std::shared_ptr<dealii::Mapping<dim> const> const & target_mapping,
386 std::vector<dealii::DoFHandler<dim> const *> const & target_dof_handlers,
387 std::vector<std::vector<VectorType *>> & target_vectors_per_dof_handler,
389{
390 // Setup a single `dealii::MatrixFree` object with multiple `dealii::DoFHandler`s.
391 MatrixFreeData<dim, Number> target_matrix_free_data;
392
393 MappingFlags mapping_flags;
394 mapping_flags.cells =
395 dealii::update_quadrature_points | dealii::update_values | dealii::update_JxW_values;
396 target_matrix_free_data.append_mapping_flags(mapping_flags);
397
398 dealii::AffineConstraints<Number> empty_constraints;
399 empty_constraints.clear();
400 empty_constraints.close();
401 for(unsigned int i = 0; i < target_dof_handlers.size(); ++i)
402 {
403 target_matrix_free_data.insert_dof_handler(target_dof_handlers[i], std::to_string(i));
404 target_matrix_free_data.insert_constraint(&empty_constraints, std::to_string(i));
405
406 ElementType element_type = get_element_type(target_dof_handlers[i]->get_triangulation());
407
408 std::shared_ptr<dealii::Quadrature<dim>> quadrature = create_quadrature<dim>(
409 element_type, target_dof_handlers[i]->get_fe().degree + data.additional_quadrature_points);
410
411 target_matrix_free_data.insert_quadrature(*quadrature, std::to_string(i));
412 }
413
414 dealii::MatrixFree<dim, Number, dealii::VectorizedArray<Number>> target_matrix_free;
415
416 target_matrix_free.reinit(*target_mapping,
417 target_matrix_free_data.get_dof_handler_vector(),
418 target_matrix_free_data.get_constraint_vector(),
419 target_matrix_free_data.get_quadrature_vector(),
420 target_matrix_free_data.data);
421
422 do_grid_to_grid_projection<dim, Number, VectorType>(source_mapping,
423 source_dof_handlers,
424 source_vectors_per_dof_handler,
425 target_dof_handlers,
426 target_matrix_free,
427 target_vectors_per_dof_handler,
428 data);
429}
430
431} // namespace GridToGridProjection
432} // namespace ExaDG
433
434#endif /* EXADG_OPERATORS_SOLUTION_PROJECTION_BETWEEN_TRIANGULATIONS_H_ */
Definition driver.cpp:33
std::shared_ptr< dealii::Quadrature< dim > > create_quadrature(ElementType const &element_type, unsigned int const n_q_points_1d)
Definition quadrature.h:40
ElementType get_element_type(dealii::Triangulation< dim > const &tria)
Definition grid_data.h:70
Definition multigrid_parameters.h:98
Definition solution_projection_between_triangulations.h:53
Definition solver_data.h:92