Reference documentation for deal.II version 9.3.3
\(\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 - 2021 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
28
30#include <boost/container/small_vector.hpp>
32
33#include <vector>
34
36
37// Forward declarations
38#ifndef DOXYGEN
39template <typename number>
40class FullMatrix;
41template <typename number>
42class Vector;
43template <typename number>
45
46template <typename Accessor>
47class TriaRawIterator;
48
49template <int, int>
50class FiniteElement;
51
52namespace internal
53{
54 namespace DoFCellAccessorImplementation
55 {
56 struct Implementation;
57 }
58
59 namespace DoFHandlerImplementation
60 {
61 struct Implementation;
62 namespace Policy
63 {
64 struct Implementation;
65 }
66 } // namespace DoFHandlerImplementation
67
68 namespace hp
69 {
70 namespace DoFHandlerImplementation
71 {
72 struct Implementation;
73 }
74 } // namespace hp
75} // namespace internal
76#endif
77
78// note: the file dof_accessor.templates.h is included at the end of
79// this file. this includes a lot of templates and thus makes
80// compilation slower, but at the same time allows for more aggressive
81// inlining and thus faster code.
82
83
84namespace internal
85{
86 namespace DoFAccessorImplementation
87 {
104 template <int structdim, int dim, int spacedim>
106 {
112 };
113
114
120 template <int dim, int spacedim>
121 struct Inheritance<dim, dim, spacedim>
122 {
128 };
129
130 struct Implementation;
131 } // namespace DoFAccessorImplementation
132} // namespace internal
133
134
135/* -------------------------------------------------------------------------- */
136
137
138
210template <int structdim, int dim, int spacedim, bool level_dof_access>
211class DoFAccessor : public ::internal::DoFAccessorImplementation::
212 Inheritance<structdim, dim, spacedim>::BaseClass
213{
214public:
219 static const unsigned int dimension = dim;
220
225 static const unsigned int space_dimension = spacedim;
226
231 using BaseClass = typename ::internal::DoFAccessorImplementation::
232 Inheritance<structdim, dimension, space_dimension>::BaseClass;
233
238
250
267 const int level,
268 const int index,
270
275 default;
276
280 DoFAccessor( // NOLINT
282 default; // NOLINT
283
287 ~DoFAccessor() = default;
288
301 template <int structdim2, int dim2, int spacedim2>
303
308 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
311
315 template <bool level_dof_access2>
317
329 delete;
330
335 operator=( // NOLINT
337 default; // NOLINT
338
348
352 template <bool level_dof_access2>
353 void
355
360 void
362
367 static bool
369
381 child(const unsigned int c) const;
382
388 typename ::internal::DoFHandlerImplementation::
389 Iterators<dim, spacedim, level_dof_access>::line_iterator
390 line(const unsigned int i) const;
391
397 typename ::internal::DoFHandlerImplementation::
398 Iterators<dim, spacedim, level_dof_access>::quad_iterator
399 quad(const unsigned int i) const;
400
446 void
447 get_dof_indices(std::vector<types::global_dof_index> &dof_indices,
448 const unsigned int fe_index =
450
457 void
459 std::vector<types::global_dof_index> &dof_indices,
460 const unsigned int fe_index =
462
466 void
468 const int level,
469 const std::vector<types::global_dof_index> &dof_indices,
470 const unsigned int fe_index = DoFHandler<dim, spacedim>::invalid_fe_index);
471
490 vertex_dof_index(const unsigned int vertex,
491 const unsigned int i,
492 const unsigned int fe_index =
494
502 const unsigned int vertex,
503 const unsigned int i,
504 const unsigned int fe_index =
506
535 dof_index(const unsigned int i,
536 const unsigned int fe_index =
538
543 mg_dof_index(const int level, const unsigned int i) const;
544
566 unsigned int
568
576 unsigned int
577 nth_active_fe_index(const unsigned int n) const;
578
585 std::set<unsigned int>
587
597 bool
598 fe_index_is_active(const unsigned int fe_index) const;
599
606 get_fe(const unsigned int fe_index) const;
607
618 "This accessor object has not been "
619 "associated with any DoFHandler object.");
652
653protected:
658
659public:
673 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
674 bool
677
681 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
682 bool
685
686protected:
690 void
692
711 void
712 set_dof_index(const unsigned int i,
713 const types::global_dof_index index,
714 const unsigned int fe_index =
716
717 void
719 const unsigned int i,
720 const types::global_dof_index index) const;
721
740 void
741 set_vertex_dof_index(const unsigned int vertex,
742 const unsigned int i,
743 const types::global_dof_index index,
744 const unsigned int fe_index =
746
747 void
749 const unsigned int vertex,
750 const unsigned int i,
751 const types::global_dof_index index,
752 const unsigned int fe_index =
754
755 // Iterator classes need to be friends because they need to access
756 // operator== and operator!=.
757 template <typename>
758 friend class TriaRawIterator;
759 template <int, int, int, bool>
760 friend class DoFAccessor;
761
762private:
763 // Make the DoFHandler class a friend so that it can call the set_xxx()
764 // functions.
765 template <int, int>
766 friend class DoFHandler;
767
768 friend struct ::internal::DoFHandlerImplementation::Policy::
769 Implementation;
770 friend struct ::internal::DoFHandlerImplementation::Implementation;
771 friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
772 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
773 friend struct ::internal::DoFAccessorImplementation::Implementation;
774};
775
776
777
785template <int spacedim, bool level_dof_access>
786class DoFAccessor<0, 1, spacedim, level_dof_access>
787 : public TriaAccessor<0, 1, spacedim>
788{
789public:
794 static const unsigned int dimension = 1;
795
800 static const unsigned int space_dimension = spacedim;
801
807
812
824
842 const Triangulation<1, spacedim> * tria,
843 const typename TriaAccessor<0, 1, spacedim>::VertexKind vertex_kind,
844 const unsigned int vertex_index,
846
853 const int = 0,
854 const int = 0,
856
869 template <int structdim2, int dim2, int spacedim2>
871
876 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
879
884
888 // NOLINTNEXTLINE OSX does not compile with noexcept
890
894 ~DoFAccessor() = default;
895
907
913
921 const DoFHandler<1, spacedim> &
923
927 template <bool level_dof_access2>
928 void
929 copy_from(const DoFAccessor<0, 1, spacedim, level_dof_access2> &a);
930
935 void
936 copy_from(const TriaAccessorBase<0, 1, spacedim> &da);
937
950 TriaIterator<DoFAccessor<0, 1, spacedim, level_dof_access>>
951 child(const unsigned int c) const;
952
959 typename ::internal::DoFHandlerImplementation::
960 Iterators<1, spacedim, level_dof_access>::line_iterator
961 line(const unsigned int i) const;
962
969 typename ::internal::DoFHandlerImplementation::
970 Iterators<1, spacedim, level_dof_access>::quad_iterator
971 quad(const unsigned int i) const;
972
1016 void
1018 std::vector<types::global_dof_index> &dof_indices,
1019 const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1020
1027 void
1029 const int level,
1030 std::vector<types::global_dof_index> &dof_indices,
1031 const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1032
1053 const unsigned int vertex,
1054 const unsigned int i,
1055 const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1056
1075 dof_index(const unsigned int i,
1076 const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1077
1096 unsigned int
1098
1106 unsigned int
1107 nth_active_fe_index(const unsigned int n) const;
1108
1117 bool
1118 fe_index_is_active(const unsigned int fe_index) const;
1119
1125 const FiniteElement<1, spacedim> &
1126 get_fe(const unsigned int fe_index) const;
1127
1138 "This accessor object has not been "
1139 "associated with any DoFHandler object.");
1172
1173protected:
1178
1182 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1183 bool
1184 operator==(
1185 const DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1186
1190 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1191 bool
1192 operator!=(
1193 const DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1194
1198 void set_dof_handler(DoFHandler<1, spacedim> *dh);
1199
1218 void
1220 const unsigned int i,
1221 const types::global_dof_index index,
1222 const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1223
1242 void
1244 const unsigned int vertex,
1245 const unsigned int i,
1246 const types::global_dof_index index,
1247 const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1248
1249 // Iterator classes need to be friends because they need to access
1250 // operator== and operator!=.
1251 template <typename>
1252 friend class TriaRawIterator;
1253
1254
1255 // Make the DoFHandler class a friend so that it can call the set_xxx()
1256 // functions.
1257 template <int, int>
1258 friend class DoFHandler;
1259
1260 friend struct ::internal::DoFHandlerImplementation::Policy::
1261 Implementation;
1262 friend struct ::internal::DoFHandlerImplementation::Implementation;
1263 friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
1264 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
1265};
1266
1267
1268
1269/* -------------------------------------------------------------------------- */
1270
1271
1292template <int structdim, int dim, int spacedim = dim>
1293class DoFInvalidAccessor : public InvalidAccessor<structdim, dim, spacedim>
1294{
1295public:
1301
1310 const int level = -1,
1311 const int index = -1,
1312 const AccessorData * local_data = 0);
1313
1322
1327 template <typename OtherAccessor>
1328 DoFInvalidAccessor(const OtherAccessor &);
1329
1336 dof_index(const unsigned int i,
1337 const unsigned int fe_index =
1339
1345 void
1346 set_dof_index(const unsigned int i,
1347 const types::global_dof_index index,
1348 const unsigned int fe_index =
1350};
1351
1352
1353
1354/* -------------------------------------------------------------------------- */
1355
1356
1368template <int dimension_, int space_dimension_, bool level_dof_access>
1369class DoFCellAccessor : public DoFAccessor<dimension_,
1370 dimension_,
1371 space_dimension_,
1372 level_dof_access>
1373{
1374public:
1378 static const unsigned int dim = dimension_;
1379
1383 static const unsigned int spacedim = space_dimension_;
1384
1385
1390
1397
1402
1408 dimension_,
1409 space_dimension_,
1410 level_dof_access>>;
1411
1423 const int level,
1424 const int index,
1425 const AccessorData *local_data);
1426
1439 template <int structdim2, int dim2, int spacedim2>
1441
1446 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1449
1455 default;
1456
1462 &&) = default; // NOLINT
1463
1467 ~DoFCellAccessor() = default;
1468
1481 delete;
1482
1487 operator=( // NOLINT
1489 &&) = default; // NOLINT
1490
1505 parent() const;
1506
1520 neighbor(const unsigned int i) const;
1521
1528 periodic_neighbor(const unsigned int i) const;
1529
1536 neighbor_or_periodic_neighbor(const unsigned int i) const;
1537
1544 child(const unsigned int i) const;
1545
1549 boost::container::small_vector<
1554
1562 face(const unsigned int i) const;
1563
1567 boost::container::small_vector<face_iterator,
1570
1578 neighbor_child_on_subface(const unsigned int face_no,
1579 const unsigned int subface_no) const;
1580
1588 periodic_neighbor_child_on_subface(const unsigned int face_no,
1589 const unsigned int subface_no) const;
1590
1619 template <class InputVector, typename number>
1620 void
1621 get_dof_values(const InputVector &values, Vector<number> &local_values) const;
1622
1640 template <class InputVector, typename ForwardIterator>
1641 void
1642 get_dof_values(const InputVector &values,
1643 ForwardIterator local_values_begin,
1644 ForwardIterator local_values_end) const;
1645
1666 template <class InputVector, typename ForwardIterator>
1667 void
1670 const InputVector & values,
1671 ForwardIterator local_values_begin,
1672 ForwardIterator local_values_end) const;
1673
1698 template <class OutputVector, typename number>
1699 void
1700 set_dof_values(const Vector<number> &local_values,
1701 OutputVector & values) const;
1702
1734 template <class InputVector, typename number>
1735 void
1736 get_interpolated_dof_values(
1737 const InputVector &values,
1738 Vector<number> & interpolated_values,
1739 const unsigned int fe_index =
1741
1794 template <class OutputVector, typename number>
1795 void
1796 set_dof_values_by_interpolation(
1797 const Vector<number> &local_values,
1798 OutputVector & values,
1799 const unsigned int fe_index =
1801
1816 template <typename number, typename OutputVector>
1817 void
1819 OutputVector & global_destination) const;
1820
1835 template <typename ForwardIterator, typename OutputVector>
1836 void
1837 distribute_local_to_global(ForwardIterator local_source_begin,
1838 ForwardIterator local_source_end,
1839 OutputVector & global_destination) const;
1840
1854 template <typename ForwardIterator, typename OutputVector>
1855 void
1858 ForwardIterator local_source_begin,
1859 ForwardIterator local_source_end,
1860 OutputVector & global_destination) const;
1861
1868 template <typename number, typename OutputMatrix>
1869 void
1871 OutputMatrix &global_destination) const;
1872
1877 template <typename number, typename OutputMatrix, typename OutputVector>
1878 void
1880 const Vector<number> & local_vector,
1881 OutputMatrix & global_matrix,
1882 OutputVector & global_vector) const;
1883
1908 void
1910 std::vector<types::global_dof_index> &dof_indices) const;
1911
1943 void
1944 get_dof_indices(std::vector<types::global_dof_index> &dof_indices) const;
1945
1950 void
1951 get_mg_dof_indices(std::vector<types::global_dof_index> &dof_indices) const;
1952
1978 get_fe() const;
1979
2005 unsigned int
2007
2034 void
2035 set_active_fe_index(const unsigned int i) const;
2044 void
2045 set_dof_indices(const std::vector<types::global_dof_index> &dof_indices);
2046
2050 void
2051 set_mg_dof_indices(const std::vector<types::global_dof_index> &dof_indices);
2052
2056 void
2057 update_cell_dof_indices_cache() const;
2058
2086
2105 unsigned int
2107
2116 void
2117 set_future_fe_index(const unsigned int i) const;
2118
2125 bool
2127
2136 void
2142private:
2143 // Make the DoFHandler class a friend so that it can call the
2144 // update_cell_dof_indices_cache() function
2145 template <int, int>
2146 friend class DoFHandler;
2147 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
2148};
2149
2150
2151template <int structdim, int dim, int spacedim, bool level_dof_access>
2152inline bool
2154{
2155 return level_dof_access;
2156}
2157
2158
2159
2160template <int structdim, int dim, int spacedim>
2161template <typename OtherAccessor>
2163 const OtherAccessor &)
2164{
2165 Assert(false,
2166 ExcMessage("You are attempting an illegal conversion between "
2167 "iterator/accessor types. The constructor you call "
2168 "only exists to make certain template constructs "
2169 "easier to write as dimension independent code but "
2170 "the conversion is not valid in the current context."));
2171}
2172
2173
2174
2176
2177// include more templates
2178#include "dof_accessor.templates.h"
2179
2180
2181#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(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
static const unsigned int space_dimension
Definition: dof_accessor.h:225
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
DoFAccessor(DoFAccessor< structdim, dim, spacedim, level_dof_access > &&)
DoFAccessor(const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &)
void set_mg_dof_indices(const int level, const std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_fe_index)
TriaIterator< DoFAccessor< structdim, dim, spacedim, level_dof_access > > child(const unsigned int c) const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int fe_index) const
types::global_dof_index vertex_dof_index(const unsigned int vertex, const unsigned int i, const unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_fe_index) const
~DoFAccessor()=default
static const unsigned int dimension
Definition: dof_accessor.h:219
void set_dof_handler(DoFHandler< dim, spacedim > *dh)
DoFAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
types::global_dof_index dof_index(const unsigned int i, const unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_fe_index) const
DoFAccessor< structdim, dim, spacedim, level_dof_access > & operator=(DoFAccessor< structdim, dim, spacedim, level_dof_access > &&)
DoFAccessor(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &)=default
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 get_dof_indices(std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_fe_index) const
void copy_from(const TriaAccessorBase< structdim, dim, spacedim > &da)
typename::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::quad_iterator quad(const unsigned int i) const
bool fe_index_is_active(const unsigned int fe_index) const
void copy_from(const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &a)
void set_dof_index(const unsigned int i, const types::global_dof_index index, const unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_fe_index) const
void set_vertex_dof_index(const unsigned int vertex, const unsigned int i, const types::global_dof_index index, const unsigned int fe_index=DoFHandler< dim, spacedim >::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
types::global_dof_index mg_vertex_dof_index(const int level, const unsigned int vertex, const unsigned int i, const unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_fe_index) const
unsigned int n_active_fe_indices() const
unsigned int nth_active_fe_index(const unsigned int n) const
static bool is_level_cell()
void get_mg_dof_indices(const int level, std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_fe_index) const
typename ::internal::DoFAccessorImplementation::Inheritance< structdim, dimension, space_dimension >::BaseClass BaseClass
Definition: dof_accessor.h:232
DoFHandler< dim, spacedim > * dof_handler
Definition: dof_accessor.h:657
DoFAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
std::set< unsigned int > get_active_fe_indices() 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 unsigned int fe_index=DoFHandler< dim, spacedim >::invalid_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
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_active_fe_index(const unsigned int i) const
void get_dof_values(const InputVector &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) 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
unsigned int active_fe_index() 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
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_future_fe_index(const unsigned int i) const
void set_dof_values(const Vector< number > &local_values, OutputVector &values) const
void set_mg_dof_indices(const std::vector< types::global_dof_index > &dof_indices)
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 get_mg_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
void distribute_local_to_global(ForwardIterator local_source_begin, ForwardIterator local_source_end, OutputVector &global_destination) const
~DoFCellAccessor()=default
unsigned int future_fe_index() const
bool future_fe_index_set() const
DoFCellAccessor< dimension_, space_dimension_, level_dof_access > & operator=(DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &&)=default
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > parent() const
DoFInvalidAccessor(const Triangulation< dim, spacedim > *parent=0, const int level=-1, const int index=-1, const AccessorData *local_data=0)
Definition: dof_accessor.cc:46
typename InvalidAccessor< structdim, dim, spacedim >::AccessorData AccessorData
typename TriaAccessorBase< structdim, dim, spacedim >::AccessorData AccessorData
Definition: vector.h:110
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition: config.h:416
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition: config.h:454
unsigned int level
Definition: grid_out.cc:4590
static ::ExceptionBase & ExcVectorNotEmpty()
#define DeclException0(Exception0)
Definition: exceptions.h:470
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcNotActive()
static ::ExceptionBase & ExcMatrixDoesNotMatch()
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
static ::ExceptionBase & ExcInvalidObject()
static ::ExceptionBase & ExcCantCompareIterators()
static ::ExceptionBase & ExcVectorDoesNotMatch()
static ::ExceptionBase & ExcMessage(std::string arg1)
Definition: hp.h:118
STL namespace.
Definition: types.h:32
unsigned int global_dof_index
Definition: types.h:76