Processing math: 25%
 Reference documentation for deal.II version 9.4.1
\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\}}
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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
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
211template <int structdim, int dim, int spacedim, bool level_dof_access>
212class DoFAccessor : public ::internal::DoFAccessorImplementation::
213 Inheritance<structdim, dim, spacedim>::BaseClass
214{
215public:
220 static constexpr unsigned int dimension = dim;
221
226 static constexpr unsigned int space_dimension = spacedim;
227
232 using BaseClass = typename ::internal::DoFAccessorImplementation::
233 Inheritance<structdim, dimension, space_dimension>::BaseClass;
234
239
251
268 const int level,
269 const int index,
271
276 default;
277
281 DoFAccessor( // NOLINT
283 default; // NOLINT
284
288 ~DoFAccessor() = default;
289
302 template <int structdim2, int dim2, int spacedim2>
304
309 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
312
316 template <bool level_dof_access2>
318
330 delete;
331
336 operator=( // NOLINT
338 default; // NOLINT
339
349
353 template <bool level_dof_access2>
354 void
356
361 void
363
368 static bool
370
382 child(const unsigned int c) const;
383
389 typename ::internal::DoFHandlerImplementation::
390 Iterators<dim, spacedim, level_dof_access>::line_iterator
391 line(const unsigned int i) const;
392
398 typename ::internal::DoFHandlerImplementation::
399 Iterators<dim, spacedim, level_dof_access>::quad_iterator
400 quad(const unsigned int i) const;
401
447 void
448 get_dof_indices(std::vector<types::global_dof_index> &dof_indices,
449 const unsigned int fe_index =
451
458 void
460 std::vector<types::global_dof_index> &dof_indices,
461 const unsigned int fe_index =
463
467 void
469 const int level,
470 const std::vector<types::global_dof_index> &dof_indices,
471 const unsigned int fe_index = DoFHandler<dim, spacedim>::invalid_fe_index);
472
495 vertex_dof_index(const unsigned int vertex,
496 const unsigned int i,
497 const unsigned int fe_index =
499
507 const unsigned int vertex,
508 const unsigned int i,
509 const unsigned int fe_index =
511
540 dof_index(const unsigned int i,
541 const unsigned int fe_index =
543
548 mg_dof_index(const int level, const unsigned int i) const;
549
571 unsigned int
573
581 unsigned int
582 nth_active_fe_index(const unsigned int n) const;
583
590 std::set<unsigned int>
592
602 bool
603 fe_index_is_active(const unsigned int fe_index) const;
604
611 get_fe(const unsigned int fe_index) const;
612
623 "This accessor object has not been "
624 "associated with any DoFHandler object.");
657
658protected:
663
664public:
678 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
679 bool
682
686 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
687 bool
690
691protected:
695 void
697
716 void
717 set_dof_index(const unsigned int i,
718 const types::global_dof_index index,
719 const unsigned int fe_index =
721
722 void
724 const unsigned int i,
725 const types::global_dof_index index) const;
726
727 void
729 const unsigned int vertex,
730 const unsigned int i,
731 const types::global_dof_index index,
732 const unsigned int fe_index =
734
735 // Iterator classes need to be friends because they need to access
736 // operator== and operator!=.
737 template <typename>
738 friend class TriaRawIterator;
739 template <int, int, int, bool>
740 friend class DoFAccessor;
741
742private:
743 // Make the DoFHandler class a friend so that it can call the set_xxx()
744 // functions.
745 template <int, int>
746 friend class DoFHandler;
747
748 friend struct ::internal::DoFHandlerImplementation::Policy::
749 Implementation;
750 friend struct ::internal::DoFHandlerImplementation::Implementation;
751 friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
752 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
753 friend struct ::internal::DoFAccessorImplementation::Implementation;
754};
755
756
757
765template <int spacedim, bool level_dof_access>
766class DoFAccessor<0, 1, spacedim, level_dof_access>
767 : public TriaAccessor<0, 1, spacedim>
768{
769public:
774 static constexpr unsigned int dimension = 1;
775
780 static constexpr unsigned int space_dimension = spacedim;
781
787
792
804
823 const typename TriaAccessor<0, 1, spacedim>::VertexKind vertex_kind,
824 const unsigned int vertex_index,
826
834 const int = 0,
835 const int = 0,
837
850 template <int structdim2, int dim2, int spacedim2>
852
857 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
860
865
869 // NOLINTNEXTLINE OSX does not compile with noexcept
871
875 ~DoFAccessor() = default;
876
888
894 default;
895
903 const DoFHandler<1, spacedim> &
905
909 template <bool level_dof_access2>
910 void
911 copy_from(const DoFAccessor<0, 1, spacedim, level_dof_access2> &a);
912
917 void
918 copy_from(const TriaAccessorBase<0, 1, spacedim> &da);
919
932 TriaIterator<DoFAccessor<0, 1, spacedim, level_dof_access>>
933 child(const unsigned int c) const;
934
941 typename ::internal::DoFHandlerImplementation::
942 Iterators<1, spacedim, level_dof_access>::line_iterator
943 line(const unsigned int i) const;
944
951 typename ::internal::DoFHandlerImplementation::
952 Iterators<1, spacedim, level_dof_access>::quad_iterator
953 quad(const unsigned int i) const;
954
998 void
1000 std::vector<types::global_dof_index> &dof_indices,
1001 const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1002
1009 void
1011 const int level,
1012 std::vector<types::global_dof_index> &dof_indices,
1013 const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1014
1033 types::global_dof_index
1035 const unsigned int vertex,
1036 const unsigned int i,
1037 const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1038
1056 types::global_dof_index
1057 dof_index(const unsigned int i,
1058 const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1059
1078 unsigned int
1080
1088 unsigned int
1089 nth_active_fe_index(const unsigned int n) const;
1090
1099 bool
1100 fe_index_is_active(const unsigned int fe_index) const;
1101
1107 const FiniteElement<1, spacedim> &
1108 get_fe(const unsigned int fe_index) const;
1109
1120 "This accessor object has not been "
1121 "associated with any DoFHandler object.");
1154
1155protected:
1160
1164 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1165 bool
1166 operator==(
1167 const DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1168
1172 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1173 bool
1174 operator!=(
1175 const DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1176
1180 void
1181 set_dof_handler(DoFHandler<1, spacedim> *dh);
1182
1201 void
1203 const unsigned int i,
1204 const types::global_dof_index index,
1205 const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1206
1207 // Iterator classes need to be friends because they need to access
1208 // operator== and operator!=.
1209 template <typename>
1210 friend class TriaRawIterator;
1211
1212
1213 // Make the DoFHandler class a friend so that it can call the set_xxx()
1214 // functions.
1215 template <int, int>
1216 friend class DoFHandler;
1217
1218 friend struct ::internal::DoFHandlerImplementation::Policy::
1220 friend struct ::internal::DoFHandlerImplementation::Implementation;
1221 friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
1222 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
1223};
1224
1225
1226
1227/* -------------------------------------------------------------------------- */
1228
1229
1250template <int structdim, int dim, int spacedim = dim>
1251class DoFInvalidAccessor : public InvalidAccessor<structdim, dim, spacedim>
1252{
1253public:
1259
1268 const int level = -1,
1269 const int index = -1,
1270 const AccessorData * local_data = 0);
1271
1280
1285 template <typename OtherAccessor>
1286 DoFInvalidAccessor(const OtherAccessor &);
1287
1294 dof_index(const unsigned int i,
1295 const unsigned int fe_index =
1297
1303 void
1304 set_dof_index(const unsigned int i,
1305 const types::global_dof_index index,
1306 const unsigned int fe_index =
1308};
1309
1310
1311
1312/* -------------------------------------------------------------------------- */
1313
1314
1326template <int dimension_, int space_dimension_, bool level_dof_access>
1327class DoFCellAccessor : public DoFAccessor<dimension_,
1328 dimension_,
1329 space_dimension_,
1330 level_dof_access>
1331{
1332public:
1336 static const unsigned int dim = dimension_;
1337
1341 static const unsigned int spacedim = space_dimension_;
1342
1343
1348
1355
1360
1366 dimension_,
1367 space_dimension_,
1368 level_dof_access>>;
1369
1381 const int level,
1382 const int index,
1383 const AccessorData *local_data);
1384
1397 template <int structdim2, int dim2, int spacedim2>
1399
1404 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1407
1413 default;
1414
1420 &&) = default; // NOLINT
1421
1425 ~DoFCellAccessor() = default;
1426
1439 delete;
1440
1445 operator=( // NOLINT
1447 &&) = default; // NOLINT
1448
1463 parent() const;
1464
1478 neighbor(const unsigned int i) const;
1479
1486 periodic_neighbor(const unsigned int i) const;
1487
1494 neighbor_or_periodic_neighbor(const unsigned int i) const;
1495
1502 child(const unsigned int i) const;
1503
1507 boost::container::small_vector<
1512
1520 face(const unsigned int i) const;
1521
1525 boost::container::small_vector<face_iterator,
1528
1536 neighbor_child_on_subface(const unsigned int face_no,
1537 const unsigned int subface_no) const;
1538
1546 periodic_neighbor_child_on_subface(const unsigned int face_no,
1547 const unsigned int subface_no) const;
1548
1577 template <class InputVector, typename number>
1578 void
1579 get_dof_values(const InputVector &values, Vector<number> &local_values) const;
1580
1598 template <class InputVector, typename ForwardIterator>
1599 void
1600 get_dof_values(const InputVector &values,
1601 ForwardIterator local_values_begin,
1602 ForwardIterator local_values_end) const;
1603
1624 template <class InputVector, typename ForwardIterator>
1625 void
1628 const InputVector & values,
1629 ForwardIterator local_values_begin,
1630 ForwardIterator local_values_end) const;
1631
1656 template <class OutputVector, typename number>
1657 void
1658 set_dof_values(const Vector<number> &local_values,
1659 OutputVector & values) const;
1660
1692 template <class InputVector, typename number>
1693 void
1694 get_interpolated_dof_values(
1695 const InputVector &values,
1696 Vector<number> & interpolated_values,
1697 const unsigned int fe_index =
1699
1761 template <class OutputVector, typename number>
1762 void
1763 set_dof_values_by_interpolation(
1764 const Vector<number> &local_values,
1765 OutputVector & values,
1766 const unsigned int fe_index =
1768 const bool perform_check = false) const;
1769
1779 template <class OutputVector, typename number>
1780 void
1781 distribute_local_to_global_by_interpolation(
1782 const Vector<number> &local_values,
1783 OutputVector & values,
1784 const unsigned int fe_index =
1786
1801 template <typename number, typename OutputVector>
1802 void
1804 OutputVector & global_destination) const;
1805
1820 template <typename ForwardIterator, typename OutputVector>
1821 void
1822 distribute_local_to_global(ForwardIterator local_source_begin,
1823 ForwardIterator local_source_end,
1824 OutputVector & global_destination) const;
1825
1839 template <typename ForwardIterator, typename OutputVector>
1840 void
1843 ForwardIterator local_source_begin,
1844 ForwardIterator local_source_end,
1845 OutputVector & global_destination) const;
1846
1853 template <typename number, typename OutputMatrix>
1854 void
1856 OutputMatrix &global_destination) const;
1857
1862 template <typename number, typename OutputMatrix, typename OutputVector>
1863 void
1865 const Vector<number> & local_vector,
1866 OutputMatrix & global_matrix,
1867 OutputVector & global_vector) const;
1868
1893 void
1895 std::vector<types::global_dof_index> &dof_indices) const;
1896
1928 void
1929 get_dof_indices(std::vector<types::global_dof_index> &dof_indices) const;
1930
1935 void
1936 get_mg_dof_indices(std::vector<types::global_dof_index> &dof_indices) const;
1937
1963 get_fe() const;
1964
1990 unsigned int
1992
2019 void
2020 set_active_fe_index(const unsigned int i) const;
2029 void
2030 set_dof_indices(const std::vector<types::global_dof_index> &dof_indices);
2031
2035 void
2036 set_mg_dof_indices(const std::vector<types::global_dof_index> &dof_indices);
2037
2041 void
2042 update_cell_dof_indices_cache() const;
2043
2071
2090 unsigned int
2092
2101 void
2102 set_future_fe_index(const unsigned int i) const;
2103
2110 bool
2112
2121 void
2127private:
2128 // Make the DoFHandler class a friend so that it can call the
2129 // update_cell_dof_indices_cache() function
2130 template <int, int>
2131 friend class DoFHandler;
2132 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
2133};
2134
2135
2136template <int structdim, int dim, int spacedim, bool level_dof_access>
2137inline bool
2139{
2140 return level_dof_access;
2141}
2142
2143
2144
2145template <int structdim, int dim, int spacedim>
2146template <typename OtherAccessor>
2148 const OtherAccessor &)
2149{
2150 Assert(false,
2151 ExcMessage("You are attempting an illegal conversion between "
2152 "iterator/accessor types. The constructor you call "
2153 "only exists to make certain template constructs "
2154 "easier to write as dimension independent code but "
2155 "the conversion is not valid in the current context."));
2156}
2157
2158
2159
2161
2162// include more templates
2163#include "dof_accessor.templates.h"
2164
2165
2166#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
Definition: dof_accessor.h:749
static constexpr unsigned int space_dimension
Definition: dof_accessor.h:226
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 constexpr unsigned int dimension
Definition: dof_accessor.h:220
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
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:233
DoFHandler< dim, spacedim > * dof_handler
Definition: dof_accessor.h:662
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:109
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition: config.h:456
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition: config.h:495
unsigned int level
Definition: grid_out.cc:4606
static ::ExceptionBase & ExcVectorNotEmpty()
#define DeclException0(Exception0)
Definition: exceptions.h:464
#define Assert(cond, exc)
Definition: exceptions.h:1473
static ::ExceptionBase & ExcNotActive()
static ::ExceptionBase & ExcMatrixDoesNotMatch()
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:487
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
const ::Triangulation< dim, spacedim > & tria