22#ifndef EXADG_OPERATORS_SOLUTION_PROJECTION_BETWEEN_TRIANGULATIONS_H_
23#define EXADG_OPERATORS_SOLUTION_PROJECTION_BETWEEN_TRIANGULATIONS_H_
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>
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>
48namespace GridToGridProjection
52struct GridToGridProjectionData
54 GridToGridProjectionData()
57 preconditioner(PreconditionerMass::PointJacobi),
59 additional_quadrature_points(1)
64 print(dealii::ConditionalOStream
const & pcout)
const
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);
72 solver_data.print(pcout);
73 print_parameter(pcout,
"Preconditioner", preconditioner);
74 amg_data.print(pcout);
77 typename dealii::Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData rpe_data;
79 PreconditionerMass preconditioner;
85 unsigned int additional_quadrature_points;
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)
99 CellIntegrator<dim, n_components, Number> fe_eval(matrix_free, dof_index, quad_index);
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() *
107 for(
unsigned int cell_batch_idx = 0; cell_batch_idx < matrix_free.n_cell_batches();
110 fe_eval.reinit(cell_batch_idx);
111 for(
const unsigned int q : fe_eval.quadrature_point_indices())
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)
117 dealii::Point<dim> p;
118 for(
unsigned int d = 0; d < dim; ++d)
120 p[d] = cell_batch_points[d][i];
122 integration_points.push_back(p);
127 return integration_points;
134template<
int dim,
int n_components,
typename Number,
typename VectorType>
136assemble_projection_rhs(
137 dealii::MatrixFree<dim, Number, dealii::VectorizedArray<Number>>
const & matrix_free,
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)
144 VectorType system_rhs;
145 matrix_free.initialize_dof_vector(system_rhs, dof_index);
147 CellIntegrator<dim, n_components, Number> fe_eval(matrix_free, dof_index, quad_index);
149 unsigned int idx_q_point = 0;
151 for(
unsigned int cell_batch_idx = 0; cell_batch_idx < matrix_free.n_cell_batches();
154 fe_eval.reinit(cell_batch_idx);
155 for(
unsigned int const q : fe_eval.quadrature_point_indices())
157 dealii::Tensor<1, n_components, dealii::VectorizedArray<Number>> tmp;
159 for(
unsigned int i = 0; i < dealii::VectorizedArray<Number>::size(); ++i)
161 typename dealii::FEPointEvaluation<n_components, dim, dim, Number>::value_type
const
162 values = values_source_in_q_points_target[idx_q_point];
168 if constexpr(n_components == 1)
174 for(
unsigned int c = 0; c < n_components; ++c)
176 tmp[c][i] = values[c];
181 fe_eval.submit_value(tmp, q);
183 fe_eval.integrate(dealii::EvaluationFlags::values);
184 fe_eval.distribute_local_to_global(system_rhs);
186 system_rhs.compress(dealii::VectorOperation::add);
195template<
int dim,
typename Number,
int n_components,
typename VectorType>
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,
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 inverse_mass_operator_data.parameters.implementation_type =
215 InverseMassOperatorData<Number>::get_optimal_inverse_mass_type(
216 target_matrix_free.get_dof_handler(dof_index).get_fe());
218 InverseMassOperator<dim, n_components, Number> inverse_mass_operator;
219 inverse_mass_operator.initialize(target_matrix_free, inverse_mass_operator_data);
221 dealii::Utilities::MPI::RemotePointEvaluation<dim> rpe_source(data.rpe_data);
225 std::vector<dealii::Point<dim>> integration_points_target =
226 collect_integration_points<dim, n_components, Number>(target_matrix_free,
230 rpe_source.reinit(integration_points_target,
231 source_dof_handler.get_triangulation(),
234 if(not rpe_source.all_points_found())
236 write_points_in_dummy_triangulation(
237 integration_points_target,
"./",
"points_all", 0, source_dof_handler.get_mpi_communicator());
239 std::vector<dealii::Point<dim>> points_not_found;
240 points_not_found.reserve(integration_points_target.size());
241 for(
unsigned int i = 0; i < integration_points_target.size(); ++i)
243 if(not rpe_source.point_found(i))
245 points_not_found.push_back(integration_points_target[i]);
249 write_points_in_dummy_triangulation(
250 points_not_found,
"./",
"points_not_found", 0, source_dof_handler.get_mpi_communicator());
253 rpe_source.all_points_found(),
255 "Could not interpolate source grid vector in target grid. "
256 "Points exported to `./points_all_points.pvtu` and `./points_not_found_points.pvtu`"));
260 for(
unsigned int i = 0; i < target_vectors.size(); ++i)
263 VectorType
const & source_vector = *source_vectors.at(i);
264 source_vector.update_ghost_values();
267 typename dealii::FEPointEvaluation<n_components, dim, dim, Number>::value_type>
const
268 values_source_in_q_points_target = dealii::VectorTools::point_values<n_components>(
269 rpe_source, source_dof_handler, source_vector, dealii::VectorTools::EvaluationFlags::avg);
272 VectorType system_rhs = assemble_projection_rhs<dim, n_components, Number, VectorType>(
273 target_matrix_free, values_source_in_q_points_target, dof_index, quad_index);
277 sol.reinit(system_rhs,
false );
279 inverse_mass_operator.apply(sol, system_rhs);
282 *target_vectors[i] = sol;
291template<
int dim,
typename Number,
typename VectorType>
293do_grid_to_grid_projection(
294 std::shared_ptr<dealii::Mapping<dim>
const>
const & source_mapping,
295 std::vector<dealii::DoFHandler<dim>
const *>
const & source_dof_handlers,
296 std::vector<std::vector<VectorType *>>
const & source_vectors_per_dof_handler,
297 std::vector<dealii::DoFHandler<dim>
const *>
const & target_dof_handlers,
298 dealii::MatrixFree<dim, Number, dealii::VectorizedArray<Number>>
const & target_matrix_free,
299 std::vector<std::vector<VectorType *>> & target_vectors_per_dof_handler,
303 AssertThrow(source_vectors_per_dof_handler.size() == source_dof_handlers.size(),
304 dealii::ExcMessage(
"First dimension of source vector of vectors "
305 "has to match source DoFHandler count."));
306 AssertThrow(target_vectors_per_dof_handler.size() == target_dof_handlers.size(),
307 dealii::ExcMessage(
"First dimension of target vector of vectors "
308 "has to match target DoFHandler count."));
309 AssertThrow(source_dof_handlers.size() == target_dof_handlers.size(),
310 dealii::ExcMessage(
"Target and source DoFHandler counts have to match"));
311 AssertThrow(source_vectors_per_dof_handler.size() > 0,
312 dealii::ExcMessage(
"Vector of source vectors empty."));
313 for(
unsigned int i = 0; i < source_vectors_per_dof_handler.size(); ++i)
315 AssertThrow(source_vectors_per_dof_handler[i].size() ==
316 target_vectors_per_dof_handler.at(i).size(),
317 dealii::ExcMessage(
"Vectors of source and target vectors need to have same size."));
321 for(
unsigned int i = 0; i < target_dof_handlers.size(); ++i)
323 unsigned int const n_components = target_dof_handlers[i]->get_fe().n_components();
324 if(n_components == 1)
326 project_vectors<dim, Number, 1 , VectorType>(
328 *source_dof_handlers.at(i),
329 source_vectors_per_dof_handler.at(i),
331 target_vectors_per_dof_handler.at(i),
336 else if(n_components == dim)
338 project_vectors<dim, Number, dim , VectorType>(
340 *source_dof_handlers.at(i),
341 source_vectors_per_dof_handler.at(i),
343 target_vectors_per_dof_handler.at(i),
348 else if(n_components == dim + 2)
350 project_vectors<dim, Number, dim + 2 , VectorType>(
352 *source_dof_handlers.at(i),
353 source_vectors_per_dof_handler.at(i),
355 target_vectors_per_dof_handler.at(i),
362 AssertThrow(n_components == 1 or n_components == dim,
363 dealii::ExcMessage(
"The requested number of components is not"
364 "supported in `grid_to_grid_projection()`."));
372template<
int dim,
typename Number,
typename VectorType>
374grid_to_grid_projection(
375 std::shared_ptr<dealii::Mapping<dim>
const>
const & source_mapping,
376 std::vector<dealii::DoFHandler<dim>
const *>
const & source_dof_handlers,
377 std::vector<std::vector<VectorType *>>
const & source_vectors_per_dof_handler,
378 std::shared_ptr<dealii::Mapping<dim>
const>
const & target_mapping,
379 std::vector<dealii::DoFHandler<dim>
const *>
const & target_dof_handlers,
380 std::vector<std::vector<VectorType *>> & target_vectors_per_dof_handler,
384 MatrixFreeData<dim, Number> target_matrix_free_data;
386 MappingFlags mapping_flags;
387 mapping_flags.cells =
388 dealii::update_quadrature_points | dealii::update_values | dealii::update_JxW_values;
389 target_matrix_free_data.append_mapping_flags(mapping_flags);
391 dealii::AffineConstraints<Number> empty_constraints;
392 empty_constraints.clear();
393 empty_constraints.close();
394 for(
unsigned int i = 0; i < target_dof_handlers.size(); ++i)
396 target_matrix_free_data.insert_dof_handler(target_dof_handlers[i], std::to_string(i));
397 target_matrix_free_data.insert_constraint(&empty_constraints, std::to_string(i));
399 ElementType element_type =
get_element_type(target_dof_handlers[i]->get_triangulation());
402 element_type, target_dof_handlers[i]->get_fe().degree + data.additional_quadrature_points);
404 target_matrix_free_data.insert_quadrature(*quadrature, std::to_string(i));
407 dealii::MatrixFree<dim, Number, dealii::VectorizedArray<Number>> target_matrix_free;
409 target_matrix_free.reinit(*target_mapping,
410 target_matrix_free_data.get_dof_handler_vector(),
411 target_matrix_free_data.get_constraint_vector(),
412 target_matrix_free_data.get_quadrature_vector(),
413 target_matrix_free_data.data);
415 do_grid_to_grid_projection<dim, Number, VectorType>(source_mapping,
417 source_vectors_per_dof_handler,
420 target_vectors_per_dof_handler,
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:34