ExaDG
Loading...
Searching...
No Matches
mapping_dof_vector.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_FUNCTIONALITIES_MESH_H_
23#define INCLUDE_FUNCTIONALITIES_MESH_H_
24
25// deal.II
26#include <deal.II/base/conditional_ostream.h>
27#include <deal.II/base/mg_level_object.h>
28#include <deal.II/base/multithread_info.h>
29#include <deal.II/dofs/dof_tools.h>
30#include <deal.II/fe/fe_nothing.h>
31#include <deal.II/fe/fe_q.h>
32#include <deal.II/fe/fe_system.h>
33#include <deal.II/fe/fe_tools.h>
34#include <deal.II/fe/fe_values.h>
35#include <deal.II/fe/mapping_q_cache.h>
36#include <deal.II/lac/la_parallel_block_vector.h>
37#include <deal.II/multigrid/mg_transfer_matrix_free.h>
38
39// ExaDG
40#include <exadg/grid/grid_data.h>
41#include <exadg/solvers_and_preconditioners/multigrid/transfer.h>
42
43namespace ExaDG
44{
52template<int dim, typename Number>
54{
55public:
56 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
57
61 MappingDoFVector(unsigned int const mapping_degree)
62 {
63 mapping_q_cache = std::make_shared<dealii::MappingQCache<dim>>(mapping_degree);
64
65 hierarchic_to_lexicographic_numbering =
66 dealii::FETools::hierarchic_to_lexicographic_numbering<dim>(mapping_degree);
67 lexicographic_to_hierarchic_numbering =
68 dealii::Utilities::invert_permutation(hierarchic_to_lexicographic_numbering);
69 }
70
75 {
76 }
77
81 std::shared_ptr<dealii::Mapping<dim> const>
83 {
84 AssertThrow(mapping_q_cache.get(),
85 dealii::ExcMessage("Mapping object mapping_q_cache is not initialized."));
86
87 return mapping_q_cache;
88 }
89
90 std::shared_ptr<dealii::MappingQCache<dim>>
91 get_mapping_q_cache() const
92 {
93 AssertThrow(mapping_q_cache.get(),
94 dealii::ExcMessage("Mapping object mapping_q_cache is not initialized."));
95
96 return mapping_q_cache;
97 }
98
104 void
105 fill_grid_coordinates_vector(VectorType & grid_coordinates,
106 dealii::DoFHandler<dim> const & dof_handler) const
107 {
108 AssertThrow(mapping_q_cache.get(),
109 dealii::ExcMessage("Mapping object mapping_q_cache is not initialized."));
110
111 // use the deformed state described by the dealii::MappingQCache object
112 fill_grid_coordinates_vector(*mapping_q_cache, grid_coordinates, dof_handler);
113 }
114
119 void
120 fill_grid_coordinates_vector(dealii::Mapping<dim> const & mapping,
121 VectorType & grid_coordinates,
122 dealii::DoFHandler<dim> const & dof_handler) const
123 {
124 if(grid_coordinates.size() != dof_handler.n_dofs())
125 {
126 dealii::IndexSet relevant_dofs_grid;
127 dealii::DoFTools::extract_locally_relevant_dofs(dof_handler, relevant_dofs_grid);
128 grid_coordinates.reinit(dof_handler.locally_owned_dofs(),
129 relevant_dofs_grid,
130 dof_handler.get_communicator());
131 }
132 else
133 {
134 grid_coordinates = 0;
135 }
136
137 AssertThrow(get_element_type(dof_handler.get_triangulation()) == ElementType::Hypercube,
138 dealii::ExcMessage("Only implemented for hypercube elements."));
139
140 dealii::FiniteElement<dim> const & fe = dof_handler.get_fe();
141
142 std::vector<std::array<unsigned int, dim>> component_to_system_index(
143 fe.base_element(0).dofs_per_cell);
144
145 if(fe.dofs_per_vertex > 0) // dealii::FE_Q
146 {
147 for(unsigned int i = 0; i < fe.dofs_per_cell; ++i)
148 {
149 component_to_system_index
150 [hierarchic_to_lexicographic_numbering[fe.system_to_component_index(i).second]]
151 [fe.system_to_component_index(i).first] = i;
152 }
153 }
154 else // dealii::FE_DGQ
155 {
156 for(unsigned int i = 0; i < fe.dofs_per_cell; ++i)
157 {
158 component_to_system_index[fe.system_to_component_index(i).second]
159 [fe.system_to_component_index(i).first] = i;
160 }
161 }
162
163 // Set up dealii::FEValues with FE_Nothing and the Gauss-Lobatto quadrature to
164 // reduce setup cost, as we only use the geometry information (this means
165 // we need to call fe_values.reinit(cell) with Triangulation::cell_iterator
166 // rather than dealii::DoFHandler::cell_iterator).
167 dealii::FE_Nothing<dim> fe_nothing;
168 dealii::FEValues<dim> fe_values(mapping,
169 fe_nothing,
170 dealii::QGaussLobatto<dim>(fe.degree + 1),
171 dealii::update_quadrature_points);
172
173 std::vector<dealii::types::global_dof_index> dof_indices(fe.dofs_per_cell);
174
175 for(auto const & cell : dof_handler.active_cell_iterators())
176 {
177 if(not cell->is_artificial())
178 {
179 fe_values.reinit(typename dealii::Triangulation<dim>::cell_iterator(cell));
180 cell->get_dof_indices(dof_indices);
181 for(unsigned int i = 0; i < fe_values.n_quadrature_points; ++i)
182 {
183 dealii::Point<dim> const point = fe_values.quadrature_point(i);
184 for(unsigned int d = 0; d < dim; ++d)
185 {
186 if(grid_coordinates.get_partitioner()->in_local_range(
187 dof_indices[component_to_system_index[i][d]]))
188 {
189 grid_coordinates(dof_indices[component_to_system_index[i][d]]) = point[d];
190 }
191 }
192 }
193 }
194 }
195
196 grid_coordinates.update_ghost_values();
197 }
198
211 void
212 initialize_mapping_from_dof_vector(std::shared_ptr<dealii::Mapping<dim> const> mapping,
213 VectorType const & displacement_vector,
214 dealii::DoFHandler<dim> const & dof_handler)
215 {
216 AssertThrow(dealii::MultithreadInfo::n_threads() == 1, dealii::ExcNotImplemented());
217
218 AssertThrow(mapping_q_cache.get(),
219 dealii::ExcMessage("Mapping object mapping_q_cache is not initialized."));
220
221 VectorType displacement_vector_ghosted;
222 if(dof_handler.n_dofs() > 0 and displacement_vector.size() == dof_handler.n_dofs())
223 {
224 dealii::IndexSet locally_relevant_dofs;
225 dealii::DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
226 displacement_vector_ghosted.reinit(dof_handler.locally_owned_dofs(),
227 locally_relevant_dofs,
228 dof_handler.get_communicator());
229 displacement_vector_ghosted.copy_locally_owned_data_from(displacement_vector);
230 displacement_vector_ghosted.update_ghost_values();
231 }
232
233 AssertThrow(get_element_type(dof_handler.get_triangulation()) == ElementType::Hypercube,
234 dealii::ExcMessage("Only implemented for hypercube elements."));
235
236 std::shared_ptr<dealii::FEValues<dim>> fe_values;
237
238 // Set up dealii::FEValues with FE_Nothing and the Gauss-Lobatto quadrature to
239 // reduce setup cost, as we only use the geometry information (this means
240 // we need to call fe_values.reinit(cell) with Triangulation::cell_iterator
241 // rather than dealii::DoFHandler::cell_iterator).
242 dealii::FE_Nothing<dim> fe_nothing;
243
244 if(mapping.get() != 0)
245 {
246 fe_values = std::make_shared<dealii::FEValues<dim>>(*mapping,
247 fe_nothing,
248 dealii::QGaussLobatto<dim>(
249 mapping_q_cache->get_degree() + 1),
250 dealii::update_quadrature_points);
251 }
252
253 // take the grid coordinates described by mapping and add deformation described by displacement
254 // vector
255 mapping_q_cache->initialize(
256 dof_handler.get_triangulation(),
257 [&](const typename dealii::Triangulation<dim>::cell_iterator & cell_tria)
258 -> std::vector<dealii::Point<dim>> {
259 unsigned int const scalar_dofs_per_cell =
260 dealii::Utilities::pow(mapping_q_cache->get_degree() + 1, dim);
261
262 std::vector<dealii::Point<dim>> grid_coordinates(scalar_dofs_per_cell);
263
264 if(mapping.get() != 0)
265 {
266 fe_values->reinit(cell_tria);
267 // extract displacement and add to original position
268 for(unsigned int i = 0; i < scalar_dofs_per_cell; ++i)
269 {
270 grid_coordinates[i] =
271 fe_values->quadrature_point(this->hierarchic_to_lexicographic_numbering[i]);
272 }
273 }
274
275 // if this function is called with an empty dof-vector, this indicates that the
276 // displacements are zero and the points do not have to be moved
277 if(dof_handler.n_dofs() > 0 and displacement_vector.size() > 0 and
278 cell_tria->is_active() and not(cell_tria->is_artificial()))
279 {
280 typename dealii::DoFHandler<dim>::cell_iterator cell(&cell_tria->get_triangulation(),
281 cell_tria->level(),
282 cell_tria->index(),
283 &dof_handler);
284
285 dealii::FiniteElement<dim> const & fe = dof_handler.get_fe();
286 AssertThrow(fe.element_multiplicity(0) == dim,
287 dealii::ExcMessage("Expected finite element with dim components."));
288
289 std::vector<dealii::types::global_dof_index> dof_indices(fe.dofs_per_cell);
290 cell->get_dof_indices(dof_indices);
291
292 for(unsigned int i = 0; i < dof_indices.size(); ++i)
293 {
294 std::pair<unsigned int, unsigned int> const id = fe.system_to_component_index(i);
295
296 if(fe.dofs_per_vertex > 0) // dealii::FE_Q
297 {
298 grid_coordinates[id.second][id.first] += displacement_vector_ghosted(dof_indices[i]);
299 }
300 else // dealii::FE_DGQ
301 {
302 grid_coordinates[this->lexicographic_to_hierarchic_numbering[id.second]][id.first] +=
303 displacement_vector_ghosted(dof_indices[i]);
304 }
305 }
306 }
307
308 return grid_coordinates;
309 });
310 }
311
312 std::vector<unsigned int> hierarchic_to_lexicographic_numbering;
313 std::vector<unsigned int> lexicographic_to_hierarchic_numbering;
314
315protected:
316 std::shared_ptr<dealii::MappingQCache<dim>> mapping_q_cache;
317};
318
319
320namespace MappingTools
321{
335template<int dim, typename Number>
336void
337initialize_coarse_mappings_from_mapping_dof_vector(
338 std::vector<std::shared_ptr<MappingDoFVector<dim, Number>>> & coarse_mappings,
339 unsigned int const degree_coarse_mappings,
340 std::shared_ptr<MappingDoFVector<dim, Number> const> const & fine_mapping,
341 dealii::Triangulation<dim> const & triangulation)
342{
343 AssertThrow(dealii::MultithreadInfo::n_threads() == 1, dealii::ExcNotImplemented());
344
345 std::shared_ptr<dealii::MappingQCache<dim>> mapping_q_cache = fine_mapping->get_mapping_q_cache();
346 AssertThrow(mapping_q_cache.get(),
347 dealii::ExcMessage("Shared pointer mapping_q_cache is invalid."));
348
349 dealii::FESystem<dim> fe(dealii::FE_Q<dim>(degree_coarse_mappings), dim);
350 dealii::DoFHandler<dim> dof_handler(triangulation);
351 dof_handler.distribute_dofs(fe);
352 dof_handler.distribute_mg_dofs();
353
354 std::shared_ptr<MappingDoFVector<dim, Number>> mapping_dof_vector_all_levels =
355 std::make_shared<MappingDoFVector<dim, Number>>(degree_coarse_mappings);
356
357 // fill a dof vector with grid coordinates of the fine level using degree_coarse_mappings
358 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
359 VectorType grid_coordinates_fine_level;
360
361 {
362 mapping_dof_vector_all_levels->fill_grid_coordinates_vector(*fine_mapping->get_mapping(),
363 grid_coordinates_fine_level,
364 dof_handler);
365 }
366
367 // project the solution onto all coarse levels of the triangulation using degree_coarse_mappings
368 dealii::MGLevelObject<VectorType> grid_coordinates_all_levels,
369 grid_coordinates_all_levels_ghosted;
370 unsigned int const n_levels = triangulation.n_global_levels();
371 grid_coordinates_all_levels.resize(0, n_levels - 1);
372 grid_coordinates_all_levels_ghosted.resize(0, n_levels - 1);
373
374 dealii::MGTransferMatrixFree<dim, Number> transfer;
375 transfer.build(dof_handler);
376 transfer.interpolate_to_mg(dof_handler, grid_coordinates_all_levels, grid_coordinates_fine_level);
377
378 // ghosting
379 for(unsigned int level = 0; level < n_levels; level++)
380 {
381 dealii::IndexSet relevant_dofs;
382 dealii::DoFTools::extract_locally_relevant_level_dofs(dof_handler, level, relevant_dofs);
383
384 grid_coordinates_all_levels_ghosted[level].reinit(
385 dof_handler.locally_owned_mg_dofs(level),
386 relevant_dofs,
387 grid_coordinates_all_levels[level].get_mpi_communicator());
388
389 grid_coordinates_all_levels_ghosted[level].copy_locally_owned_data_from(
390 grid_coordinates_all_levels[level]);
391
392 grid_coordinates_all_levels_ghosted[level].update_ghost_values();
393 }
394
395 // Call the initialize() function of dealii::MappingQCache, which initializes the mapping for all
396 // levels according to grid_coordinates_all_levels_ghosted.
397 mapping_dof_vector_all_levels->get_mapping_q_cache()->initialize(
398 dof_handler.get_triangulation(),
399 [&](const typename dealii::Triangulation<dim>::cell_iterator & cell_tria)
400 -> std::vector<dealii::Point<dim>> {
401 unsigned int const level = cell_tria->level();
402
403 typename dealii::DoFHandler<dim>::cell_iterator cell(&cell_tria->get_triangulation(),
404 level,
405 cell_tria->index(),
406 &dof_handler);
407
408 AssertThrow(fe.element_multiplicity(0) == dim,
409 dealii::ExcMessage("Expected finite element with dim components."));
410
411 unsigned int const scalar_dofs_per_cell = dealii::Utilities::pow(fe.degree + 1, dim);
412
413 std::vector<dealii::Point<dim>> grid_coordinates(scalar_dofs_per_cell);
414
415 if(cell->level_subdomain_id() != dealii::numbers::artificial_subdomain_id)
416 {
417 std::vector<dealii::types::global_dof_index> dof_indices(fe.dofs_per_cell);
418 cell->get_mg_dof_indices(dof_indices);
419
420 for(unsigned int i = 0; i < dof_indices.size(); ++i)
421 {
422 std::pair<unsigned int, unsigned int> const id = fe.system_to_component_index(i);
423
424 if(fe.dofs_per_vertex > 0) // dealii::FE_Q
425 {
426 grid_coordinates[id.second][id.first] =
427 grid_coordinates_all_levels_ghosted[level](dof_indices[i]);
428 }
429 else // dealii::FE_DGQ
430 {
431 grid_coordinates[mapping_dof_vector_all_levels
432 ->lexicographic_to_hierarchic_numbering[id.second]][id.first] =
433 grid_coordinates_all_levels_ghosted[level](dof_indices[i]);
434 }
435 }
436 }
437
438 return grid_coordinates;
439 });
440
441
442 AssertThrow(
443 coarse_mappings.size() == n_levels - 1,
444 dealii::ExcMessage(
445 "coarse_mappings does not have correct size relative to the number of levels of the triangulation."));
446
447 // Finally, let all coarse grid mappings point to the same MappingDoFVector object. Using the
448 // same Mapping object for all multigrid h-levels is some form of legacy code. The class
449 // dealii::MatrixFree can internally extract the coarse-level mapping information (provided
450 // through the fine-level Mapping object).
451 for(unsigned int h_level = 0; h_level < coarse_mappings.size(); ++h_level)
452 {
453 coarse_mappings[h_level] = mapping_dof_vector_all_levels;
454 }
455}
456
472template<int dim, typename Number>
473void
474initialize_coarse_mappings_from_mapping_dof_vector(
475 std::vector<std::shared_ptr<MappingDoFVector<dim, Number>>> & coarse_mappings,
476 unsigned int const degree_coarse_mappings,
477 std::shared_ptr<MappingDoFVector<dim, Number> const> const & fine_mapping,
478 std::shared_ptr<dealii::Triangulation<dim> const> const & fine_triangulation,
479 std::vector<std::shared_ptr<dealii::Triangulation<dim> const>> const & coarse_triangulations)
480{
481 std::shared_ptr<dealii::MappingQCache<dim>> mapping_q_cache = fine_mapping->get_mapping_q_cache();
482 AssertThrow(mapping_q_cache.get(),
483 dealii::ExcMessage("Shared pointer mapping_q_cache is invalid."));
484
485 // setup dof-handlers and constraints for all levels using degree_coarse_mappings
486 dealii::FESystem<dim> fe(dealii::FE_Q<dim>(degree_coarse_mappings), dim);
487 unsigned int const n_h_levels = coarse_triangulations.size() + 1;
488 std::vector<dealii::DoFHandler<dim>> dof_handlers_all_levels(n_h_levels);
489 std::vector<dealii::AffineConstraints<Number>> constraints_all_levels(n_h_levels);
490 for(unsigned int h_level = 0; h_level < n_h_levels; ++h_level)
491 {
492 if(h_level == n_h_levels - 1)
493 dof_handlers_all_levels[h_level].reinit(*fine_triangulation);
494 else
495 dof_handlers_all_levels[h_level].reinit(*coarse_triangulations[h_level]);
496 dof_handlers_all_levels[h_level].distribute_dofs(fe);
497 // constraints are irrelevant for interpolation
498 constraints_all_levels[h_level].close();
499 }
500
501 // fill a dof vector with grid coordinates of the fine level using degree_coarse_mappings
502 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
503 VectorType grid_coordinates_fine_level;
504
505 {
506 std::shared_ptr<MappingDoFVector<dim, Number>> mapping_dof_vector_fine_level =
507 std::make_shared<MappingDoFVector<dim, Number>>(degree_coarse_mappings);
508
509 auto const & dof_handler_fine_level = dof_handlers_all_levels[n_h_levels - 1];
510 mapping_dof_vector_fine_level->fill_grid_coordinates_vector(*fine_mapping->get_mapping(),
511 grid_coordinates_fine_level,
512 dof_handler_fine_level);
513 }
514
515 // create transfer objects
516 dealii::MGLevelObject<dealii::MGTwoLevelTransfer<dim, VectorType>> transfers(0, n_h_levels - 1);
517 for(unsigned int h_level = 1; h_level < n_h_levels; ++h_level)
518 {
519 transfers[h_level].reinit_geometric_transfer(dof_handlers_all_levels[h_level],
520 dof_handlers_all_levels[h_level - 1],
521 constraints_all_levels[h_level],
522 constraints_all_levels[h_level - 1]);
523 }
524
525 // a function that initializes the dof-vector for a given level and dof_handler
526 const std::function<void(unsigned int const, VectorType &)> initialize_dof_vector =
527 [&](unsigned int const h_level, VectorType & vector) {
528 dealii::IndexSet locally_relevant_dofs;
529 dealii::DoFTools::extract_locally_relevant_dofs(dof_handlers_all_levels[h_level],
530 locally_relevant_dofs);
531 vector.reinit(dof_handlers_all_levels[h_level].locally_owned_dofs(),
532 locally_relevant_dofs,
533 dof_handlers_all_levels[h_level].get_communicator());
534 };
535
536 dealii::MGTransferGlobalCoarsening<dim, VectorType> mg_transfer_global_coarsening(
537 transfers, initialize_dof_vector);
538
539 // Transfer grid coordinates to coarser h-levels.
540 // The dealii::DoFHandler object will not be used for global coarsening.
541 dealii::DoFHandler<dim> dof_handler_dummy;
542 dealii::MGLevelObject<VectorType> grid_coordinates_all_levels(0, n_h_levels - 1);
543 mg_transfer_global_coarsening.interpolate_to_mg(dof_handler_dummy,
544 grid_coordinates_all_levels,
545 grid_coordinates_fine_level);
546
547 // initialize mapping for all coarse h-levels using the dof-vectors with grid coordinates
548 AssertThrow(coarse_mappings.size() == n_h_levels - 1,
549 dealii::ExcMessage(
550 "coarse_mappings does not have correct size relative to coarse_triangulations."));
551
552 for(unsigned int h_level = 0; h_level < coarse_mappings.size(); ++h_level)
553 {
554 coarse_mappings[h_level] =
555 std::make_shared<MappingDoFVector<dim, Number>>(degree_coarse_mappings);
556
557 // grid_coordinates_all_levels describes absolute coordinates -> use an uninitialized mapping
558 // in order to interpret the grid coordinates vector as absolute coordinates and not as
559 // displacements.
560 std::shared_ptr<dealii::Mapping<dim> const> mapping_dummy;
561 coarse_mappings[h_level]->initialize_mapping_from_dof_vector(
562 mapping_dummy, grid_coordinates_all_levels[h_level], dof_handlers_all_levels[h_level]);
563 }
564}
565
581template<int dim, typename Number>
582void
583initialize_coarse_mappings(
584 std::vector<std::shared_ptr<MappingDoFVector<dim, Number>>> & coarse_mappings,
585 unsigned int const degree_coarse_mappings,
586 std::shared_ptr<MappingDoFVector<dim, Number> const> const & fine_mapping,
587 std::shared_ptr<dealii::Triangulation<dim> const> const & fine_triangulation,
588 std::vector<std::shared_ptr<dealii::Triangulation<dim> const>> const & coarse_triangulations)
589{
590 if(fine_mapping.get())
591 {
592 if(coarse_triangulations.size() > 0)
593 {
594 MappingTools::initialize_coarse_mappings_from_mapping_dof_vector<dim, Number>(
595 coarse_mappings,
596 degree_coarse_mappings,
597 fine_mapping,
598 fine_triangulation,
599 coarse_triangulations);
600 }
601 else
602 {
603 MappingTools::initialize_coarse_mappings_from_mapping_dof_vector<dim, Number>(
604 coarse_mappings, degree_coarse_mappings, fine_mapping, *fine_triangulation);
605 }
606 }
607}
608
609
610} // namespace MappingTools
611
612} // namespace ExaDG
613
614#endif /* INCLUDE_FUNCTIONALITIES_MESH_H_ */
Definition mapping_dof_vector.h:54
MappingDoFVector(unsigned int const mapping_degree)
Definition mapping_dof_vector.h:61
void initialize_mapping_from_dof_vector(std::shared_ptr< dealii::Mapping< dim > const > mapping, VectorType const &displacement_vector, dealii::DoFHandler< dim > const &dof_handler)
Definition mapping_dof_vector.h:212
virtual ~MappingDoFVector()
Definition mapping_dof_vector.h:74
void fill_grid_coordinates_vector(dealii::Mapping< dim > const &mapping, VectorType &grid_coordinates, dealii::DoFHandler< dim > const &dof_handler) const
Definition mapping_dof_vector.h:120
void fill_grid_coordinates_vector(VectorType &grid_coordinates, dealii::DoFHandler< dim > const &dof_handler) const
Definition mapping_dof_vector.h:105
std::shared_ptr< dealii::Mapping< dim > const > get_mapping() const
Definition mapping_dof_vector.h:82
Definition driver.cpp:33
ElementType get_element_type(dealii::Triangulation< dim > const &tria)
Definition grid_data.h:70