22#ifndef EXADG_SOLVERS_AND_PRECONDITIONERS_MULTIGRID_CONSTRAINTS_H_
23#define EXADG_SOLVERS_AND_PRECONDITIONERS_MULTIGRID_CONSTRAINTS_H_
26#include <deal.II/multigrid/mg_constrained_dofs.h>
30namespace ConstraintUtil
34template<
int dim,
typename Number>
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)
45 unsigned int const dofs_per_face = face1->get_fe(0).dofs_per_face;
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);
50 face1->get_mg_dof_indices(target_level, dofs_1, 0);
51 face2->get_mg_dof_indices(target_level, dofs_2, 0);
53 for(
unsigned int i = 0; i < dofs_per_face; ++i)
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])))
59 constraints.add_line(dofs_2[i]);
61 constraints.add_entry(dofs_2[i], dofs_1[i], 1.);
65 else if(face1->has_children() and face2->has_children())
68 for(
unsigned int c = 0; c < face1->n_children(); ++c)
70 add_periodicity_constraints<dim>(
71 level - 1, target_level, face1->child(c), face2->child(c), constraints);
77template<
int dim,
typename Number>
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)
87 for(
auto & it : periodic_face_pairs_level0)
90 typename dealii::DoFHandler<dim>::cell_iterator cell1(&dof_handler.get_triangulation(),
94 typename dealii::DoFHandler<dim>::cell_iterator cell0(&dof_handler.get_triangulation(),
101 add_periodicity_constraints<dim, Number>(level,
103 cell1->face(it.face_idx[1]),
104 cell0->face(it.face_idx[0]),
105 affine_constraints_own);
109template<
int dim,
typename Number>
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)
123 affine_constraints_own.close();
127 affine_constraints_own.clear();
130 if(level == dealii::numbers::invalid_unsigned_int)
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);
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);
144 add_periodicity_constraints<dim, Number>(dof_handler,
147 affine_constraints_own);
150 if(mg_constrained_dofs.have_boundary_indices())
151 affine_constraints_own.add_lines(mg_constrained_dofs.get_boundary_indices(level));
157 if(operator_is_singular and affine_constraints_own.can_store_line(0))
161 dealii::types::global_dof_index line_index = 0;
164 auto const * lines = affine_constraints_own.get_constraint_entries(line_index);
168 affine_constraints_own.add_line(line_index);
171 if(mg_constrained_dofs.have_boundary_indices() &&
172 level != dealii::numbers::invalid_unsigned_int)
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);
185 Assert(lines->size() == 1 and std::abs((*lines)[0].second - 1.) < 1e-15,
186 dealii::ExcMessage(
"Periodic index expected, bailing out"));
188 line_index = (*lines)[0].first;
193 affine_constraints_own.close();