Reference documentation for deal.II version GIT b065060f03 2023-05-30 23:50:02+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.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2022 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_dof_accessor_h
17 #define dealii_dof_accessor_h
18 
19 
20 #include <deal.II/base/config.h>
21 
24 
26 
27 #include <boost/container/small_vector.hpp>
28 
29 #include <set>
30 #include <vector>
31 
33 
34 // Forward declarations
35 #ifndef DOXYGEN
36 template <typename number>
37 class FullMatrix;
38 template <typename number>
39 class Vector;
40 template <typename number>
41 class AffineConstraints;
42 
43 template <typename Accessor>
44 class TriaRawIterator;
45 
46 template <int, int>
47 class FiniteElement;
48 
49 namespace internal
50 {
51  namespace DoFCellAccessorImplementation
52  {
53  struct Implementation;
54  }
55 
56  namespace DoFHandlerImplementation
57  {
58  struct Implementation;
59  namespace Policy
60  {
61  struct Implementation;
62  }
63  } // namespace DoFHandlerImplementation
64 
65  namespace hp
66  {
67  namespace DoFHandlerImplementation
68  {
69  struct Implementation;
70  }
71  } // namespace hp
72 } // namespace internal
73 #endif
74 
75 // note: the file dof_accessor.templates.h is included at the end of
76 // this file. this includes a lot of templates and thus makes
77 // compilation slower, but at the same time allows for more aggressive
78 // inlining and thus faster code.
79 
80 
81 namespace internal
82 {
83  namespace DoFAccessorImplementation
84  {
101  template <int structdim, int dim, int spacedim>
102  struct Inheritance
103  {
109  };
110 
111 
117  template <int dim, int spacedim>
118  struct Inheritance<dim, dim, spacedim>
119  {
125  };
126 
127  struct Implementation;
128  } // namespace DoFAccessorImplementation
129 } // namespace internal
130 
131 
132 /* -------------------------------------------------------------------------- */
133 
134 
135 
208 template <int structdim, int dim, int spacedim, bool level_dof_access>
209 class DoFAccessor : public ::internal::DoFAccessorImplementation::
210  Inheritance<structdim, dim, spacedim>::BaseClass
211 {
212 public:
217  static constexpr unsigned int dimension = dim;
218 
223  static constexpr unsigned int space_dimension = spacedim;
224 
229  using BaseClass = typename ::internal::DoFAccessorImplementation::
230  Inheritance<structdim, dimension, space_dimension>::BaseClass;
231 
236 
248 
265  const int level,
266  const int index,
268 
273  default;
274 
278  DoFAccessor( // NOLINT
280  default; // NOLINT
281 
285  ~DoFAccessor() = default;
286 
299  template <int structdim2, int dim2, int spacedim2>
301 
306  template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
309 
313  template <bool level_dof_access2>
315 
327  delete;
328 
333  operator=( // NOLINT
335  default; // NOLINT
336 
346 
350  template <bool level_dof_access2>
351  void
353 
358  void
360 
365  static bool
367 
379  child(const unsigned int c) const;
380 
386  typename ::internal::DoFHandlerImplementation::
387  Iterators<dim, spacedim, level_dof_access>::line_iterator
388  line(const unsigned int i) const;
389 
395  typename ::internal::DoFHandlerImplementation::
396  Iterators<dim, spacedim, level_dof_access>::quad_iterator
397  quad(const unsigned int i) const;
398 
444  void
446  std::vector<types::global_dof_index> &dof_indices,
448 
455  void
457  const int level,
458  std::vector<types::global_dof_index> &dof_indices,
460 
464  void
466  const int level,
467  const std::vector<types::global_dof_index> &dof_indices,
469 
493  const unsigned int vertex,
494  const unsigned int i,
496 
504  const int level,
505  const unsigned int vertex,
506  const unsigned int i,
508 
537  dof_index(const unsigned int i,
539 
544  mg_dof_index(const int level, const unsigned int i) const;
545 
567  unsigned int
569 
578  nth_active_fe_index(const unsigned int n) const;
579 
586  std::set<types::fe_index>
588 
598  bool
600 
608 
619  "This accessor object has not been "
620  "associated with any DoFHandler object.");
653 
654 protected:
659 
660 public:
674  template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
675  bool
678 
682  template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
683  bool
686 
687 protected:
691  void
693 
712  void
714  const unsigned int i,
715  const types::global_dof_index index,
717 
718  void
720  const unsigned int i,
721  const types::global_dof_index index) const;
722 
723  void
725  const int level,
726  const unsigned int vertex,
727  const unsigned int i,
728  const types::global_dof_index index,
730 
731  // Iterator classes need to be friends because they need to access
732  // operator== and operator!=.
733  template <typename>
734  friend class TriaRawIterator;
735  template <int, int, int, bool>
736  friend class DoFAccessor;
737 
738 private:
739  // Make the DoFHandler class a friend so that it can call the set_xxx()
740  // functions.
741  friend class DoFHandler<dim, spacedim>;
742 
743  friend struct ::internal::DoFHandlerImplementation::Policy::
744  Implementation;
745  friend struct ::internal::DoFHandlerImplementation::Implementation;
746  friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
747  friend struct ::internal::DoFCellAccessorImplementation::Implementation;
748  friend struct ::internal::DoFAccessorImplementation::Implementation;
749 };
750 
751 
752 
760 template <int spacedim, bool level_dof_access>
761 class DoFAccessor<0, 1, spacedim, level_dof_access>
762  : public TriaAccessor<0, 1, spacedim>
763 {
764 public:
769  static constexpr unsigned int dimension = 1;
770 
775  static constexpr unsigned int space_dimension = spacedim;
776 
782 
787 
799 
818  const typename TriaAccessor<0, 1, spacedim>::VertexKind vertex_kind,
819  const unsigned int vertex_index,
821 
829  const int = 0,
830  const int = 0,
832 
845  template <int structdim2, int dim2, int spacedim2>
847 
852  template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
855 
860 
864  // NOLINTNEXTLINE OSX does not compile with noexcept
866 
870  ~DoFAccessor() = default;
871 
883 
889  default;
890 
898  const DoFHandler<1, spacedim> &
900 
904  template <bool level_dof_access2>
905  void
906  copy_from(const DoFAccessor<0, 1, spacedim, level_dof_access2> &a);
907 
912  void
913  copy_from(const TriaAccessorBase<0, 1, spacedim> &da);
914 
927  TriaIterator<DoFAccessor<0, 1, spacedim, level_dof_access>>
928  child(const unsigned int c) const;
929 
936  typename ::internal::DoFHandlerImplementation::
937  Iterators<1, spacedim, level_dof_access>::line_iterator
938  line(const unsigned int i) const;
939 
946  typename ::internal::DoFHandlerImplementation::
947  Iterators<1, spacedim, level_dof_access>::quad_iterator
948  quad(const unsigned int i) const;
949 
993  void
995  std::vector<types::global_dof_index> &dof_indices,
997 
1004  void
1006  const int level,
1007  std::vector<types::global_dof_index> &dof_indices,
1008  const types::fe_index fe_index = numbers::invalid_fe_index) const;
1009 
1030  const unsigned int vertex,
1031  const unsigned int i,
1032  const types::fe_index fe_index = numbers::invalid_fe_index) const;
1033 
1052  dof_index(const unsigned int i,
1053  const types::fe_index fe_index = numbers::invalid_fe_index) const;
1054 
1073  unsigned int
1075 
1083  types::fe_index
1084  nth_active_fe_index(const unsigned int n) const;
1085 
1094  bool
1096 
1102  const FiniteElement<1, spacedim> &
1104 
1115  "This accessor object has not been "
1116  "associated with any DoFHandler object.");
1149 
1150 protected:
1154  DoFHandler<1, spacedim> *dof_handler;
1155 
1159  template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1160  bool
1161  operator==(
1162  const DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1163 
1167  template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1168  bool
1169  operator!=(
1170  const DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1171 
1175  void
1176  set_dof_handler(DoFHandler<1, spacedim> *dh);
1177 
1196  void
1198  const unsigned int i,
1199  const types::global_dof_index index,
1200  const types::fe_index fe_index = numbers::invalid_fe_index) const;
1201 
1202  // Iterator classes need to be friends because they need to access
1203  // operator== and operator!=.
1204  template <typename>
1205  friend class TriaRawIterator;
1206 
1207 
1208  // Make the DoFHandler class a friend so that it can call the set_xxx()
1209  // functions.
1210  friend class DoFHandler<1, spacedim>;
1211 
1212  friend struct ::internal::DoFHandlerImplementation::Policy::
1213  Implementation;
1214  friend struct ::internal::DoFHandlerImplementation::Implementation;
1215  friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
1216  friend struct ::internal::DoFCellAccessorImplementation::Implementation;
1217 };
1218 
1219 
1220 
1221 /* -------------------------------------------------------------------------- */
1222 
1223 
1244 template <int structdim, int dim, int spacedim = dim>
1245 class DoFInvalidAccessor : public InvalidAccessor<structdim, dim, spacedim>
1246 {
1247 public:
1253 
1261  DoFInvalidAccessor(const void * parent = nullptr,
1262  const int level = -1,
1263  const int index = -1,
1264  const AccessorData *local_data = nullptr);
1265 
1274 
1279  template <typename OtherAccessor>
1280  DoFInvalidAccessor(const OtherAccessor &);
1281 
1288  dof_index(const unsigned int i,
1289  const types::fe_index fe_index =
1291 
1297  void
1298  set_dof_index(
1299  const unsigned int i,
1300  const types::global_dof_index index,
1302 };
1303 
1304 
1305 
1306 /* -------------------------------------------------------------------------- */
1307 
1308 
1320 template <int dimension_, int space_dimension_, bool level_dof_access>
1321 class DoFCellAccessor : public DoFAccessor<dimension_,
1322  dimension_,
1323  space_dimension_,
1324  level_dof_access>
1325 {
1326 public:
1330  static const unsigned int dim = dimension_;
1331 
1335  static const unsigned int spacedim = space_dimension_;
1336 
1337 
1342 
1347  using BaseClass =
1349 
1354 
1359  using face_iterator = TriaIterator<DoFAccessor<dimension_ - 1,
1360  dimension_,
1361  space_dimension_,
1362  level_dof_access>>;
1363 
1375  const int level,
1376  const int index,
1377  const AccessorData *local_data);
1378 
1391  template <int structdim2, int dim2, int spacedim2>
1393 
1398  template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1401 
1407  default;
1408 
1412  DoFCellAccessor( // NOLINT
1414  &&) = default; // NOLINT
1415 
1419  ~DoFCellAccessor() = default;
1420 
1433  delete;
1434 
1439  operator=( // NOLINT
1441  &&) = default; // NOLINT
1442 
1457  parent() const;
1458 
1472  neighbor(const unsigned int i) const;
1473 
1480  periodic_neighbor(const unsigned int i) const;
1481 
1488  neighbor_or_periodic_neighbor(const unsigned int i) const;
1489 
1496  child(const unsigned int i) const;
1497 
1501  boost::container::small_vector<
1502  TriaIterator<
1506 
1514  face(const unsigned int i) const;
1515 
1519  boost::container::small_vector<face_iterator,
1522 
1530  neighbor_child_on_subface(const unsigned int face_no,
1531  const unsigned int subface_no) const;
1532 
1540  periodic_neighbor_child_on_subface(const unsigned int face_no,
1541  const unsigned int subface_no) const;
1542 
1571  template <class InputVector, typename number>
1572  void
1573  get_dof_values(const InputVector &values, Vector<number> &local_values) const;
1574 
1592  template <class InputVector, typename ForwardIterator>
1593  void
1594  get_dof_values(const InputVector &values,
1595  ForwardIterator local_values_begin,
1596  ForwardIterator local_values_end) const;
1597 
1618  template <class InputVector, typename ForwardIterator>
1619  void
1622  const InputVector & values,
1623  ForwardIterator local_values_begin,
1624  ForwardIterator local_values_end) const;
1625 
1650  template <class OutputVector, typename number>
1651  void
1652  set_dof_values(const Vector<number> &local_values,
1653  OutputVector & values) const;
1654 
1686  template <class InputVector, typename number>
1687  void
1688  get_interpolated_dof_values(
1689  const InputVector & values,
1690  Vector<number> & interpolated_values,
1692 
1754  template <class OutputVector, typename number>
1755  void
1756  set_dof_values_by_interpolation(
1757  const Vector<number> &local_values,
1758  OutputVector & values,
1760  const bool perform_check = false) const;
1761 
1771  template <class OutputVector, typename number>
1772  void
1773  distribute_local_to_global_by_interpolation(
1774  const Vector<number> &local_values,
1775  OutputVector & values,
1777 
1792  template <typename number, typename OutputVector>
1793  void
1795  OutputVector & global_destination) const;
1796 
1811  template <typename ForwardIterator, typename OutputVector>
1812  void
1813  distribute_local_to_global(ForwardIterator local_source_begin,
1814  ForwardIterator local_source_end,
1815  OutputVector & global_destination) const;
1816 
1830  template <typename ForwardIterator, typename OutputVector>
1831  void
1834  ForwardIterator local_source_begin,
1835  ForwardIterator local_source_end,
1836  OutputVector & global_destination) const;
1837 
1844  template <typename number, typename OutputMatrix>
1845  void
1847  OutputMatrix &global_destination) const;
1848 
1853  template <typename number, typename OutputMatrix, typename OutputVector>
1854  void
1856  const Vector<number> & local_vector,
1857  OutputMatrix & global_matrix,
1858  OutputVector & global_vector) const;
1859 
1884  void
1886  std::vector<types::global_dof_index> &dof_indices) const;
1887 
1919  void
1920  get_dof_indices(std::vector<types::global_dof_index> &dof_indices) const;
1921 
1926  void
1927  get_mg_dof_indices(std::vector<types::global_dof_index> &dof_indices) const;
1928 
1954  get_fe() const;
1955 
1983 
2010  void
2020  void
2021  set_dof_indices(const std::vector<types::global_dof_index> &dof_indices);
2022 
2026  void
2027  set_mg_dof_indices(const std::vector<types::global_dof_index> &dof_indices);
2028 
2055  get_future_fe() const;
2056 
2077 
2086  void
2088 
2095  bool
2097 
2106  void
2112 private:
2113  friend struct ::internal::DoFCellAccessorImplementation::Implementation;
2114 };
2115 
2116 
2117 template <int structdim, int dim, int spacedim, bool level_dof_access>
2118 inline bool
2120 {
2121  return level_dof_access;
2122 }
2123 
2124 
2125 
2126 template <int structdim, int dim, int spacedim>
2127 template <typename OtherAccessor>
2129  const OtherAccessor &)
2130 {
2131  Assert(false,
2132  ExcMessage("You are attempting an illegal conversion between "
2133  "iterator/accessor types. The constructor you call "
2134  "only exists to make certain template constructs "
2135  "easier to write as dimension independent code but "
2136  "the conversion is not valid in the current context."));
2137 }
2138 
2139 
2140 
2142 
2143 // include more templates
2144 #include <deal.II/dofs/dof_accessor.templates.h>
2145 
2146 
2147 #endif
DoFAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
DoFAccessor(const Triangulation< 1, spacedim > *, const int=0, const int=0, const DoFHandler< 1, spacedim > *dof_handler=0)
DoFAccessor(DoFAccessor< 0, 1, spacedim, level_dof_access > &&)=default
DoFAccessor(const Triangulation< 1, spacedim > *tria, const typename TriaAccessor< 0, 1, spacedim >::VertexKind vertex_kind, const unsigned int vertex_index, const DoFHandler< 1, spacedim > *dof_handler)
DoFAccessor< 0, 1, spacedim, level_dof_access > & operator=(DoFAccessor< 0, 1, spacedim, level_dof_access > &&) noexcept=default
DoFAccessor(const DoFAccessor< 0, 1, spacedim, level_dof_access > &)=default
DoFAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
DoFAccessor< 0, 1, spacedim, level_dof_access > & operator=(const DoFAccessor< 0, 1, spacedim, level_dof_access > &da)=delete
static constexpr unsigned int space_dimension
Definition: dof_accessor.h:223
void set_mg_dof_index(const int level, const unsigned int i, const types::global_dof_index index) const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index fe_index) const
DoFAccessor< structdim, dim, spacedim, level_dof_access > & operator=(DoFAccessor< structdim, dim, spacedim, level_dof_access > &&)
types::fe_index nth_active_fe_index(const unsigned int n) const
DoFAccessor(DoFAccessor< structdim, dim, spacedim, level_dof_access > &&)
DoFAccessor(const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &)
const DoFHandler< dim, spacedim > & get_dof_handler() const
void set_mg_vertex_dof_index(const int level, const unsigned int vertex, const unsigned int i, const types::global_dof_index index, const types::fe_index fe_index=numbers::invalid_fe_index) const
~DoFAccessor()=default
static constexpr unsigned int dimension
Definition: dof_accessor.h:217
void set_dof_handler(DoFHandler< dim, spacedim > *dh)
DoFAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
typename ::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::quad_iterator quad(const unsigned int i) const
std::set< types::fe_index > get_active_fe_indices() const
DoFAccessor(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &)=default
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
types::global_dof_index mg_vertex_dof_index(const int level, const unsigned int vertex, const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
TriaIterator< DoFAccessor< structdim, dim, spacedim, level_dof_access > > child(const unsigned int c) const
DoFAccessor(const Triangulation< dim, spacedim > *tria, const int level, const int index, const DoFHandler< dim, spacedim > *dof_handler)
bool operator!=(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &) const
void copy_from(const TriaAccessorBase< structdim, dim, spacedim > &da)
types::global_dof_index dof_index(const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
void copy_from(const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &a)
typename ::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::line_iterator line(const unsigned int i) const
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index) const
bool operator==(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &) const
types::global_dof_index mg_dof_index(const int level, const unsigned int i) const
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
unsigned int n_active_fe_indices() const
DoFAccessor< structdim, dim, spacedim, level_dof_access > & operator=(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &da)=delete
static bool is_level_cell()
typename ::internal::DoFAccessorImplementation::Inheritance< structdim, dimension, space_dimension >::BaseClass BaseClass
Definition: dof_accessor.h:230
types::global_dof_index vertex_dof_index(const unsigned int vertex, const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
DoFHandler< dim, spacedim > * dof_handler
Definition: dof_accessor.h:658
DoFAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
bool fe_index_is_active(const types::fe_index fe_index) const
boost::container::small_vector< TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > >, GeometryInfo< dimension_ >::max_children_per_cell > child_iterators() const
DoFCellAccessor< dimension_, space_dimension_, level_dof_access > & operator=(DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &&)=default
DoFCellAccessor(DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &&)=default
void get_active_or_mg_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
void get_dof_values(const InputVector &values, Vector< number > &local_values) const
face_iterator face(const unsigned int i) const
types::fe_index future_fe_index() const
DoFCellAccessor(const Triangulation< dimension_, space_dimension_ > *tria, const int level, const int index, const AccessorData *local_data)
void get_dof_values(const InputVector &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) const
void set_future_fe_index(const types::fe_index i) const
void distribute_local_to_global(const Vector< number > &local_source, OutputVector &global_destination) const
void get_dof_values(const AffineConstraints< typename InputVector::value_type > &constraints, const InputVector &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) const
void distribute_local_to_global(const FullMatrix< number > &local_source, OutputMatrix &global_destination) const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > child(const unsigned int i) const
boost::container::small_vector< face_iterator, GeometryInfo< dimension_ >::faces_per_cell > face_iterators() const
void set_active_fe_index(const types::fe_index i) const
const FiniteElement< dimension_, space_dimension_ > & get_future_fe() const
DoFCellAccessor< dimension_, space_dimension_, level_dof_access > & operator=(const DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &da)=delete
DoFCellAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > neighbor(const unsigned int i) const
DoFCellAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
void clear_future_fe_index() const
void distribute_local_to_global(const FullMatrix< number > &local_matrix, const Vector< number > &local_vector, OutputMatrix &global_matrix, OutputVector &global_vector) const
void distribute_local_to_global(const AffineConstraints< typename OutputVector::value_type > &constraints, ForwardIterator local_source_begin, ForwardIterator local_source_end, OutputVector &global_destination) const
void set_dof_values(const Vector< number > &local_values, OutputVector &values) const
DoFCellAccessor(const DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &)=default
const FiniteElement< dimension_, space_dimension_ > & get_fe() const
void distribute_local_to_global(ForwardIterator local_source_begin, ForwardIterator local_source_end, OutputVector &global_destination) const
~DoFCellAccessor()=default
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > parent() const
bool future_fe_index_set() const
types::fe_index active_fe_index() const
DoFInvalidAccessor(const void *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
Definition: dof_accessor.cc:44
typename InvalidAccessor< structdim, dim, spacedim >::AccessorData AccessorData
Definition: vector.h:109
#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
#define DeclException0(Exception0)
Definition: exceptions.h:465
static ::ExceptionBase & ExcVectorDoesNotMatch()
#define Assert(cond, exc)
Definition: exceptions.h:1614
static ::ExceptionBase & ExcNotActive()
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:488
static ::ExceptionBase & ExcInvalidObject()
static ::ExceptionBase & ExcMatrixDoesNotMatch()
static ::ExceptionBase & ExcCantCompareIterators()
static ::ExceptionBase & ExcVectorNotEmpty()
static ::ExceptionBase & ExcMessage(std::string arg1)
Definition: hp.h:118
const types::fe_index invalid_fe_index
Definition: types.h:228
Definition: types.h:33
unsigned int global_dof_index
Definition: types.h:82
unsigned short int fe_index
Definition: types.h:60
const ::Triangulation< dim, spacedim > & tria