Reference documentation for deal.II version 9.6.0
\(\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\}}\)
Loading...
Searching...
No Matches
dof_accessor.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1998 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_dof_accessor_h
16#define dealii_dof_accessor_h
17
18
19#include <deal.II/base/config.h>
20
23
25
27
28#include <boost/container/small_vector.hpp>
29
30#include <set>
31#include <vector>
32
34
35// Forward declarations
36#ifndef DOXYGEN
37template <typename number>
38class FullMatrix;
39template <typename number>
40class Vector;
41template <typename number>
43
44template <typename Accessor>
45class TriaRawIterator;
46
47template <int, int>
48class FiniteElement;
49
50namespace internal
51{
52 namespace DoFCellAccessorImplementation
53 {
54 struct Implementation;
55 }
56
57 namespace DoFHandlerImplementation
58 {
59 struct Implementation;
60 namespace Policy
61 {
62 struct Implementation;
63 }
64 } // namespace DoFHandlerImplementation
65
66 namespace hp
67 {
68 namespace DoFHandlerImplementation
69 {
70 struct Implementation;
71 }
72 } // namespace hp
73} // namespace internal
74#endif
75
76// note: the file dof_accessor.templates.h is included at the end of
77// this file. this includes a lot of templates and thus makes
78// compilation slower, but at the same time allows for more aggressive
79// inlining and thus faster code.
80
81
82namespace internal
83{
84 namespace DoFAccessorImplementation
85 {
102 template <int structdim, int dim, int spacedim>
111
112
118 template <int dim, int spacedim>
119 struct Inheritance<dim, dim, spacedim>
120 {
126 };
127
128 struct Implementation;
129 } // namespace DoFAccessorImplementation
130} // namespace internal
131
132
133/* -------------------------------------------------------------------------- */
134
135
136
209template <int structdim, int dim, int spacedim, bool level_dof_access>
210class DoFAccessor : public ::internal::DoFAccessorImplementation::
211 Inheritance<structdim, dim, spacedim>::BaseClass
212{
213public:
218 static constexpr unsigned int dimension = dim;
219
224 static constexpr unsigned int space_dimension = spacedim;
225
230 using BaseClass = typename ::internal::DoFAccessorImplementation::
231 Inheritance<structdim, dimension, space_dimension>::BaseClass;
232
237
249
266 const int level,
267 const int index,
269
274 default;
275
279 DoFAccessor( // NOLINT
281 default; // NOLINT
282
286 ~DoFAccessor() = default;
287
300 template <int structdim2, int dim2, int spacedim2>
302
307 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
310
314 template <bool level_dof_access2>
316
328 delete;
329
334 operator=( // NOLINT
336 default; // NOLINT
337
347
351 template <bool level_dof_access2>
352 void
354
359 void
361
366 static bool
368
380 child(const unsigned int c) const;
381
387 typename ::internal::DoFHandlerImplementation::
388 Iterators<dim, spacedim, level_dof_access>::line_iterator
389 line(const unsigned int i) const;
390
396 typename ::internal::DoFHandlerImplementation::
397 Iterators<dim, spacedim, level_dof_access>::quad_iterator
398 quad(const unsigned int i) const;
399
445 void
447 std::vector<types::global_dof_index> &dof_indices,
448 const types::fe_index fe_index = numbers::invalid_fe_index) const;
449
456 void
458 const int level,
459 std::vector<types::global_dof_index> &dof_indices,
460 const types::fe_index fe_index = numbers::invalid_fe_index) const;
461
465 void
467 const int level,
468 const std::vector<types::global_dof_index> &dof_indices,
470
494 const unsigned int vertex,
495 const unsigned int i,
496 const types::fe_index fe_index = numbers::invalid_fe_index) const;
497
505 const int level,
506 const unsigned int vertex,
507 const unsigned int i,
508 const types::fe_index fe_index = numbers::invalid_fe_index) const;
509
538 dof_index(const unsigned int i,
539 const types::fe_index fe_index = numbers::invalid_fe_index) const;
540
545 mg_dof_index(const int level, const unsigned int i) const;
546
568 unsigned int
570
579 nth_active_fe_index(const unsigned int n) const;
580
587 std::set<types::fe_index>
589
599 bool
600 fe_index_is_active(const types::fe_index fe_index) const;
601
608 get_fe(const types::fe_index fe_index) const;
609
620 "This accessor object has not been "
621 "associated with any DoFHandler object.");
654
655protected:
660
661public:
675 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
676 bool
679
683 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
684 bool
687
688protected:
692 void
694
713 void
715 const unsigned int i,
716 const types::global_dof_index index,
717 const types::fe_index fe_index = numbers::invalid_fe_index) const;
718
719 void
721 const unsigned int i,
722 const types::global_dof_index index) const;
723
724 void
726 const int level,
727 const unsigned int vertex,
728 const unsigned int i,
729 const types::global_dof_index index,
730 const types::fe_index fe_index = numbers::invalid_fe_index) const;
731
732 // Iterator classes need to be friends because they need to access
733 // operator== and operator!=.
734 template <typename>
735 friend class TriaRawIterator;
736 template <int, int, int, bool>
737 friend class DoFAccessor;
738
739private:
740 // Make the DoFHandler class a friend so that it can call the set_xxx()
741 // functions.
742 friend class DoFHandler<dim, spacedim>;
743
744 friend struct ::internal::DoFHandlerImplementation::Policy::
745 Implementation;
746 friend struct ::internal::DoFHandlerImplementation::Implementation;
747 friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
748 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
749 friend struct ::internal::DoFAccessorImplementation::Implementation;
750};
751
752
753
761template <int spacedim, bool level_dof_access>
762class DoFAccessor<0, 1, spacedim, level_dof_access>
763 : public TriaAccessor<0, 1, spacedim>
764{
765public:
770 static constexpr unsigned int dimension = 1;
771
776 static constexpr unsigned int space_dimension = spacedim;
777
783
788
800
818 const Triangulation<1, spacedim> *tria,
819 const typename TriaAccessor<0, 1, spacedim>::VertexKind vertex_kind,
820 const unsigned int vertex_index,
822
830 const int = 0,
831 const int = 0,
833
846 template <int structdim2, int dim2, int spacedim2>
848
853 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
856
861
865 // NOLINTNEXTLINE OSX does not compile with noexcept
867
871 ~DoFAccessor() = default;
872
884
890 default;
891
899 const DoFHandler<1, spacedim> &
901
905 template <bool level_dof_access2>
906 void
907 copy_from(const DoFAccessor<0, 1, spacedim, level_dof_access2> &a);
908
913 void
914 copy_from(const TriaAccessorBase<0, 1, spacedim> &da);
915
928 TriaIterator<DoFAccessor<0, 1, spacedim, level_dof_access>>
929 child(const unsigned int c) const;
930
937 typename ::internal::DoFHandlerImplementation::
938 Iterators<1, spacedim, level_dof_access>::line_iterator
939 line(const unsigned int i) const;
940
947 typename ::internal::DoFHandlerImplementation::
948 Iterators<1, spacedim, level_dof_access>::quad_iterator
949 quad(const unsigned int i) const;
950
994 void
996 std::vector<types::global_dof_index> &dof_indices,
997 const types::fe_index fe_index = numbers::invalid_fe_index) const;
998
1005 void
1007 const int level,
1008 std::vector<types::global_dof_index> &dof_indices,
1009 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1010
1029 types::global_dof_index
1031 const unsigned int vertex,
1032 const unsigned int i,
1033 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1034
1052 types::global_dof_index
1053 dof_index(const unsigned int i,
1054 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1055
1074 unsigned int
1076
1084 types::fe_index
1085 nth_active_fe_index(const unsigned int n) const;
1086
1095 bool
1096 fe_index_is_active(const types::fe_index fe_index) const;
1097
1103 const FiniteElement<1, spacedim> &
1104 get_fe(const types::fe_index fe_index) const;
1105
1116 "This accessor object has not been "
1117 "associated with any DoFHandler object.");
1150
1151protected:
1156
1160 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1161 bool
1162 operator==(
1163 const DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1164
1168 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1169 bool
1170 operator!=(
1171 const DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1172
1176 void
1177 set_dof_handler(DoFHandler<1, spacedim> *dh);
1178
1197 void
1199 const unsigned int i,
1200 const types::global_dof_index index,
1201 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1202
1203 // Iterator classes need to be friends because they need to access
1204 // operator== and operator!=.
1205 template <typename>
1206 friend class TriaRawIterator;
1207
1208
1209 // Make the DoFHandler class a friend so that it can call the set_xxx()
1210 // functions.
1211 friend class DoFHandler<1, spacedim>;
1212
1213 friend struct ::internal::DoFHandlerImplementation::Policy::
1215 friend struct ::internal::DoFHandlerImplementation::Implementation;
1216 friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
1217 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
1218};
1219
1220
1221
1222/* -------------------------------------------------------------------------- */
1223
1224
1245template <int structdim, int dim, int spacedim = dim>
1246class DoFInvalidAccessor : public InvalidAccessor<structdim, dim, spacedim>
1247{
1248public:
1254
1262 DoFInvalidAccessor(const void *parent = nullptr,
1263 const int level = -1,
1264 const int index = -1,
1265 const AccessorData *local_data = nullptr);
1266
1275
1280 template <typename OtherAccessor>
1281 DoFInvalidAccessor(const OtherAccessor &);
1282
1289 dof_index(const unsigned int i,
1290 const types::fe_index fe_index =
1292
1298 void
1300 const unsigned int i,
1301 const types::global_dof_index index,
1302 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1303};
1304
1305
1306
1307/* -------------------------------------------------------------------------- */
1308
1309
1321template <int dimension_, int space_dimension_, bool level_dof_access>
1322class DoFCellAccessor : public DoFAccessor<dimension_,
1323 dimension_,
1324 space_dimension_,
1325 level_dof_access>
1326{
1327public:
1331 static const unsigned int dim = dimension_;
1332
1336 static const unsigned int spacedim = space_dimension_;
1337
1338
1343
1350
1355
1361 dimension_,
1362 space_dimension_,
1363 level_dof_access>>;
1364
1376 const int level,
1377 const int index,
1378 const AccessorData *local_data);
1379
1392 template <int structdim2, int dim2, int spacedim2>
1394
1399 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1402
1408 default;
1409
1415 &&) = default; // NOLINT
1416
1420 ~DoFCellAccessor() = default;
1421
1434 delete;
1435
1440 operator=( // NOLINT
1442 &&) = default; // NOLINT
1443
1458 parent() const;
1459
1473 neighbor(const unsigned int i) const;
1474
1481 periodic_neighbor(const unsigned int i) const;
1482
1489 neighbor_or_periodic_neighbor(const unsigned int i) const;
1490
1497 child(const unsigned int i) const;
1498
1502 boost::container::small_vector<
1507
1515 face(const unsigned int i) const;
1516
1520 boost::container::small_vector<face_iterator,
1523
1531 neighbor_child_on_subface(const unsigned int face_no,
1532 const unsigned int subface_no) const;
1533
1541 periodic_neighbor_child_on_subface(const unsigned int face_no,
1542 const unsigned int subface_no) const;
1543
1572 template <class InputVector, typename number>
1573 void
1574 get_dof_values(const InputVector &values, Vector<number> &local_values) const;
1575
1593 template <typename Number, typename ForwardIterator>
1594 void
1596 ForwardIterator local_values_begin,
1597 ForwardIterator local_values_end) const;
1598
1619 template <class InputVector, typename ForwardIterator>
1620 void
1623 const InputVector &values,
1624 ForwardIterator local_values_begin,
1625 ForwardIterator local_values_end) const;
1626
1651 template <class OutputVector, typename number>
1652 void
1653 set_dof_values(const Vector<number> &local_values,
1654 OutputVector &values) const;
1655
1687 template <typename Number>
1688 void
1689 get_interpolated_dof_values(
1690 const ReadVector<Number> &values,
1691 Vector<Number> &interpolated_values,
1692 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1693
1755 template <class OutputVector, typename number>
1756 void
1757 set_dof_values_by_interpolation(
1758 const Vector<number> &local_values,
1759 OutputVector &values,
1761 const bool perform_check = false) const;
1762
1772 template <class OutputVector, typename number>
1773 void
1774 distribute_local_to_global_by_interpolation(
1775 const Vector<number> &local_values,
1776 OutputVector &values,
1777 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1778
1793 template <typename number, typename OutputVector>
1794 void
1796 OutputVector &global_destination) const;
1797
1812 template <typename ForwardIterator, typename OutputVector>
1813 void
1814 distribute_local_to_global(ForwardIterator local_source_begin,
1815 ForwardIterator local_source_end,
1816 OutputVector &global_destination) const;
1817
1831 template <typename ForwardIterator, typename OutputVector>
1832 void
1835 ForwardIterator local_source_begin,
1836 ForwardIterator local_source_end,
1837 OutputVector &global_destination) const;
1838
1845 template <typename number, typename OutputMatrix>
1846 void
1848 OutputMatrix &global_destination) const;
1849
1854 template <typename number, typename OutputMatrix, typename OutputVector>
1855 void
1857 const Vector<number> &local_vector,
1858 OutputMatrix &global_matrix,
1859 OutputVector &global_vector) const;
1860
1885 void
1887 std::vector<types::global_dof_index> &dof_indices) const;
1888
1920 void
1921 get_dof_indices(std::vector<types::global_dof_index> &dof_indices) const;
1922
1927 void
1928 get_mg_dof_indices(std::vector<types::global_dof_index> &dof_indices) const;
1929
1955 get_fe() const;
1956
1984
2011 void
2021 void
2022 set_dof_indices(const std::vector<types::global_dof_index> &dof_indices);
2023
2027 void
2028 set_mg_dof_indices(const std::vector<types::global_dof_index> &dof_indices);
2029
2057
2078
2087 void
2089
2096 bool
2098
2107 void
2113private:
2114 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
2115};
2116
2117
2118template <int structdim, int dim, int spacedim, bool level_dof_access>
2119inline bool
2124
2125
2126
2127template <int structdim, int dim, int spacedim>
2128template <typename OtherAccessor>
2130 const OtherAccessor &)
2131{
2132 Assert(false,
2133 ExcMessage("You are attempting an illegal conversion between "
2134 "iterator/accessor types. The constructor you call "
2135 "only exists to make certain template constructs "
2136 "easier to write as dimension independent code but "
2137 "the conversion is not valid in the current context."));
2138}
2139
2140
2141
2143
2144// include more templates
2145#include <deal.II/dofs/dof_accessor.templates.h>
2146
2147
2148#endif
DoFAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
DoFAccessor< 0, 1, spacedim, level_dof_access > & operator=(const DoFAccessor< 0, 1, spacedim, level_dof_access > &da)=delete
DoFAccessor< 0, 1, spacedim, level_dof_access > & operator=(DoFAccessor< 0, 1, spacedim, level_dof_access > &&) noexcept=default
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(const DoFAccessor< 0, 1, spacedim, level_dof_access > &)=default
DoFAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
typename::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::line_iterator line(const unsigned int i) const
friend struct ::internal::DoFHandlerImplementation::Policy::Implementation
static constexpr unsigned int space_dimension
const DoFHandler< dim, spacedim > & get_dof_handler() const
void set_mg_dof_index(const int level, const unsigned int i, const types::global_dof_index index) const
DoFAccessor< structdim, dim, spacedim, level_dof_access > & operator=(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &da)=delete
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 > &)
TriaIterator< DoFAccessor< structdim, dim, spacedim, level_dof_access > > child(const unsigned int c) 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
void set_dof_handler(DoFHandler< dim, spacedim > *dh)
DoFAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
DoFAccessor< structdim, dim, spacedim, level_dof_access > & operator=(DoFAccessor< structdim, dim, spacedim, level_dof_access > &&)
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
typename ::internal::DoFAccessorImplementation:: Inheritance< structdim, dimension, space_dimension >::BaseClass BaseClass
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
typename::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::quad_iterator quad(const unsigned int i) const
void copy_from(const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &a)
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index 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
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
static bool is_level_cell()
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
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(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< dimension_, space_dimension_, level_dof_access > & operator=(const DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &da)=delete
const FiniteElement< dimension_, space_dimension_ > & get_fe() const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > child(const unsigned int i) const
DoFCellAccessor(const Triangulation< dimension_, space_dimension_ > *tria, const int level, const int index, const AccessorData *local_data)
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
const FiniteElement< dimension_, space_dimension_ > & get_future_fe() 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 > > neighbor(const unsigned int i) const
void set_active_fe_index(const types::fe_index i) const
DoFCellAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
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
boost::container::small_vector< face_iterator, GeometryInfo< dimension_ >::faces_per_cell > face_iterators() const
DoFCellAccessor(const DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &)=default
void distribute_local_to_global(ForwardIterator local_source_begin, ForwardIterator local_source_end, OutputVector &global_destination) const
void get_dof_values(const ReadVector< Number > &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) const
~DoFCellAccessor()=default
bool future_fe_index_set() const
DoFCellAccessor< dimension_, space_dimension_, level_dof_access > & operator=(DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &&)=default
types::fe_index active_fe_index() const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > parent() const
typename InvalidAccessor< structdim, dim, spacedim >::AccessorData AccessorData
DoFInvalidAccessor(const void *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
unsigned int level
Definition grid_out.cc:4626
static ::ExceptionBase & ExcVectorNotEmpty()
#define DeclException0(Exception0)
Definition exceptions.h:471
#define Assert(cond, exc)
static ::ExceptionBase & ExcNotActive()
static ::ExceptionBase & ExcMatrixDoesNotMatch()
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:494
static ::ExceptionBase & ExcInvalidObject()
static ::ExceptionBase & ExcCantCompareIterators()
static ::ExceptionBase & ExcVectorDoesNotMatch()
static ::ExceptionBase & ExcMessage(std::string arg1)
Definition hp.h:117
const types::fe_index invalid_fe_index
Definition types.h:243
STL namespace.
Definition types.h:32
unsigned short int fe_index
Definition types.h:59