Reference documentation for deal.II version GIT 3d869ba6cd 2022-05-16 18:15:01+00:00
\(\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 
23 
25 
27 
28 #include <set>
29 #include <vector>
30 
32 
33 // Forward declaration
34 #ifndef DOXYGEN
35 template <int dim, int spacedim>
36 class DoFHandler;
37 #endif
38 
39 
47 {
48 public:
49  using size_dof = std::vector<std::set<types::global_dof_index>>::size_type;
75  template <int dim, int spacedim>
76  void
78  const MGLevelObject<IndexSet> & level_relevant_dofs =
80 
91  template <int dim, int spacedim>
92  void
94  const DoFHandler<dim, spacedim> & dof,
95  const std::set<types::boundary_id> &boundary_ids,
96  const ComponentMask & component_mask = ComponentMask());
97 
104  template <int dim, int spacedim>
105  void
107  const unsigned int level,
108  const IndexSet & boundary_indices);
109 
124  void
125  add_user_constraints(const unsigned int level,
126  const AffineConstraints<double> &constraints_on_level);
127 
142  template <int dim, int spacedim>
143  void
145  const types::boundary_id bid,
146  const unsigned int first_vector_component);
147 
151  void
152  clear();
153 
157  bool
158  is_boundary_index(const unsigned int level,
159  const types::global_dof_index index) const;
160 
164  bool
165  at_refinement_edge(const unsigned int level,
166  const types::global_dof_index index) const;
167 
168 
175  bool
176  is_interface_matrix_entry(const unsigned int level,
177  const types::global_dof_index i,
178  const types::global_dof_index j) const;
179 
186  const IndexSet &
187  get_boundary_indices(const unsigned int level) const;
188 
189 
194  const IndexSet &
195  get_refinement_edge_indices(unsigned int level) const;
196 
197 
201  bool
202  have_boundary_indices() const;
203 
209  get_level_constraints(const unsigned int level) const;
210 
219  get_user_constraint_matrix(const unsigned int level) const;
220 
221 private:
225  std::vector<IndexSet> boundary_indices;
226 
231  std::vector<IndexSet> refinement_edge_indices;
232 
237  std::vector<AffineConstraints<double>> level_constraints;
238 
242  std::vector<AffineConstraints<double>> user_constraints;
243 };
244 
245 
246 template <int dim, int spacedim>
247 inline void
249  const DoFHandler<dim, spacedim> &dof,
250  const MGLevelObject<IndexSet> & level_relevant_dofs)
251 {
252  boundary_indices.clear();
253  refinement_edge_indices.clear();
254  level_constraints.clear();
255  user_constraints.clear();
256 
257  const unsigned int nlevels = dof.get_triangulation().n_global_levels();
258  const unsigned int min_level = level_relevant_dofs.min_level();
259  const unsigned int max_level = (level_relevant_dofs.max_level() == 0) ?
260  nlevels - 1 :
261  level_relevant_dofs.max_level();
262  const bool user_level_dofs =
263  (level_relevant_dofs.max_level() == 0) ? false : true;
264 
265  // At this point level_constraint and refinement_edge_indices are empty.
266  refinement_edge_indices.resize(nlevels);
267  level_constraints.resize(nlevels);
268  user_constraints.resize(nlevels);
269  for (unsigned int l = min_level; l <= max_level; ++l)
270  {
271  if (user_level_dofs)
272  {
273  level_constraints[l].reinit(level_relevant_dofs[l]);
274  }
275  else
276  {
277  const IndexSet relevant_dofs =
279  level_constraints[l].reinit(relevant_dofs);
280  }
281 
282  // Loop through relevant cells and faces finding those which are periodic
283  // neighbors.
284  typename DoFHandler<dim, spacedim>::cell_iterator cell = dof.begin(l),
285  endc = dof.end(l);
286  for (; cell != endc; ++cell)
287  if (cell->level_subdomain_id() != numbers::artificial_subdomain_id)
288  {
289  for (auto f : cell->face_indices())
290  if (cell->has_periodic_neighbor(f) &&
291  cell->periodic_neighbor(f)->level() == cell->level())
292  {
293  if (cell->is_locally_owned_on_level())
294  {
295  Assert(
296  cell->periodic_neighbor(f)->level_subdomain_id() !=
298  ExcMessage(
299  "Periodic neighbor of a locally owned cell must either be owned or ghost."));
300  }
301  // Cell is a level-ghost and its neighbor is a
302  // level-artificial cell nothing to do here
303  else if (cell->periodic_neighbor(f)->level_subdomain_id() ==
305  {
306  Assert(cell->is_locally_owned_on_level() == false,
307  ExcInternalError());
308  continue;
309  }
310 
311  const unsigned int dofs_per_face =
312  dof.get_fe(0).n_dofs_per_face(f);
313  std::vector<types::global_dof_index> dofs_1(dofs_per_face);
314  std::vector<types::global_dof_index> dofs_2(dofs_per_face);
315 
316  cell->periodic_neighbor(f)
317  ->face(cell->periodic_neighbor_face_no(f))
318  ->get_mg_dof_indices(l, dofs_1, 0);
319  cell->face(f)->get_mg_dof_indices(l, dofs_2, 0);
320  // Store periodicity information in the level
321  // AffineConstraints object. Skip DoFs for which we've
322  // previously entered periodicity constraints already; this
323  // can happen, for example, for a vertex dof at a periodic
324  // boundary that we visit from more than one cell
325  for (unsigned int i = 0; i < dofs_per_face; ++i)
326  if (level_constraints[l].can_store_line(dofs_2[i]) &&
327  level_constraints[l].can_store_line(dofs_1[i]) &&
328  !level_constraints[l].is_constrained(dofs_2[i]) &&
329  !level_constraints[l].is_constrained(dofs_1[i]))
330  {
331  level_constraints[l].add_line(dofs_2[i]);
332  level_constraints[l].add_entry(dofs_2[i],
333  dofs_1[i],
334  1.);
335  }
336  }
337  }
338  level_constraints[l].close();
339 
340  // Initialize with empty IndexSet of correct size
342  }
343 
345 }
346 
347 
348 template <int dim, int spacedim>
349 inline void
351  const DoFHandler<dim, spacedim> & dof,
352  const std::set<types::boundary_id> &boundary_ids,
353  const ComponentMask & component_mask)
354 {
355  // allocate an IndexSet for each global level. Contents will be
356  // overwritten inside make_boundary_list.
357  const unsigned int n_levels = dof.get_triangulation().n_global_levels();
358  Assert(boundary_indices.size() == 0 || boundary_indices.size() == n_levels,
359  ExcInternalError());
360  boundary_indices.resize(n_levels);
361 
363  boundary_ids,
365  component_mask);
366 }
367 
368 
369 
370 template <int dim, int spacedim>
371 inline void
373  const unsigned int level,
374  const IndexSet &level_boundary_indices)
375 {
376  const unsigned int n_levels = dof.get_triangulation().n_global_levels();
377  if (boundary_indices.size() == 0)
378  {
379  boundary_indices.resize(n_levels);
380  for (unsigned int i = 0; i < n_levels; ++i)
381  boundary_indices[i] = IndexSet(dof.n_dofs(i));
382  }
383  AssertDimension(boundary_indices.size(), n_levels);
384  boundary_indices[level].add_indices(level_boundary_indices);
385 }
386 
387 
388 
389 template <int dim, int spacedim>
390 inline void
392  const DoFHandler<dim, spacedim> &dof,
393  const types::boundary_id bid,
394  const unsigned int first_vector_component)
395 {
396  // For a given boundary id, find which vector component is on the boundary
397  // and set a zero boundary constraint for those degrees of freedom.
398  const unsigned int n_components = dof.get_fe_collection().n_components();
399  AssertIndexRange(first_vector_component + dim - 1, n_components);
400 
401  ComponentMask comp_mask(n_components, false);
402 
403 
405  face = dof.get_triangulation().begin_face(),
406  endf = dof.get_triangulation().end_face();
407  for (; face != endf; ++face)
408  if (face->at_boundary() && face->boundary_id() == bid)
409  for (unsigned int d = 0; d < dim; ++d)
410  {
411  Tensor<1, dim, double> unit_vec;
412  unit_vec[d] = 1.0;
413 
414  const Tensor<1, dim> normal_vec =
415  face->get_manifold().normal_vector(face, face->center());
416 
417  if (std::abs(std::abs(unit_vec * normal_vec) - 1.0) < 1e-10)
418  comp_mask.set(d + first_vector_component, true);
419  else
420  Assert(
421  std::abs(unit_vec * normal_vec) < 1e-10,
422  ExcMessage(
423  "We can currently only support no normal flux conditions "
424  "for a specific boundary id if all faces are normal to the "
425  "x, y, or z axis."));
426  }
427 
428  Assert(comp_mask.n_selected_components() == 1,
429  ExcMessage(
430  "We can currently only support no normal flux conditions "
431  "for a specific boundary id if all faces are facing in the "
432  "same direction, i.e., a boundary normal to the x-axis must "
433  "have a different boundary id than a boundary normal to the "
434  "y- or z-axis and so on. If the mesh here was produced using "
435  "GridGenerator::..., setting colorize=true during mesh generation "
436  "and calling make_no_normal_flux_constraints() for each no normal "
437  "flux boundary will fulfill the condition."));
438 
439  this->make_zero_boundary_constraints(dof, {bid}, comp_mask);
440 }
441 
442 
443 inline void
445  const unsigned int level,
446  const AffineConstraints<double> &constraints_on_level)
447 {
449 
450  // Get the relevant DoFs from level_constraints if
451  // the user constraint matrix has not been initialized
452  if (user_constraints[level].get_local_lines().size() == 0)
453  user_constraints[level].reinit(level_constraints[level].get_local_lines());
454 
455  user_constraints[level].merge(
456  constraints_on_level,
458  user_constraints[level].close();
459 }
460 
461 
462 inline void
464 {
465  boundary_indices.clear();
466  refinement_edge_indices.clear();
467  user_constraints.clear();
468 }
469 
470 
471 inline bool
473  const types::global_dof_index index) const
474 {
475  if (boundary_indices.size() == 0)
476  return false;
477 
479  return boundary_indices[level].is_element(index);
480 }
481 
482 inline bool
484  const types::global_dof_index index) const
485 {
487 
488  return refinement_edge_indices[level].is_element(index);
489 }
490 
491 inline bool
493  const unsigned int level,
494  const types::global_dof_index i,
495  const types::global_dof_index j) const
496 {
497  const IndexSet &interface_dofs_on_level =
498  this->get_refinement_edge_indices(level);
499 
500  return interface_dofs_on_level.is_element(i) // at_refinement_edge(i)
501  && !interface_dofs_on_level.is_element(j) // !at_refinement_edge(j)
502  && !this->is_boundary_index(level, i) // !on_boundary(i)
503  && !this->is_boundary_index(level, j); // !on_boundary(j)
504 }
505 
506 
507 
508 inline const IndexSet &
510 {
512  return boundary_indices[level];
513 }
514 
515 
516 
517 inline const IndexSet &
519 {
522 }
523 
524 
525 
526 inline bool
528 {
529  return boundary_indices.size() != 0;
530 }
531 
532 
533 
534 inline const AffineConstraints<double> &
536 {
538  return level_constraints[level];
539 }
540 
541 
542 
543 inline const AffineConstraints<double> &
545 {
547  return user_constraints[level];
548 }
549 
550 
551 
553 
554 #endif
void set(const unsigned int index, const bool value)
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
cell_iterator end() const
const Triangulation< dim, spacedim > & get_triangulation() const
types::global_dof_index n_dofs() const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
cell_iterator begin(const unsigned int level=0) const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
bool is_element(const size_type index) const
Definition: index_set.h:1767
std::vector< IndexSet > refinement_edge_indices
std::vector< AffineConstraints< double > > user_constraints
void make_no_normal_flux_constraints(const DoFHandler< dim, spacedim > &dof, const types::boundary_id bid, const unsigned int first_vector_component)
bool is_interface_matrix_entry(const unsigned int level, const types::global_dof_index i, const types::global_dof_index j) const
std::vector< IndexSet > boundary_indices
void add_boundary_indices(const DoFHandler< dim, spacedim > &dof, const unsigned int level, const IndexSet &boundary_indices)
bool at_refinement_edge(const unsigned int level, const types::global_dof_index index) const
void initialize(const DoFHandler< dim, spacedim > &dof, const MGLevelObject< IndexSet > &level_relevant_dofs=MGLevelObject< IndexSet >())
bool is_boundary_index(const unsigned int level, const types::global_dof_index index) const
std::vector< AffineConstraints< double > > level_constraints
bool have_boundary_indices() const
void add_user_constraints(const unsigned int level, const AffineConstraints< double > &constraints_on_level)
const IndexSet & get_refinement_edge_indices(unsigned int level) const
const AffineConstraints< double > & get_user_constraint_matrix(const unsigned int level) const
std::vector< std::set< types::global_dof_index > >::size_type size_dof
void make_zero_boundary_constraints(const DoFHandler< dim, spacedim > &dof, const std::set< types::boundary_id > &boundary_ids, const ComponentMask &component_mask=ComponentMask())
const AffineConstraints< double > & get_level_constraints(const unsigned int level) const
const IndexSet & get_boundary_indices(const unsigned int level) const
unsigned int max_level() const
unsigned int min_level() const
face_iterator end_face() const
virtual unsigned int n_global_levels() const
face_iterator begin_face() const
unsigned int n_components() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:416
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:417
unsigned int level
Definition: grid_out.cc:4606
static ::ExceptionBase & ExcInternalError()
#define Assert(cond, exc)
Definition: exceptions.h:1473
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1667
#define AssertIndexRange(index, range)
Definition: exceptions.h:1732
static ::ExceptionBase & ExcMessage(std::string arg1)
IndexSet extract_locally_relevant_level_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level)
Definition: dof_tools.cc:1194
types::global_dof_index size_type
Definition: cuda_kernels.h:45
void extract_inner_interface_dofs(const DoFHandler< dim, spacedim > &mg_dof_handler, std::vector< IndexSet > &interface_dofs)
Definition: mg_tools.cc:1478
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:1272
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
const types::subdomain_id artificial_subdomain_id
Definition: types.h:298
unsigned int global_dof_index
Definition: types.h:76
unsigned int boundary_id
Definition: types.h:129