Reference documentation for deal.II version Git 9e617570d7 2020-09-22 20:39:52 -0400
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
mg_constrained_dofs.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2020 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 
22 
24 
26 
27 #include <set>
28 #include <vector>
29 
31 
32 // Forward declaration
33 #ifndef DOXYGEN
34 template <int dim, int spacedim>
35 class DoFHandler;
36 #endif
37 
38 
46 {
47 public:
48  using size_dof = std::vector<std::set<types::global_dof_index>>::size_type;
68  template <int dim, int spacedim>
69  void
71 
82  template <int dim, int spacedim>
83  void
85  const DoFHandler<dim, spacedim> & dof,
86  const std::set<types::boundary_id> &boundary_ids,
87  const ComponentMask & component_mask = ComponentMask());
88 
103  void
104  add_user_constraints(const unsigned int level,
105  const AffineConstraints<double> &constraints_on_level);
106 
121  template <int dim, int spacedim>
122  void
124  const types::boundary_id bid,
125  const unsigned int first_vector_component);
126 
130  void
131  clear();
132 
136  bool
137  is_boundary_index(const unsigned int level,
138  const types::global_dof_index index) const;
139 
143  bool
144  at_refinement_edge(const unsigned int level,
145  const types::global_dof_index index) const;
146 
147 
154  bool
155  is_interface_matrix_entry(const unsigned int level,
156  const types::global_dof_index i,
157  const types::global_dof_index j) const;
158 
165  const IndexSet &
166  get_boundary_indices(const unsigned int level) const;
167 
168 
173  const IndexSet &
174  get_refinement_edge_indices(unsigned int level) const;
175 
176 
180  bool
181  have_boundary_indices() const;
182 
188  get_level_constraints(const unsigned int level) const;
189 
198  get_user_constraint_matrix(const unsigned int level) const;
199 
200 private:
204  std::vector<IndexSet> boundary_indices;
205 
210  std::vector<IndexSet> refinement_edge_indices;
211 
216  std::vector<AffineConstraints<double>> level_constraints;
217 
221  std::vector<AffineConstraints<double>> user_constraints;
222 };
223 
224 
225 template <int dim, int spacedim>
226 inline void
228 {
229  boundary_indices.clear();
230  refinement_edge_indices.clear();
231  level_constraints.clear();
232  user_constraints.clear();
233 
234  const unsigned int nlevels = dof.get_triangulation().n_global_levels();
235 
236  // At this point level_constraint and refinement_edge_indices are empty.
237  refinement_edge_indices.resize(nlevels);
238  level_constraints.resize(nlevels);
239  user_constraints.resize(nlevels);
240  for (unsigned int l = 0; l < nlevels; ++l)
241  {
242  IndexSet relevant_dofs;
244  level_constraints[l].reinit(relevant_dofs);
245 
246  // Loop through relevant cells and faces finding those which are periodic
247  // neighbors.
248  typename DoFHandler<dim, spacedim>::cell_iterator cell = dof.begin(l),
249  endc = dof.end(l);
250  for (; cell != endc; ++cell)
251  if (cell->level_subdomain_id() != numbers::artificial_subdomain_id)
252  {
253  for (auto f : cell->face_indices())
254  if (cell->has_periodic_neighbor(f) &&
255  cell->periodic_neighbor(f)->level() == cell->level())
256  {
257  if (cell->is_locally_owned_on_level())
258  {
259  Assert(
260  cell->periodic_neighbor(f)->level_subdomain_id() !=
262  ExcMessage(
263  "Periodic neighbor of a locally owned cell must either be owned or ghost."));
264  }
265  // Cell is a level-ghost and its neighbor is a
266  // level-artificial cell nothing to do here
267  else if (cell->periodic_neighbor(f)->level_subdomain_id() ==
269  {
270  Assert(cell->is_locally_owned_on_level() == false,
271  ExcInternalError());
272  continue;
273  }
274 
275  const unsigned int dofs_per_face =
276  cell->face(f)->get_fe(0).n_dofs_per_face();
277  std::vector<types::global_dof_index> dofs_1(dofs_per_face);
278  std::vector<types::global_dof_index> dofs_2(dofs_per_face);
279 
280  cell->periodic_neighbor(f)
281  ->face(cell->periodic_neighbor_face_no(f))
282  ->get_mg_dof_indices(l, dofs_1, 0);
283  cell->face(f)->get_mg_dof_indices(l, dofs_2, 0);
284  // Store periodicity information in the level
285  // AffineConstraints object. Skip DoFs for which we've
286  // previously entered periodicity constraints already; this
287  // can happen, for example, for a vertex dof at a periodic
288  // boundary that we visit from more than one cell
289  for (unsigned int i = 0; i < dofs_per_face; ++i)
290  if (level_constraints[l].can_store_line(dofs_2[i]) &&
291  level_constraints[l].can_store_line(dofs_1[i]) &&
292  !level_constraints[l].is_constrained(dofs_2[i]) &&
293  !level_constraints[l].is_constrained(dofs_1[i]))
294  {
295  level_constraints[l].add_line(dofs_2[i]);
296  level_constraints[l].add_entry(dofs_2[i],
297  dofs_1[i],
298  1.);
299  }
300  }
301  }
302  level_constraints[l].close();
303 
304  // Initialize with empty IndexSet of correct size
306  }
307 
309 }
310 
311 
312 
313 template <int dim, int spacedim>
314 inline void
316  const DoFHandler<dim, spacedim> & dof,
317  const std::set<types::boundary_id> &boundary_ids,
318  const ComponentMask & component_mask)
319 {
320  // allocate an IndexSet for each global level. Contents will be
321  // overwritten inside make_boundary_list.
322  const unsigned int n_levels = dof.get_triangulation().n_global_levels();
323  Assert(boundary_indices.size() == 0 || boundary_indices.size() == n_levels,
324  ExcInternalError());
325  boundary_indices.resize(n_levels);
326 
328  boundary_ids,
330  component_mask);
331 }
332 
333 
334 template <int dim, int spacedim>
335 inline void
337  const DoFHandler<dim, spacedim> &dof,
338  const types::boundary_id bid,
339  const unsigned int first_vector_component)
340 {
341  // For a given boundary id, find which vector component is on the boundary
342  // and set a zero boundary constraint for those degrees of freedom.
343  const unsigned int n_components = dof.get_fe_collection().n_components();
344  AssertIndexRange(first_vector_component + dim - 1, n_components);
345 
346  ComponentMask comp_mask(n_components, false);
347 
348 
350  face = dof.get_triangulation().begin_face(),
351  endf = dof.get_triangulation().end_face();
352  for (; face != endf; ++face)
353  if (face->at_boundary() && face->boundary_id() == bid)
354  for (unsigned int d = 0; d < dim; ++d)
355  {
356  Tensor<1, dim, double> unit_vec;
357  unit_vec[d] = 1.0;
358 
359  const Tensor<1, dim> normal_vec =
360  face->get_manifold().normal_vector(face, face->center());
361 
362  if (std::abs(std::abs(unit_vec * normal_vec) - 1.0) < 1e-10)
363  comp_mask.set(d + first_vector_component, true);
364  else
365  Assert(
366  std::abs(unit_vec * normal_vec) < 1e-10,
367  ExcMessage(
368  "We can currently only support no normal flux conditions "
369  "for a specific boundary id if all faces are normal to the "
370  "x, y, or z axis."));
371  }
372 
373  Assert(comp_mask.n_selected_components() == 1,
374  ExcMessage(
375  "We can currently only support no normal flux conditions "
376  "for a specific boundary id if all faces are facing in the "
377  "same direction, i.e., a boundary normal to the x-axis must "
378  "have a different boundary id than a boundary normal to the "
379  "y- or z-axis and so on. If the mesh here was produced using "
380  "GridGenerator::..., setting colorize=true during mesh generation "
381  "and calling make_no_normal_flux_constraints() for each no normal "
382  "flux boundary will fulfill the condition."));
383 
384  this->make_zero_boundary_constraints(dof, {bid}, comp_mask);
385 }
386 
387 
388 inline void
390  const unsigned int level,
391  const AffineConstraints<double> &constraints_on_level)
392 {
393  AssertIndexRange(level, user_constraints.size());
394 
395  // Get the relevant DoFs from level_constraints if
396  // the user constraint matrix has not been initialized
397  if (user_constraints[level].get_local_lines().size() == 0)
398  user_constraints[level].reinit(level_constraints[level].get_local_lines());
399 
400  user_constraints[level].merge(
401  constraints_on_level,
403  user_constraints[level].close();
404 }
405 
406 
407 inline void
409 {
410  boundary_indices.clear();
411  refinement_edge_indices.clear();
412  user_constraints.clear();
413 }
414 
415 
416 inline bool
418  const types::global_dof_index index) const
419 {
420  if (boundary_indices.size() == 0)
421  return false;
422 
423  AssertIndexRange(level, boundary_indices.size());
424  return boundary_indices[level].is_element(index);
425 }
426 
427 inline bool
429  const types::global_dof_index index) const
430 {
432 
433  return refinement_edge_indices[level].is_element(index);
434 }
435 
436 inline bool
438  const unsigned int level,
439  const types::global_dof_index i,
440  const types::global_dof_index j) const
441 {
442  const IndexSet &interface_dofs_on_level =
443  this->get_refinement_edge_indices(level);
444 
445  return interface_dofs_on_level.is_element(i) // at_refinement_edge(i)
446  && !interface_dofs_on_level.is_element(j) // !at_refinement_edge(j)
447  && !this->is_boundary_index(level, i) // !on_boundary(i)
448  && !this->is_boundary_index(level, j); // !on_boundary(j)
449 }
450 
451 
452 
453 inline const IndexSet &
455 {
456  AssertIndexRange(level, boundary_indices.size());
457  return boundary_indices[level];
458 }
459 
460 
461 
462 inline const IndexSet &
464 {
467 }
468 
469 
470 
471 inline bool
473 {
474  return boundary_indices.size() != 0;
475 }
476 
477 
478 
479 inline const AffineConstraints<double> &
481 {
482  AssertIndexRange(level, level_constraints.size());
483  return level_constraints[level];
484 }
485 
486 
487 
488 inline const AffineConstraints<double> &
490 {
491  AssertIndexRange(level, user_constraints.size());
492  return user_constraints[level];
493 }
494 
495 
496 
498 
499 #endif
const Triangulation< dim, spacedim > & get_triangulation() const
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:618
cell_iterator begin(const unsigned int level=0) const
types::global_dof_index size_type
Definition: cuda_kernels.h:45
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
bool is_boundary_index(const unsigned int level, const types::global_dof_index index) const
void make_no_normal_flux_constraints(const DoFHandler< dim, spacedim > &dof, const types::boundary_id bid, const unsigned int first_vector_component)
cell_iterator end() const
#define AssertIndexRange(index, range)
Definition: exceptions.h:1636
std::vector< IndexSet > boundary_indices
const AffineConstraints< double > & get_level_constraints(const unsigned int level) const
std::vector< std::set< types::global_dof_index > >::size_type size_dof
void set(const unsigned int index, const bool value)
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
void extract_locally_relevant_level_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level, IndexSet &dof_set)
Definition: dof_tools.cc:1180
bool is_interface_matrix_entry(const unsigned int level, const types::global_dof_index i, const types::global_dof_index j) const
static ::ExceptionBase & ExcMessage(std::string arg1)
const IndexSet & get_boundary_indices(const unsigned int level) const
void add_user_constraints(const unsigned int level, const AffineConstraints< double > &constraints_on_level)
#define Assert(cond, exc)
Definition: exceptions.h:1411
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:363
unsigned int level
Definition: grid_out.cc:4338
types::global_dof_index n_dofs() const
const IndexSet & get_refinement_edge_indices(unsigned int level) const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
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:293
bool at_refinement_edge(const unsigned int level, const types::global_dof_index index) const
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:1252
Definition: tensor.h:448
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:362
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
void extract_inner_interface_dofs(const DoFHandler< dim, spacedim > &mg_dof_handler, std::vector< IndexSet > &interface_dofs)
Definition: mg_tools.cc:1461
bool is_element(const size_type index) const
Definition: index_set.h:1763
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const
Definition: manifold.cc:237
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:371
bool have_boundary_indices() const
std::vector< AffineConstraints< double > > level_constraints
const AffineConstraints< double > & get_user_constraint_matrix(const unsigned int level) const
const Manifold< dim, spacedim > & get_manifold(const types::manifold_id number) const
Definition: tria.cc:9482
std::vector< AffineConstraints< double > > user_constraints
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
static ::ExceptionBase & ExcInternalError()
void initialize(const DoFHandler< dim, spacedim > &dof)