Reference documentation for deal.II version GIT 6a72d26406 2023-06-07 13:05: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\}}\)
dof_accessor.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 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 
19 
20 #include <deal.II/fe/fe.h>
21 
23 #include <deal.II/grid/tria_iterator.templates.h>
24 
25 #include <vector>
26 
28 
29 /* --------------------- Static variables: DoFAccessor --------------------- */
30 
31 template <int structdim, int dim, int spacedim, bool level_dof_access>
32 const unsigned int
34 
35 template <int structdim, int dim, int spacedim, bool level_dof_access>
36 const unsigned int
38 
39 
40 
41 /* ---------------------- Functions: DoFInvalidAccessor -------------------- */
42 
43 template <int structdim, int dim, int spacedim>
45  const void *,
46  const int,
47  const int,
48  const AccessorData *)
49 {
50  Assert(false,
51  ExcMessage("You are attempting an illegal conversion between "
52  "iterator/accessor types. The constructor you call "
53  "only exists to make certain template constructs "
54  "easier to write as dimension independent code but "
55  "the conversion is not valid in the current context."));
56 }
57 
58 
59 
60 template <int structdim, int dim, int spacedim>
62  const DoFInvalidAccessor &i)
63  : InvalidAccessor<structdim, dim, spacedim>(
64  static_cast<const InvalidAccessor<structdim, dim, spacedim> &>(i))
65 {
66  Assert(false,
67  ExcMessage("You are attempting an illegal conversion between "
68  "iterator/accessor types. The constructor you call "
69  "only exists to make certain template constructs "
70  "easier to write as dimension independent code but "
71  "the conversion is not valid in the current context."));
72 }
73 
74 
75 
76 template <int structdim, int dim, int spacedim>
79  const unsigned int,
80  const types::fe_index) const
81 {
82  Assert(false, ExcInternalError());
83  return 0;
84 }
85 
86 
87 
88 template <int structdim, int dim, int spacedim>
89 void
91  const unsigned int,
93  const types::fe_index) const
94 {
95  Assert(false, ExcInternalError());
96 }
97 
98 
99 
100 /*------------------------- Functions: DoFCellAccessor -----------------------*/
101 
102 
103 
104 template <int dim, int spacedim, bool lda>
105 void
107  const std::vector<types::global_dof_index> &local_dof_indices)
108 {
109  Assert(static_cast<unsigned int>(this->present_level) <
110  this->dof_handler->object_dof_indices.size(),
111  ExcMessage("DoFHandler not initialized"));
112 
113  Assert(this->dof_handler != nullptr, typename BaseClass::ExcInvalidObject());
114 
115  internal::DoFAccessorImplementation::Implementation::
116  template set_dof_indices<dim, spacedim, lda, dim>(*this,
117  local_dof_indices,
118  this->active_fe_index());
119 }
120 
121 
122 
123 template <int dim, int spacedim, bool lda>
126  const unsigned int face,
127  const unsigned int subface) const
128 {
132  this->dof_handler);
133 }
134 
135 
136 
137 template <int dim, int spacedim, bool lda>
140  const unsigned int face,
141  const unsigned int subface) const
142 {
145  subface);
147  this->dof_handler);
148 }
149 
150 
151 
152 template <int dim, int spacedim, bool lda>
155  const unsigned int face) const
156 {
160  this->dof_handler);
161 }
162 
163 
164 
165 template <int dim, int spacedim, bool lda>
168  const unsigned int face) const
169 {
173  this->dof_handler);
174 }
175 
176 
177 
178 template <int structdim, int dim, int spacedim, bool level_dof_access>
179 inline void
181  std::vector<types::global_dof_index> &dof_indices,
182  const types::fe_index fe_index_) const
183 {
184  Assert(this->dof_handler != nullptr, ExcInvalidObject());
185 
186  const auto fe_index =
187  internal::DoFAccessorImplementation::get_fe_index_or_default(*this,
188  fe_index_);
189 
190  Assert(static_cast<unsigned int>(this->level()) <
191  this->dof_handler->object_dof_indices.size(),
192  ExcMessage(
193  "The DoFHandler to which this accessor points has not "
194  "been initialized, i.e., it doesn't appear that DoF indices "
195  "have been distributed on it."));
196 
197  // this function really only makes sense if either a) there are degrees of
198  // freedom defined on the present object, or b) the object is non-active
199  // objects but all degrees of freedom are located on vertices, since
200  // otherwise there are degrees of freedom on sub-objects which are not
201  // allocated for this non-active thing
202  Assert(this->fe_index_is_active(fe_index) ||
203  (this->dof_handler->get_fe(fe_index).n_dofs_per_cell() ==
204  this->n_vertices() *
205  this->dof_handler->get_fe(fe_index).n_dofs_per_vertex()),
206  ExcInternalError());
207 
208  // now do the actual work
209  ::internal::DoFAccessorImplementation::Implementation::get_dof_indices(
210  *this, dof_indices, fe_index);
211 }
212 
213 
214 
215 template <int structdim, int dim, int spacedim, bool level_dof_access>
216 inline void
218  const int level,
219  std::vector<types::global_dof_index> &dof_indices,
220  const types::fe_index fe_index_) const
221 {
222  Assert(this->dof_handler != nullptr, ExcInvalidObject());
223  Assert(this->dof_handler->mg_vertex_dofs.size() > 0,
224  ExcMessage("Multigrid DoF indices can only be accessed after "
225  "DoFHandler::distribute_mg_dofs() has been called!"));
226 
227  const auto fe_index =
228  internal::DoFAccessorImplementation::get_fe_index_or_default(*this,
229  fe_index_);
230 
231  internal::DoFAccessorImplementation::Implementation::get_mg_dof_indices(
232  *this, level, dof_indices, fe_index);
233 }
234 
235 
236 
237 template <int structdim, int dim, int spacedim, bool level_dof_access>
238 inline void
240  const int level,
241  const std::vector<types::global_dof_index> &dof_indices,
242  const types::fe_index fe_index_)
243 {
244  Assert(this->dof_handler != nullptr, ExcInvalidObject());
245 
246  const auto fe_index =
247  internal::DoFAccessorImplementation::get_fe_index_or_default(*this,
248  fe_index_);
249 
250  internal::DoFAccessorImplementation::Implementation::set_mg_dof_indices(
251  *this, level, dof_indices, fe_index);
252 }
253 
254 
255 
256 namespace internal
257 {
258  namespace DoFAccessorImplementation
259  {
260  template <int dim, int spacedim, bool level_dof_access>
261  void
263  const ::DoFCellAccessor<dim, spacedim, level_dof_access> &accessor,
264  boost::container::small_vector<types::global_dof_index, 27> &dof_indices,
265  const unsigned int fe_index)
266  {
267  Implementation::process_dof_indices(
268  accessor,
269  dof_indices,
270  fe_index,
271  Implementation::DoFIndexProcessor<dim, spacedim>(),
272  [](auto stored_index, auto dof_ptr) { *dof_ptr = stored_index; },
273  false);
274  }
275  } // namespace DoFAccessorImplementation
276 } // namespace internal
277 
278 
279 
280 template <int dimension_, int space_dimension_, bool level_dof_access>
281 inline void
283  get_dof_indices(std::vector<types::global_dof_index> &dof_indices) const
284 {
285  Assert(this->is_active(),
286  ExcMessage("get_dof_indices() only works on active cells."));
287  Assert(this->is_artificial() == false,
288  ExcMessage("Can't ask for DoF indices on artificial cells."));
289  AssertDimension(dof_indices.size(), this->get_fe().n_dofs_per_cell());
290 
291  ::internal::DoFAccessorImplementation::Implementation::get_dof_indices(
292  *this, dof_indices, this->active_fe_index());
293 }
294 
295 
296 
297 template <int dimension_, int space_dimension_, bool level_dof_access>
298 inline void
300  get_mg_dof_indices(std::vector<types::global_dof_index> &dof_indices) const
301 {
302  Assert(this->dof_handler->mg_vertex_dofs.size() > 0,
303  ExcMessage("Multigrid DoF indices can only be accessed after "
304  "DoFHandler::distribute_mg_dofs() has been called!"));
306  get_mg_dof_indices(this->level(), dof_indices);
307 }
308 
309 
310 
311 template <int dimension_, int space_dimension_, bool level_dof_access>
312 inline void
314  set_mg_dof_indices(const std::vector<types::global_dof_index> &dof_indices)
315 {
316  Assert(this->dof_handler->mg_vertex_dofs.size() > 0,
317  ExcMessage("Multigrid DoF indices can only be accessed after "
318  "DoFHandler::distribute_mg_dofs() has been called!"));
320  set_mg_dof_indices(this->level(), dof_indices);
321 }
322 
323 // --------------------------------------------------------------------------
324 // explicit instantiations
325 #include "dof_accessor.inst"
326 
TriaIterator< CellAccessor< dim, spacedim > > periodic_neighbor(const unsigned int i) const
TriaIterator< CellAccessor< dim, spacedim > > neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
TriaIterator< CellAccessor< dim, spacedim > > periodic_neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
TriaIterator< CellAccessor< dim, spacedim > > neighbor_or_periodic_neighbor(const unsigned int i) const
void set_mg_dof_indices(const int level, const std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index)
void get_mg_dof_indices(const int level, std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index) const
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index) const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > periodic_neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > neighbor_or_periodic_neighbor(const unsigned int i) const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > periodic_neighbor(const unsigned int i) const
void set_dof_indices(const std::vector< types::global_dof_index > &dof_indices)
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
void set_mg_dof_indices(const std::vector< types::global_dof_index > &dof_indices)
void get_mg_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
DoFInvalidAccessor(const void *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
Definition: dof_accessor.cc:44
types::global_dof_index dof_index(const unsigned int i, const types::fe_index fe_index=DoFHandler< dim, spacedim >::default_fe_index) const
Definition: dof_accessor.cc:78
typename InvalidAccessor< structdim, dim, spacedim >::AccessorData AccessorData
void set_dof_index(const unsigned int i, const types::global_dof_index index, const types::fe_index fe_index=numbers::invalid_fe_index) const
Definition: dof_accessor.cc:90
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:474
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:475
unsigned int level
Definition: grid_out.cc:4618
static ::ExceptionBase & ExcInternalError()
#define Assert(cond, exc)
Definition: exceptions.h:1614
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1787
static ::ExceptionBase & ExcMessage(std::string arg1)
void get_cell_dof_indices(const ::DoFCellAccessor< dim, spacedim, level_dof_access > &accessor, boost::container::small_vector< types::global_dof_index, 27 > &dof_indices, const unsigned int fe_index)
unsigned int global_dof_index
Definition: types.h:82
unsigned short int fe_index
Definition: types.h:60