ExaDG
Loading...
Searching...
No Matches
dof_tools_extension.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2022 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_FLUID_STRUCTURE_INTERACTION_PRECICE_DOF_TOOLS_EXTENSION_H_
23#define EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_DOF_TOOLS_EXTENSION_H_
24
25// deal.II
26#include <deal.II/dofs/dof_handler.h>
27#include <deal.II/dofs/dof_tools.h>
28#include <deal.II/fe/fe.h>
29#include <deal.II/fe/mapping_q_generic.h>
30
31DEAL_II_NAMESPACE_OPEN
32
33namespace DoFTools
34{
39template<int dim, int spacedim>
40void
41map_boundary_dofs_to_support_points(
42 Mapping<dim, spacedim> const & mapping,
43 DoFHandler<dim, spacedim> const & dof_handler,
44 std::map<types::global_dof_index, Point<spacedim>> & support_points,
45 ComponentMask const & in_mask,
46 types::boundary_id const boundary_id)
47{
48 FiniteElement<dim, spacedim> const & fe = dof_handler.get_fe();
49 // check whether every fe in the collection has support points
50 Assert(fe.has_support_points(), typename FiniteElement<dim>::ExcFEHasNoSupportPoints());
51
52 Quadrature<dim - 1> const quad(fe.get_unit_face_support_points());
53
54 // Take care of components
55 ComponentMask const mask =
56 (in_mask.size() == 0 ? ComponentMask(fe.n_components(), true) : in_mask);
57
58 // Now loop over all cells and enquire the support points on each
59 // of these. we use dummy quadrature formulas where the quadrature
60 // points are located at the unit support points to enquire the
61 // location of the support points in real space.
62 //
63 // The weights of the quadrature rule have been set to invalid
64 // values by the used constructor.
65 FEFaceValues<dim, spacedim> fe_values(mapping, fe, quad, update_quadrature_points);
66
67 std::vector<types::global_dof_index> local_dof_indices;
68 for(auto const & cell : dof_handler.active_cell_iterators())
69 if(cell->is_locally_owned())
70 for(auto const & face : cell->face_iterators())
71 if(face->at_boundary() == true and face->boundary_id() == boundary_id)
72 // only work on locally relevant cells
73 {
74 fe_values.reinit(cell, face);
75
76 local_dof_indices.resize(fe.dofs_per_face);
77 face->get_dof_indices(local_dof_indices);
78
79 std::vector<Point<spacedim>> const & points = fe_values.get_quadrature_points();
80
81 for(unsigned int i = 0; i < fe.n_dofs_per_face(); ++i)
82 {
83 unsigned int const dof_comp = fe.face_system_to_component_index(i).first;
84
85 // insert the values into the map if it is a valid component
86 if(mask[dof_comp])
87 support_points[local_dof_indices[i]] = points[i];
88 }
89 }
90}
91} // namespace DoFTools
92
93DEAL_II_NAMESPACE_CLOSE
94
95#endif /* EXADG_FLUID_STRUCTURE_INTERACTION_PRECICE_DOF_TOOLS_EXTENSION_H_ */