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/matrix_free/matrix_free.h>
32#include <deal.II/numerics/vector_tools.h>
35#include <exadg/matrix_free/integrators.h>
36#include <exadg/matrix_free/matrix_free_data.h>
37#include <exadg/operators/inverse_mass_operator.h>
38#include <exadg/operators/inverse_mass_parameters.h>
39#include <exadg/operators/mapping_flags.h>
40#include <exadg/operators/quadrature.h>
41#include <exadg/postprocessor/write_output.h>
45namespace GridToGridProjection
49struct GridToGridProjectionData
51 GridToGridProjectionData()
54 preconditioner(PreconditionerMass::PointJacobi),
55 inverse_mass_type(InverseMassType::MatrixfreeOperator),
56 additional_quadrature_points(1),
62 print(dealii::ConditionalOStream
const & pcout)
const
64 print_parameter(pcout,
"RPE tolerance", rpe_data.tolerance);
65 print_parameter(pcout,
"RPE enforce unique mapping", rpe_data.enforce_unique_mapping);
66 print_parameter(pcout,
"RPE rtree level", rpe_data.rtree_level);
69 typename dealii::Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData::AdditionalData
72 PreconditionerMass preconditioner;
73 InverseMassType inverse_mass_type;
78 unsigned int additional_quadrature_points;
88template<
int dim,
int n_components,
typename Number>
89std::vector<dealii::Point<dim>>
90collect_integration_points(
91 dealii::MatrixFree<dim, Number, dealii::VectorizedArray<Number>>
const & matrix_free,
92 unsigned int const dof_index,
93 unsigned int const quad_index)
95 CellIntegrator<dim, n_components, Number> fe_eval(matrix_free, dof_index, quad_index);
98 std::vector<dealii::Point<dim>> integration_points;
99 integration_points.reserve(
100 matrix_free.get_dof_handler(dof_index).get_triangulation().n_active_cells() *
103 for(
unsigned int cell_batch_idx = 0; cell_batch_idx < matrix_free.n_cell_batches();
106 fe_eval.reinit(cell_batch_idx);
107 for(
const unsigned int q : fe_eval.quadrature_point_indices())
109 dealii::Point<dim, dealii::VectorizedArray<Number>>
const cell_batch_points =
110 fe_eval.quadrature_point(q);
111 for(
unsigned int i = 0; i < dealii::VectorizedArray<Number>::size(); ++i)
113 dealii::Point<dim> p;
114 for(
unsigned int d = 0; d < dim; ++d)
116 p[d] = cell_batch_points[d][i];
118 integration_points.push_back(p);
123 return integration_points;
130template<
int dim,
int n_components,
typename Number,
typename VectorType>
132assemble_projection_rhs(
133 dealii::MatrixFree<dim, Number, dealii::VectorizedArray<Number>>
const & matrix_free,
135 typename dealii::FEPointEvaluation<n_components, dim, dim, Number>::value_type>
const &
136 values_source_in_q_points_target,
137 unsigned int const dof_index,
138 unsigned int const quad_index)
140 VectorType system_rhs;
141 matrix_free.initialize_dof_vector(system_rhs, dof_index);
143 CellIntegrator<dim, n_components, Number> fe_eval(matrix_free, dof_index, quad_index);
145 unsigned int idx_q_point = 0;
147 for(
unsigned int cell_batch_idx = 0; cell_batch_idx < matrix_free.n_cell_batches();
150 fe_eval.reinit(cell_batch_idx);
151 for(
unsigned int const q : fe_eval.quadrature_point_indices())
153 dealii::Tensor<1, n_components, dealii::VectorizedArray<Number>> tmp;
155 for(
unsigned int i = 0; i < dealii::VectorizedArray<Number>::size(); ++i)
157 typename dealii::FEPointEvaluation<n_components, dim, dim, Number>::value_type
const
158 values = values_source_in_q_points_target[idx_q_point];
164 if constexpr(n_components == 1)
170 for(
unsigned int c = 0; c < n_components; ++c)
172 tmp[c][i] = values[c];
177 fe_eval.submit_value(tmp, q);
179 fe_eval.integrate(dealii::EvaluationFlags::values);
180 fe_eval.distribute_local_to_global(system_rhs);
182 system_rhs.compress(dealii::VectorOperation::add);
191template<
int dim,
typename Number,
int n_components,
typename VectorType>
194 std::vector<VectorType *>
const & source_vectors,
195 dealii::DoFHandler<dim>
const & source_dof_handler,
196 std::shared_ptr<dealii::Mapping<dim>
const>
const & source_mapping,
197 std::vector<VectorType *>
const & target_vectors,
198 dealii::MatrixFree<dim, Number, dealii::VectorizedArray<Number>>
const & target_matrix_free,
199 unsigned int const dof_index,
200 unsigned int const quad_index,
204 InverseMassOperatorData<Number> inverse_mass_operator_data;
205 inverse_mass_operator_data.dof_index = dof_index;
206 inverse_mass_operator_data.quad_index = quad_index;
207 inverse_mass_operator_data.parameters.preconditioner = PreconditionerMass::PointJacobi;
208 inverse_mass_operator_data.parameters.solver_data = data.solver_data;
209 inverse_mass_operator_data.parameters.implementation_type = data.inverse_mass_type;
211 InverseMassOperator<dim, n_components, Number> inverse_mass_operator;
212 inverse_mass_operator.initialize(target_matrix_free, inverse_mass_operator_data);
214 dealii::Utilities::MPI::RemotePointEvaluation<dim> rpe_source(data.rpe_data);
218 std::vector<dealii::Point<dim>> integration_points_target =
219 collect_integration_points<dim, n_components, Number>(target_matrix_free,
223 rpe_source.reinit(integration_points_target,
224 source_dof_handler.get_triangulation(),
227 if(not rpe_source.all_points_found())
229 write_points_in_dummy_triangulation(
230 integration_points_target,
"./",
"all_points", 0, source_dof_handler.get_mpi_communicator());
232 std::vector<dealii::Point<dim>> points_not_found;
233 points_not_found.reserve(integration_points_target.size());
234 for(
unsigned int i = 0; i < integration_points_target.size(); ++i)
236 if(not rpe_source.point_found(i))
238 points_not_found.push_back(integration_points_target[i]);
242 write_points_in_dummy_triangulation(
243 points_not_found,
"./",
"points_not_found", 0, source_dof_handler.get_mpi_communicator());
245 AssertThrow(rpe_source.all_points_found(),
247 "Could not interpolate source grid vector in target grid. "
248 "Points exported to `./all_points.pvtu` and `./points_not_found.pvtu`"));
252 for(
unsigned int i = 0; i < target_vectors.size(); ++i)
255 VectorType
const & source_vector = *source_vectors.at(i);
256 source_vector.update_ghost_values();
259 typename dealii::FEPointEvaluation<n_components, dim, dim, Number>::value_type>
const
260 values_source_in_q_points_target = dealii::VectorTools::point_values<n_components>(
261 rpe_source, source_dof_handler, source_vector, dealii::VectorTools::EvaluationFlags::avg);
264 VectorType system_rhs = assemble_projection_rhs<dim, n_components, Number, VectorType>(
265 target_matrix_free, values_source_in_q_points_target, dof_index, quad_index);
269 sol.reinit(system_rhs,
false );
271 inverse_mass_operator.apply(sol, system_rhs);
274 *target_vectors[i] = sol;
283template<
int dim,
typename Number,
typename VectorType>
285do_grid_to_grid_projection(
286 std::vector<std::vector<VectorType *>>
const & source_vectors_per_dof_handler,
287 std::vector<dealii::DoFHandler<dim>
const *>
const & source_dof_handlers,
288 std::shared_ptr<dealii::Mapping<dim>
const>
const & source_mapping,
289 std::vector<std::vector<VectorType *>> & target_vectors_per_dof_handler,
290 std::vector<dealii::DoFHandler<dim>
const *>
const & target_dof_handlers,
291 dealii::MatrixFree<dim, Number, dealii::VectorizedArray<Number>>
const & matrix_free,
295 AssertThrow(source_vectors_per_dof_handler.size() == source_dof_handlers.size(),
296 dealii::ExcMessage(
"First dimension of source vector of vectors "
297 "has to match source DoFHandler count."));
298 AssertThrow(target_vectors_per_dof_handler.size() == target_dof_handlers.size(),
299 dealii::ExcMessage(
"First dimension of target vector of vectors "
300 "has to match target DoFHandler count."));
301 AssertThrow(source_dof_handlers.size() == target_dof_handlers.size(),
302 dealii::ExcMessage(
"Target and source DoFHandler counts have to match"));
303 AssertThrow(source_vectors_per_dof_handler.size() > 0,
304 dealii::ExcMessage(
"Vector of source vectors empty."));
305 for(
unsigned int i = 0; i < source_vectors_per_dof_handler.size(); ++i)
307 AssertThrow(source_vectors_per_dof_handler[i].size() ==
308 target_vectors_per_dof_handler.at(i).size(),
309 dealii::ExcMessage(
"Vectors of source and target vectors need to have same size."));
313 for(
unsigned int i = 0; i < target_dof_handlers.size(); ++i)
315 unsigned int const n_components = target_dof_handlers[i]->get_fe().n_components();
316 if(n_components == 1)
318 project_vectors<dim, Number, 1 , VectorType>(
319 source_vectors_per_dof_handler.at(i),
320 *source_dof_handlers.at(i),
322 target_vectors_per_dof_handler.at(i),
328 else if(n_components == dim)
330 project_vectors<dim, Number, dim , VectorType>(
331 source_vectors_per_dof_handler.at(i),
332 *source_dof_handlers.at(i),
334 target_vectors_per_dof_handler.at(i),
340 else if(n_components == dim + 2)
342 project_vectors<dim, Number, dim + 2 , VectorType>(
343 source_vectors_per_dof_handler.at(i),
344 *source_dof_handlers.at(i),
346 target_vectors_per_dof_handler.at(i),
354 AssertThrow(n_components == 1 or n_components == dim,
355 dealii::ExcMessage(
"The requested number of components is not"
356 "supported in `grid_to_grid_projection()`."));
364template<
int dim,
typename Number,
typename VectorType>
366grid_to_grid_projection(
367 std::vector<std::vector<VectorType *>>
const & source_vectors_per_dof_handler,
368 std::vector<dealii::DoFHandler<dim>
const *>
const & source_dof_handlers,
369 std::shared_ptr<dealii::Mapping<dim>
const>
const & source_mapping,
370 std::vector<std::vector<VectorType *>> & target_vectors_per_dof_handler,
371 std::vector<dealii::DoFHandler<dim>
const *>
const & target_dof_handlers,
372 std::shared_ptr<dealii::Mapping<dim>
const>
const & target_mapping,
376 MatrixFreeData<dim, Number> matrix_free_data;
378 MappingFlags mapping_flags;
379 mapping_flags.cells =
380 dealii::update_quadrature_points | dealii::update_values | dealii::update_JxW_values;
381 matrix_free_data.append_mapping_flags(mapping_flags);
383 dealii::AffineConstraints<Number> empty_constraints;
384 empty_constraints.clear();
385 empty_constraints.close();
386 for(
unsigned int i = 0; i < target_dof_handlers.size(); ++i)
388 matrix_free_data.insert_dof_handler(target_dof_handlers[i], std::to_string(i));
389 matrix_free_data.insert_constraint(&empty_constraints, std::to_string(i));
391 ElementType element_type =
get_element_type(target_dof_handlers[i]->get_triangulation());
394 element_type, target_dof_handlers[i]->get_fe().degree + data.additional_quadrature_points);
396 matrix_free_data.insert_quadrature(*quadrature, std::to_string(i));
399 dealii::MatrixFree<dim, Number, dealii::VectorizedArray<Number>> matrix_free;
400 matrix_free.reinit(*target_mapping,
401 matrix_free_data.get_dof_handler_vector(),
402 matrix_free_data.get_constraint_vector(),
403 matrix_free_data.get_quadrature_vector(),
404 matrix_free_data.data);
406 do_grid_to_grid_projection<dim, Number, VectorType>(source_vectors_per_dof_handler,
409 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 solution_projection_between_triangulations.h:50
Definition solver_data.h:34