Reference documentation for deal.II version 9.1.0-pre
mg_constrained_dofs.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_mg_constrained_dofs_h
17 #define dealii_mg_constrained_dofs_h
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/base/subscriptor.h>
22 
23 #include <deal.II/lac/affine_constraints.h>
24 
25 #include <deal.II/multigrid/mg_tools.h>
26 
27 #include <set>
28 #include <vector>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 template <int dim, int spacedim>
33 class DoFHandler;
34 
35 
43 {
44 public:
45  using size_dof = std::vector<std::set<types::global_dof_index>>::size_type;
65  template <int dim, int spacedim>
66  void
68 
79  template <int dim, int spacedim>
80  DEAL_II_DEPRECATED void
82  const std::map<types::boundary_id, const Function<spacedim> *>
83  & function_map,
84  const ComponentMask &component_mask = ComponentMask());
85 
96  template <int dim, int spacedim>
97  void
99  const DoFHandler<dim, spacedim> & dof,
100  const std::set<types::boundary_id> &boundary_ids,
101  const ComponentMask & component_mask = ComponentMask());
102 
106  void
107  clear();
108 
112  bool
113  is_boundary_index(const unsigned int level,
114  const types::global_dof_index index) const;
115 
119  bool
120  at_refinement_edge(const unsigned int level,
121  const types::global_dof_index index) const;
122 
123 
130  bool
131  is_interface_matrix_entry(const unsigned int level,
132  const types::global_dof_index i,
133  const types::global_dof_index j) const;
134 
141  const IndexSet &
142  get_boundary_indices(const unsigned int level) const;
143 
144 
149  const IndexSet &
150  get_refinement_edge_indices(unsigned int level) const;
151 
152 
156  bool
157  have_boundary_indices() const;
158 
164  get_level_constraint_matrix(const unsigned int level) const;
165 
166 private:
170  std::vector<IndexSet> boundary_indices;
171 
176  std::vector<IndexSet> refinement_edge_indices;
177 
182  std::vector<AffineConstraints<double>> level_constraints;
183 };
184 
185 
186 template <int dim, int spacedim>
187 inline void
189 {
190  boundary_indices.clear();
191  refinement_edge_indices.clear();
192  level_constraints.clear();
193 
194  const unsigned int nlevels = dof.get_triangulation().n_global_levels();
195 
196  // At this point level_constraint and refinement_edge_indices are empty.
197  level_constraints.resize(nlevels);
198  refinement_edge_indices.resize(nlevels);
199  for (unsigned int l = 0; l < nlevels; ++l)
200  {
201  IndexSet relevant_dofs;
202  DoFTools::extract_locally_relevant_level_dofs(dof, l, relevant_dofs);
203  level_constraints[l].reinit(relevant_dofs);
204 
205  // Loop through relevant cells and faces finding those which are periodic
206  // neighbors.
207  typename DoFHandler<dim, spacedim>::cell_iterator cell = dof.begin(l),
208  endc = dof.end(l);
209  for (; cell != endc; ++cell)
210  if (cell->level_subdomain_id() != numbers::artificial_subdomain_id)
211  {
212  for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
213  if (cell->has_periodic_neighbor(f))
214  {
215  if (cell->is_locally_owned_on_level())
216  {
217  Assert(
218  cell->periodic_neighbor(f)->level_subdomain_id() !=
220  ExcMessage(
221  "Periodic neighbor of a locally owned cell must either be owned or ghost."));
222  }
223  // Cell is a level-ghost and its neighbor is a
224  // level-artificial cell nothing to do here
225  else if (cell->periodic_neighbor(f)->level_subdomain_id() ==
227  {
228  Assert(cell->is_locally_owned_on_level() == false,
229  ExcInternalError());
230  continue;
231  }
232 
233  const unsigned int dofs_per_face =
234  cell->face(f)->get_fe(0).dofs_per_face;
235  std::vector<types::global_dof_index> dofs_1(dofs_per_face);
236  std::vector<types::global_dof_index> dofs_2(dofs_per_face);
237 
238  cell->periodic_neighbor(f)
239  ->face(cell->periodic_neighbor_face_no(f))
240  ->get_mg_dof_indices(l, dofs_1, 0);
241  cell->face(f)->get_mg_dof_indices(l, dofs_2, 0);
242  // Store periodicity information in the level constraint
243  // matrix Skip DoFs for which we've previously entered
244  // periodicity constraints already; this can happen, for
245  // example, for a vertex dof at a periodic boundary that we
246  // visit from more than one cell
247  for (unsigned int i = 0; i < dofs_per_face; ++i)
248  if (level_constraints[l].can_store_line(dofs_2[i]) &&
249  level_constraints[l].can_store_line(dofs_1[i]) &&
250  !level_constraints[l].is_constrained(dofs_2[i]) &&
251  !level_constraints[l].is_constrained(dofs_1[i]))
252  {
253  level_constraints[l].add_line(dofs_2[i]);
254  level_constraints[l].add_entry(dofs_2[i],
255  dofs_1[i],
256  1.);
257  }
258  }
259  }
260  level_constraints[l].close();
261 
262  // Initialize with empty IndexSet of correct size
264  }
265 
267 }
268 
269 
270 template <int dim, int spacedim>
271 inline void
273  const DoFHandler<dim, spacedim> & dof,
274  const std::map<types::boundary_id, const Function<spacedim> *> &function_map,
275  const ComponentMask &component_mask)
276 {
277  initialize(dof);
278 
279  // allocate an IndexSet for each global level. Contents will be
280  // overwritten inside make_boundary_list.
281  const unsigned int n_levels = dof.get_triangulation().n_global_levels();
282  // At this point boundary_indices is empty.
283  boundary_indices.resize(n_levels);
284 
286  function_map,
288  component_mask);
289 }
290 
291 
292 template <int dim, int spacedim>
293 inline void
295  const DoFHandler<dim, spacedim> & dof,
296  const std::set<types::boundary_id> &boundary_ids,
297  const ComponentMask & component_mask)
298 {
299  // allocate an IndexSet for each global level. Contents will be
300  // overwritten inside make_boundary_list.
301  const unsigned int n_levels = dof.get_triangulation().n_global_levels();
302  Assert(boundary_indices.size() == 0 || boundary_indices.size() == n_levels,
303  ExcInternalError());
304  boundary_indices.resize(n_levels);
305 
307  boundary_ids,
309  component_mask);
310 }
311 
312 
313 inline void
315 {
316  boundary_indices.clear();
317  refinement_edge_indices.clear();
318 }
319 
320 
321 inline bool
322 MGConstrainedDoFs::is_boundary_index(const unsigned int level,
323  const types::global_dof_index index) const
324 {
325  if (boundary_indices.size() == 0)
326  return false;
327 
328  AssertIndexRange(level, boundary_indices.size());
329  return boundary_indices[level].is_element(index);
330 }
331 
332 inline bool
333 MGConstrainedDoFs::at_refinement_edge(const unsigned int level,
334  const types::global_dof_index index) const
335 {
337 
338  return refinement_edge_indices[level].is_element(index);
339 }
340 
341 inline bool
343  const unsigned int level,
344  const types::global_dof_index i,
345  const types::global_dof_index j) const
346 {
347  const IndexSet &interface_dofs_on_level =
348  this->get_refinement_edge_indices(level);
349 
350  return interface_dofs_on_level.is_element(i) // at_refinement_edge(i)
351  && !interface_dofs_on_level.is_element(j) // !at_refinement_edge(j)
352  && !this->is_boundary_index(level, i) // !on_boundary(i)
353  && !this->is_boundary_index(level, j); // !on_boundary(j)
354 }
355 
356 
357 
358 inline const IndexSet &
359 MGConstrainedDoFs::get_boundary_indices(const unsigned int level) const
360 {
361  AssertIndexRange(level, boundary_indices.size());
362  return boundary_indices[level];
363 }
364 
365 
366 
367 inline const IndexSet &
369 {
371  return refinement_edge_indices[level];
372 }
373 
374 
375 
376 inline bool
378 {
379  return boundary_indices.size() != 0;
380 }
381 
382 
383 
384 inline const AffineConstraints<double> &
385 MGConstrainedDoFs::get_level_constraint_matrix(const unsigned int level) const
386 {
387  AssertIndexRange(level, level_constraints.size());
388  return level_constraints[level];
389 }
390 
391 
392 
393 DEAL_II_NAMESPACE_CLOSE
394 
395 #endif
bool is_boundary_index(const unsigned int level, const types::global_dof_index index) const
bool is_interface_matrix_entry(const unsigned int level, const types::global_dof_index i, const types::global_dof_index j) const
bool at_refinement_edge(const unsigned int level, const types::global_dof_index index) const
const AffineConstraints< double > & get_level_constraint_matrix(const unsigned int level) const
const IndexSet & get_boundary_indices(const unsigned int level) const
#define AssertIndexRange(index, range)
Definition: exceptions.h:1407
std::vector< IndexSet > boundary_indices
const Triangulation< dim, spacedim > & get_triangulation() const
void extract_locally_relevant_level_dofs(const DoFHandlerType &dof_handler, const unsigned int level, IndexSet &dof_set)
Definition: dof_tools.cc:1185
unsigned long long int global_dof_index
Definition: types.h:72
static::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1227
types::global_dof_index n_dofs() const
cell_iterator end() const
Definition: dof_handler.cc:959
void make_zero_boundary_constraints(const DoFHandler< dim, spacedim > &dof, const std::set< types::boundary_id > &boundary_ids, const ComponentMask &component_mask=ComponentMask())
std::vector< IndexSet > refinement_edge_indices
const types::subdomain_id artificial_subdomain_id
Definition: types.h:272
cell_iterator begin(const unsigned int level=0) const
Definition: dof_handler.cc:930
const IndexSet & get_refinement_edge_indices(unsigned int level) const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
bool have_boundary_indices() const
void extract_inner_interface_dofs(const DoFHandler< dim, spacedim > &mg_dof_handler, std::vector< IndexSet > &interface_dofs)
Definition: mg_tools.cc:1536
bool is_element(const size_type index) const
Definition: index_set.h:1676
void make_boundary_list(const DoFHandler< dim, spacedim > &mg_dof, const std::map< types::boundary_id, const Function< spacedim > * > &function_map, std::vector< std::set< types::global_dof_index >> &boundary_indices, const ComponentMask &component_mask=ComponentMask())
Definition: mg_tools.cc:1255
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:268
std::vector< AffineConstraints< double > > level_constraints
unsigned int boundary_id
Definition: types.h:111
static::ExceptionBase & ExcInternalError()
void initialize(const DoFHandler< dim, spacedim > &dof)