Reference documentation for deal.II version 8.4.2
mg_constrained_dofs.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2015 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;
54  template <int dim, int spacedim>
55  void initialize (const DoFHandler<dim,spacedim> &dof);
56 
65  template <int dim, int spacedim>
66  void initialize(const DoFHandler<dim,spacedim> &dof,
67  const typename FunctionMap<dim>::type &function_map,
68  const ComponentMask &component_mask = ComponentMask());
69 
73  void clear();
74 
78  bool is_boundary_index (const unsigned int level,
79  const types::global_dof_index index) const;
80 
84  bool at_refinement_edge (const unsigned int level,
85  const types::global_dof_index index) const;
86 
93  const IndexSet &
94  get_boundary_indices (const unsigned int level) const;
95 
96 
101  const IndexSet &
102  get_refinement_edge_indices (unsigned int level) const;
103 
104 
108  bool have_boundary_indices () const;
109 
110 private:
111 
115  std::vector<IndexSet> boundary_indices;
116 
121  std::vector<IndexSet> refinement_edge_indices;
122 };
123 
124 
125 template <int dim, int spacedim>
126 inline
127 void
129 {
130  boundary_indices.clear();
131 
132  const unsigned int nlevels = dof.get_triangulation().n_global_levels();
133 
134  refinement_edge_indices.resize(nlevels);
135  for (unsigned int l=0; l<nlevels; ++l)
137 
139 }
140 
141 
142 template <int dim, int spacedim>
143 inline
144 void
146  const typename FunctionMap<dim>::type &function_map,
147  const ComponentMask &component_mask)
148 {
149  initialize (dof);
150 
151  // allocate an IndexSet for each global level. Contents will be
152  // overwritten inside make_boundary_list.
153  const unsigned int n_levels = dof.get_triangulation().n_global_levels();
154  boundary_indices.resize(n_levels);
155 
157  function_map,
159  component_mask);
160 }
161 
162 
163 inline
164 void
166 {
167  boundary_indices.clear();
168  refinement_edge_indices.clear();
169 }
170 
171 
172 inline
173 bool
174 MGConstrainedDoFs::is_boundary_index (const unsigned int level,
175  const types::global_dof_index index) const
176 {
177  if (boundary_indices.size() == 0)
178  return false;
179 
180  AssertIndexRange(level, boundary_indices.size());
181  return boundary_indices[level].is_element(index);
182 }
183 
184 inline
185 bool
186 MGConstrainedDoFs::at_refinement_edge (const unsigned int level,
187  const types::global_dof_index index) const
188 {
190 
191  return refinement_edge_indices[level].is_element(index);
192 }
193 
194 
195 
196 
197 inline
198 const IndexSet &
199 MGConstrainedDoFs::get_boundary_indices (const unsigned int level) const
200 {
201  AssertIndexRange(level, boundary_indices.size());
202  return boundary_indices[level];
203 }
204 
205 
206 
207 inline
208 const IndexSet &
210 {
212  return refinement_edge_indices[level];
213 }
214 
215 
216 
217 
218 inline
219 bool
221 {
222  return boundary_indices.size()!=0;
223 }
224 
225 
226 DEAL_II_NAMESPACE_CLOSE
227 
228 #endif
bool is_boundary_index(const unsigned int level, const types::global_dof_index index) const
#define AssertIndexRange(index, range)
Definition: exceptions.h:1081
std::vector< IndexSet > boundary_indices
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:1177
const IndexSet & get_boundary_indices(const unsigned int level) const
unsigned int global_dof_index
Definition: types.h:88
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
std::vector< IndexSet > refinement_edge_indices
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:1485
const Triangulation< dim, spacedim > & get_triangulation() const
bool have_boundary_indices() const
void initialize(const DoFHandler< dim, spacedim > &dof)