Reference documentation for deal.II version 9.0.0
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 at
12 // the top level of the deal.II distribution.
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 #include <deal.II/base/subscriptor.h>
21 
22 #include <deal.II/multigrid/mg_tools.h>
23 
24 #include <vector>
25 #include <set>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 template <int dim, int spacedim> class DoFHandler;
30 template <int dim, typename Number> struct FunctionMap;
31 
32 
40 {
41 public:
42 
43  typedef std::vector<std::set<types::global_dof_index> >::size_type size_dof;
63  template <int dim, int spacedim>
64  void initialize (const DoFHandler<dim,spacedim> &dof);
65 
76  template <int dim, int spacedim>
77  DEAL_II_DEPRECATED
78  void initialize(const DoFHandler<dim,spacedim> &dof,
79  const typename FunctionMap<dim>::type &function_map,
80  const ComponentMask &component_mask = ComponentMask());
81 
92  template <int dim, int spacedim>
94  const std::set<types::boundary_id> &boundary_ids,
95  const ComponentMask &component_mask = ComponentMask());
96 
100  void clear();
101 
105  bool is_boundary_index (const unsigned int level,
106  const types::global_dof_index index) const;
107 
111  bool at_refinement_edge (const unsigned int level,
112  const types::global_dof_index index) const;
113 
114 
121  bool is_interface_matrix_entry (const unsigned int level,
122  const types::global_dof_index i,
123  const types::global_dof_index j) const;
124 
131  const IndexSet &
132  get_boundary_indices (const unsigned int level) const;
133 
134 
139  const IndexSet &
140  get_refinement_edge_indices (unsigned int level) const;
141 
142 
146  bool have_boundary_indices () const;
147 
152  const ConstraintMatrix &
153  get_level_constraint_matrix (const unsigned int level) const;
154 
155 private:
156 
160  std::vector<IndexSet> boundary_indices;
161 
166  std::vector<IndexSet> refinement_edge_indices;
167 
172  std::vector<ConstraintMatrix> level_constraints;
173 
174 };
175 
176 
177 template <int dim, int spacedim>
178 inline
179 void
181 {
182  boundary_indices.clear();
183  refinement_edge_indices.clear();
184  level_constraints.clear();
185 
186  const unsigned int nlevels = dof.get_triangulation().n_global_levels();
187 
188  //At this point level_constraint and refinement_edge_indices are empty.
189  level_constraints.resize(nlevels);
190  refinement_edge_indices.resize(nlevels);
191  for (unsigned int l=0; l<nlevels; ++l)
192  {
193  IndexSet relevant_dofs;
195  level_constraints[l].reinit(relevant_dofs);
196 
197  // Loop through relevant cells and faces finding those which are periodic neighbors.
199  cell = dof.begin(l),
200  endc = dof.end(l);
201  for (; cell!=endc; ++cell)
202  if (cell->level_subdomain_id() != numbers::artificial_subdomain_id)
203  {
204  for (unsigned int f = 0;
205  f < GeometryInfo<dim>::faces_per_cell;
206  ++f)
207  if (cell->has_periodic_neighbor(f))
208  {
209  if (cell->is_locally_owned_on_level())
210  {
211  Assert(cell->periodic_neighbor(f)->level_subdomain_id() != numbers::artificial_subdomain_id,
212  ExcMessage ("Periodic neighbor of a locally owned cell must either be owned or ghost."));
213  }
214  // Cell is a level-ghost and its neighbor is a level-artificial cell
215  // nothing to do here
216  else if (cell->periodic_neighbor(f)->level_subdomain_id() == numbers::artificial_subdomain_id)
217  {
218  Assert (cell->is_locally_owned_on_level() == false, ExcInternalError());
219  continue;
220  }
221 
222  const unsigned int dofs_per_face = cell->face(f)->get_fe(0).dofs_per_face;
223  std::vector<types::global_dof_index> dofs_1(dofs_per_face);
224  std::vector<types::global_dof_index> dofs_2(dofs_per_face);
225 
226  cell->periodic_neighbor(f)->face(cell->periodic_neighbor_face_no(f))->get_mg_dof_indices(l, dofs_1, 0);
227  cell->face(f)->get_mg_dof_indices(l, dofs_2, 0);
228  // Store periodicity information in the level constraint matrix
229  // Skip DoFs for which we've previously entered periodicity constraints already;
230  // this can happen, for example, for a vertex dof at a periodic boundary that we
231  // visit from more than one cell
232  for (unsigned int i=0; i<dofs_per_face; ++i)
233  if (level_constraints[l].can_store_line(dofs_2[i]) &&
234  level_constraints[l].can_store_line(dofs_1[i]) &&
235  !level_constraints[l].is_constrained(dofs_2[i]) &&
236  !level_constraints[l].is_constrained(dofs_1[i]))
237  {
238  level_constraints[l].add_line(dofs_2[i]);
239  level_constraints[l].add_entry(dofs_2[i], dofs_1[i], 1.);
240  }
241  }
242  }
243  level_constraints[l].close();
244 
245  // Initialize with empty IndexSet of correct size
247  }
248 
250 }
251 
252 
253 template <int dim, int spacedim>
254 inline
255 void
257  const typename FunctionMap<dim>::type &function_map,
258  const ComponentMask &component_mask)
259 {
260  initialize (dof);
261 
262  // allocate an IndexSet for each global level. Contents will be
263  // overwritten inside make_boundary_list.
264  const unsigned int n_levels = dof.get_triangulation().n_global_levels();
265  //At this point boundary_indices is empty.
266  boundary_indices.resize(n_levels);
267 
269  function_map,
271  component_mask);
272 }
273 
274 
275 template <int dim, int spacedim>
276 inline
277 void
279  const std::set<types::boundary_id> &boundary_ids,
280  const ComponentMask &component_mask)
281 {
282  // allocate an IndexSet for each global level. Contents will be
283  // overwritten inside make_boundary_list.
284  const unsigned int n_levels = dof.get_triangulation().n_global_levels();
285  Assert(boundary_indices.size() == 0 || boundary_indices.size() == n_levels,
286  ExcInternalError());
287  boundary_indices.resize(n_levels);
288 
290  boundary_ids,
292  component_mask);
293 }
294 
295 
296 inline
297 void
299 {
300  boundary_indices.clear();
301  refinement_edge_indices.clear();
302 }
303 
304 
305 inline
306 bool
307 MGConstrainedDoFs::is_boundary_index (const unsigned int level,
308  const types::global_dof_index index) const
309 {
310  if (boundary_indices.size() == 0)
311  return false;
312 
313  AssertIndexRange(level, boundary_indices.size());
314  return boundary_indices[level].is_element(index);
315 }
316 
317 inline
318 bool
319 MGConstrainedDoFs::at_refinement_edge (const unsigned int level,
320  const types::global_dof_index index) const
321 {
323 
324  return refinement_edge_indices[level].is_element(index);
325 }
326 
327 inline
328 bool
330  const types::global_dof_index i,
331  const types::global_dof_index j) const
332 {
333  const IndexSet &interface_dofs_on_level
334  = this->get_refinement_edge_indices(level);
335 
336  return interface_dofs_on_level.is_element(i) // at_refinement_edge(i)
337  &&
338  !interface_dofs_on_level.is_element(j) // !at_refinement_edge(j)
339  &&
340  !this->is_boundary_index(level, i) // !on_boundary(i)
341  &&
342  !this->is_boundary_index(level, j); // !on_boundary(j)
343 }
344 
345 
346 
347 
348 inline
349 const IndexSet &
350 MGConstrainedDoFs::get_boundary_indices (const unsigned int level) const
351 {
352  AssertIndexRange(level, boundary_indices.size());
353  return boundary_indices[level];
354 }
355 
356 
357 
358 inline
359 const IndexSet &
361 {
363  return refinement_edge_indices[level];
364 }
365 
366 
367 
368 
369 inline
370 bool
372 {
373  return boundary_indices.size()!=0;
374 }
375 
376 
377 
378 
379 inline
380 const ConstraintMatrix &
381 MGConstrainedDoFs::get_level_constraint_matrix (const unsigned int level) const
382 {
383  AssertIndexRange(level, level_constraints.size());
384  return level_constraints[level];
385 }
386 
387 
388 
389 DEAL_II_NAMESPACE_CLOSE
390 
391 #endif
cell_iterator begin(const unsigned int level=0) const
Definition: dof_handler.cc:723
bool is_boundary_index(const unsigned int level, const types::global_dof_index index) const
cell_iterator end() const
Definition: dof_handler.cc:751
#define AssertIndexRange(index, range)
Definition: exceptions.h:1284
std::vector< IndexSet > boundary_indices
void extract_locally_relevant_level_dofs(const DoFHandlerType &dof_handler, const unsigned int level, IndexSet &dof_set)
Definition: dof_tools.cc:1127
std::vector< ConstraintMatrix > level_constraints
bool is_interface_matrix_entry(const unsigned int level, const types::global_dof_index i, const types::global_dof_index j) const
void make_boundary_list(const DoFHandler< dim, spacedim > &mg_dof, const typename FunctionMap< dim >::type &function_map, std::vector< std::set< types::global_dof_index > > &boundary_indices, const ComponentMask &component_mask=ComponentMask())
Definition: mg_tools.cc:1208
static ::ExceptionBase & ExcMessage(std::string arg1)
const IndexSet & get_boundary_indices(const unsigned int level) const
ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:257
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:1142
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
const ConstraintMatrix & get_level_constraint_matrix(const unsigned int level) const
types::global_dof_index n_dofs() const
std::map< types::boundary_id, const Function< dim, Number > * > type
Definition: function_map.h:81
const IndexSet & get_refinement_edge_indices(unsigned int level) const
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:264
bool at_refinement_edge(const unsigned int level, const types::global_dof_index index) const
void extract_inner_interface_dofs(const DoFHandler< dim, spacedim > &mg_dof_handler, std::vector< IndexSet > &interface_dofs)
Definition: mg_tools.cc:1466
const Triangulation< dim, spacedim > & get_triangulation() const
bool is_element(const size_type index) const
Definition: index_set.h:1645
bool have_boundary_indices() const
static ::ExceptionBase & ExcInternalError()
void initialize(const DoFHandler< dim, spacedim > &dof)