ExaDG
Loading...
Searching...
No Matches
constraints.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 EXADG_SOLVERS_AND_PRECONDITIONERS_MULTIGRID_CONSTRAINTS_H_
23#define EXADG_SOLVERS_AND_PRECONDITIONERS_MULTIGRID_CONSTRAINTS_H_
24
25// deal.II
26#include <deal.II/multigrid/mg_constrained_dofs.h>
27
28namespace ExaDG
29{
30namespace ConstraintUtil
31{
32namespace // anonymous namespace
33{
34template<int dim, typename Number>
35void
36add_periodicity_constraints(unsigned int const level,
37 unsigned int const target_level,
38 typename dealii::DoFHandler<dim>::face_iterator const face1,
39 typename dealii::DoFHandler<dim>::face_iterator const face2,
40 dealii::AffineConstraints<Number> & constraints)
41{
42 if(level == 0)
43 {
44 // level of interest has been reached
45 unsigned int const dofs_per_face = face1->get_fe(0).dofs_per_face;
46
47 std::vector<dealii::types::global_dof_index> dofs_1(dofs_per_face);
48 std::vector<dealii::types::global_dof_index> dofs_2(dofs_per_face);
49
50 face1->get_mg_dof_indices(target_level, dofs_1, 0);
51 face2->get_mg_dof_indices(target_level, dofs_2, 0);
52
53 for(unsigned int i = 0; i < dofs_per_face; ++i)
54 {
55 if(constraints.can_store_line(dofs_2[i]) and constraints.can_store_line(dofs_1[i]) and
56 not(constraints.is_constrained(dofs_2[i])))
57 {
58 // constraint dof and ...
59 constraints.add_line(dofs_2[i]);
60 // specify type of constraint: equality (dof_2[i]=dof_1[j]*1.0)
61 constraints.add_entry(dofs_2[i], dofs_1[i], 1.);
62 }
63 }
64 }
65 else if(face1->has_children() and face2->has_children())
66 {
67 // recursively visit all subfaces
68 for(unsigned int c = 0; c < face1->n_children(); ++c)
69 {
70 add_periodicity_constraints<dim>(
71 level - 1, target_level, face1->child(c), face2->child(c), constraints);
72 }
73 }
74}
75
76
77template<int dim, typename Number>
78void
79add_periodicity_constraints(
80 dealii::DoFHandler<dim> const & dof_handler,
81 unsigned int const level,
82 std::vector<typename dealii::GridTools::PeriodicFacePair<
83 typename dealii::Triangulation<dim>::cell_iterator>> const & periodic_face_pairs_level0,
84 dealii::AffineConstraints<Number> & affine_constraints_own)
85{
86 // loop over all periodic face pairs of level 0
87 for(auto & it : periodic_face_pairs_level0)
88 {
89 // get reference to the cells on level 0 sharing the periodic face
90 typename dealii::DoFHandler<dim>::cell_iterator cell1(&dof_handler.get_triangulation(),
91 0,
92 it.cell[1]->index(),
93 &dof_handler);
94 typename dealii::DoFHandler<dim>::cell_iterator cell0(&dof_handler.get_triangulation(),
95 0,
96 it.cell[0]->index(),
97 &dof_handler);
98
99 // get reference to periodic faces on level and add recursively their
100 // subfaces on the given level
101 add_periodicity_constraints<dim, Number>(level,
102 level,
103 cell1->face(it.face_idx[1]),
104 cell0->face(it.face_idx[0]),
105 affine_constraints_own);
106 }
107}
108
109template<int dim, typename Number>
110void
111add_constraints(bool is_dg,
112 bool operator_is_singular,
113 dealii::DoFHandler<dim> const & dof_handler,
114 dealii::AffineConstraints<Number> & affine_constraints_own,
115 dealii::MGConstrainedDoFs const & mg_constrained_dofs,
116 std::vector<dealii::GridTools::PeriodicFacePair<
117 typename dealii::Triangulation<dim>::cell_iterator>> const & periodic_face_pairs,
118 unsigned int const level)
119{
120 if(is_dg)
121 {
122 // for DG: nothing to do
123 affine_constraints_own.close();
124 return;
125 }
126 // 0) clear old content (to be on the safe side)
127 affine_constraints_own.clear();
128
129 // ... and set local dofs
130 if(level == dealii::numbers::invalid_unsigned_int)
131 {
132 dealii::IndexSet const relevant_dofs =
133 dealii::DoFTools::extract_locally_relevant_dofs(dof_handler);
134 affine_constraints_own.reinit(dof_handler.locally_owned_dofs(), relevant_dofs);
135 }
136 else
137 {
138 dealii::IndexSet const relevant_dofs =
139 dealii::DoFTools::extract_locally_relevant_level_dofs(dof_handler, level);
140 affine_constraints_own.reinit(dof_handler.locally_owned_mg_dofs(level), relevant_dofs);
141 }
142
143 // 1) add periodic BCs
144 add_periodicity_constraints<dim, Number>(dof_handler,
145 level,
146 periodic_face_pairs,
147 affine_constraints_own);
148
149 // 2) add Dirichlet BCs
150 if(mg_constrained_dofs.have_boundary_indices())
151 affine_constraints_own.add_lines(mg_constrained_dofs.get_boundary_indices(level));
152
153 // constrain zeroth DoF in continuous case (the mean value constraint will
154 // be applied in the DG case). In case we have interface matrices, there are
155 // Dirichlet constraints on parts of the boundary and no such transformation
156 // is required.
157 if(operator_is_singular and affine_constraints_own.can_store_line(0))
158 {
159 // if dof 0 is constrained, it must be a periodic dof, so we take the
160 // value on the other side
161 dealii::types::global_dof_index line_index = 0;
162 while(true)
163 {
164 auto const * lines = affine_constraints_own.get_constraint_entries(line_index);
165
166 if(lines == 0)
167 {
168 affine_constraints_own.add_line(line_index);
169 // add the constraint back to the dealii::MGConstrainedDoFs field. This
170 // is potentially dangerous but we know what we are doing... ;-)
171 if(mg_constrained_dofs.have_boundary_indices() &&
172 level != dealii::numbers::invalid_unsigned_int)
173 {
174 if(mg_constrained_dofs.get_boundary_indices(level).size() != dof_handler.n_dofs(level))
175 const_cast<dealii::IndexSet &>(mg_constrained_dofs.get_boundary_indices(level))
176 .set_size(dof_handler.n_dofs(level));
177 const_cast<dealii::IndexSet &>(mg_constrained_dofs.get_boundary_indices(level))
178 .add_index(line_index);
179 }
180
181 break;
182 }
183 else
184 {
185 Assert(lines->size() == 1 and std::abs((*lines)[0].second - 1.) < 1e-15,
186 dealii::ExcMessage("Periodic index expected, bailing out"));
187
188 line_index = (*lines)[0].first;
189 }
190 }
191 }
192
193 affine_constraints_own.close();
194}
195
196} // anonymous namespace
197} // namespace ConstraintUtil
198} // namespace ExaDG
199
200#endif /* EXADG_SOLVERS_AND_PRECONDITIONERS_MULTIGRID_CONSTRAINTS_H_ */
Definition driver.cpp:33