Reference documentation for deal.II version 8.5.1
mg_constrained_dofs.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2016 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 
67  template <int dim, int spacedim>
68  void initialize(const DoFHandler<dim,spacedim> &dof,
69  const typename FunctionMap<dim>::type &function_map,
70  const ComponentMask &component_mask = ComponentMask()) DEAL_II_DEPRECATED;
71 
82  template <int dim, int spacedim>
84  const std::set<types::boundary_id> &boundary_ids,
85  const ComponentMask &component_mask = ComponentMask());
86 
90  void clear();
91 
95  bool is_boundary_index (const unsigned int level,
96  const types::global_dof_index index) const;
97 
101  bool at_refinement_edge (const unsigned int level,
102  const types::global_dof_index index) const;
103 
110  const IndexSet &
111  get_boundary_indices (const unsigned int level) const;
112 
113 
118  const IndexSet &
119  get_refinement_edge_indices (unsigned int level) const;
120 
121 
125  bool have_boundary_indices () const;
126 
127 private:
128 
132  std::vector<IndexSet> boundary_indices;
133 
138  std::vector<IndexSet> refinement_edge_indices;
139 };
140 
141 
142 template <int dim, int spacedim>
143 inline
144 void
146 {
147  boundary_indices.clear();
148  refinement_edge_indices.clear();
149 
150  const unsigned int nlevels = dof.get_triangulation().n_global_levels();
151 
152  //At this point refinement_edge_indices is empty.
153  refinement_edge_indices.resize(nlevels);
154  for (unsigned int l=0; l<nlevels; ++l)
156 
158 }
159 
160 
161 template <int dim, int spacedim>
162 inline
163 void
165  const typename FunctionMap<dim>::type &function_map,
166  const ComponentMask &component_mask)
167 {
168  initialize (dof);
169 
170  // allocate an IndexSet for each global level. Contents will be
171  // overwritten inside make_boundary_list.
172  const unsigned int n_levels = dof.get_triangulation().n_global_levels();
173  //At this point boundary_indices is empty.
174  boundary_indices.resize(n_levels);
175 
177  function_map,
179  component_mask);
180 }
181 
182 
183 template <int dim, int spacedim>
184 inline
185 void
187  const std::set<types::boundary_id> &boundary_ids,
188  const ComponentMask &component_mask)
189 {
190  // allocate an IndexSet for each global level. Contents will be
191  // overwritten inside make_boundary_list.
192  const unsigned int n_levels = dof.get_triangulation().n_global_levels();
193  Assert(boundary_indices.size() == 0 || boundary_indices.size() == n_levels,
194  ExcInternalError());
195  boundary_indices.resize(n_levels);
196 
198  boundary_ids,
200  component_mask);
201 }
202 
203 
204 inline
205 void
207 {
208  boundary_indices.clear();
209  refinement_edge_indices.clear();
210 }
211 
212 
213 inline
214 bool
215 MGConstrainedDoFs::is_boundary_index (const unsigned int level,
216  const types::global_dof_index index) const
217 {
218  if (boundary_indices.size() == 0)
219  return false;
220 
221  AssertIndexRange(level, boundary_indices.size());
222  return boundary_indices[level].is_element(index);
223 }
224 
225 inline
226 bool
227 MGConstrainedDoFs::at_refinement_edge (const unsigned int level,
228  const types::global_dof_index index) const
229 {
231 
232  return refinement_edge_indices[level].is_element(index);
233 }
234 
235 
236 
237 
238 inline
239 const IndexSet &
240 MGConstrainedDoFs::get_boundary_indices (const unsigned int level) const
241 {
242  AssertIndexRange(level, boundary_indices.size());
243  return boundary_indices[level];
244 }
245 
246 
247 
248 inline
249 const IndexSet &
251 {
253  return refinement_edge_indices[level];
254 }
255 
256 
257 
258 
259 inline
260 bool
262 {
263  return boundary_indices.size()!=0;
264 }
265 
266 
267 DEAL_II_NAMESPACE_CLOSE
268 
269 #endif
bool is_boundary_index(const unsigned int level, const types::global_dof_index index) const
#define AssertIndexRange(index, range)
Definition: exceptions.h:1170
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:1173
const IndexSet & get_boundary_indices(const unsigned int level) const
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:313
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
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:1431
const Triangulation< dim, spacedim > & get_triangulation() const
bool have_boundary_indices() const
static ::ExceptionBase & ExcInternalError()
void initialize(const DoFHandler< dim, spacedim > &dof)