Reference documentation for deal.II version GIT d7aca55de5 2022-08-10 12: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 
28 #include <boost/container/small_vector.hpp>
30 
31 #include <vector>
32 
34 
35 // Forward declarations
36 #ifndef DOXYGEN
37 template <typename number>
38 class FullMatrix;
39 template <typename number>
40 class Vector;
41 template <typename number>
42 class AffineConstraints;
43 
44 template <typename Accessor>
45 class TriaRawIterator;
46 
47 template <int, int>
48 class FiniteElement;
49 
50 namespace 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 
82 namespace internal
83 {
84  namespace DoFAccessorImplementation
85  {
102  template <int structdim, int dim, int spacedim>
103  struct Inheritance
104  {
110  };
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 
209 template <int structdim, int dim, int spacedim, bool level_dof_access>
210 class DoFAccessor : public ::internal::DoFAccessorImplementation::
211  Inheritance<structdim, dim, spacedim>::BaseClass
212 {
213 public:
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
446  get_dof_indices(std::vector<types::global_dof_index> &dof_indices,
447  const unsigned int fe_index =
449 
456  void
457  get_mg_dof_indices(const int level,
458  std::vector<types::global_dof_index> &dof_indices,
459  const unsigned int fe_index =
461 
465  void
467  const int level,
468  const std::vector<types::global_dof_index> &dof_indices,
469  const unsigned int fe_index = DoFHandler<dim, spacedim>::invalid_fe_index);
470 
493  vertex_dof_index(const unsigned int vertex,
494  const unsigned int i,
495  const unsigned int fe_index =
497 
505  const unsigned int vertex,
506  const unsigned int i,
507  const unsigned int fe_index =
509 
538  dof_index(const unsigned int i,
539  const unsigned int fe_index =
541 
546  mg_dof_index(const int level, const unsigned int i) const;
547 
569  unsigned int
571 
579  unsigned int
580  nth_active_fe_index(const unsigned int n) const;
581 
588  std::set<unsigned int>
590 
600  bool
601  fe_index_is_active(const unsigned int fe_index) const;
602 
609  get_fe(const unsigned int fe_index) const;
610 
621  "This accessor object has not been "
622  "associated with any DoFHandler object.");
655 
656 protected:
661 
662 public:
676  template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
677  bool
680 
684  template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
685  bool
688 
689 protected:
693  void
695 
714  void
715  set_dof_index(const unsigned int i,
716  const types::global_dof_index index,
717  const unsigned int fe_index =
719 
720  void
722  const unsigned int i,
723  const types::global_dof_index index) const;
724 
725  void
727  const unsigned int vertex,
728  const unsigned int i,
729  const types::global_dof_index index,
730  const unsigned int fe_index =
732 
733  // Iterator classes need to be friends because they need to access
734  // operator== and operator!=.
735  template <typename>
736  friend class TriaRawIterator;
737  template <int, int, int, bool>
738  friend class DoFAccessor;
739 
740 private:
741  // Make the DoFHandler class a friend so that it can call the set_xxx()
742  // functions.
743  template <int, int>
744  friend class DoFHandler;
745 
746  friend struct ::internal::DoFHandlerImplementation::Policy::
747  Implementation;
748  friend struct ::internal::DoFHandlerImplementation::Implementation;
749  friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
750  friend struct ::internal::DoFCellAccessorImplementation::Implementation;
751  friend struct ::internal::DoFAccessorImplementation::Implementation;
752 };
753 
754 
755 
763 template <int spacedim, bool level_dof_access>
764 class DoFAccessor<0, 1, spacedim, level_dof_access>
765  : public TriaAccessor<0, 1, spacedim>
766 {
767 public:
772  static constexpr unsigned int dimension = 1;
773 
778  static constexpr unsigned int space_dimension = spacedim;
779 
785 
790 
802 
821  const typename TriaAccessor<0, 1, spacedim>::VertexKind vertex_kind,
822  const unsigned int vertex_index,
824 
832  const int = 0,
833  const int = 0,
835 
848  template <int structdim2, int dim2, int spacedim2>
850 
855  template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
858 
863 
867  // NOLINTNEXTLINE OSX does not compile with noexcept
869 
873  ~DoFAccessor() = default;
874 
886 
892  default;
893 
901  const DoFHandler<1, spacedim> &
903 
907  template <bool level_dof_access2>
908  void
909  copy_from(const DoFAccessor<0, 1, spacedim, level_dof_access2> &a);
910 
915  void
916  copy_from(const TriaAccessorBase<0, 1, spacedim> &da);
917 
930  TriaIterator<DoFAccessor<0, 1, spacedim, level_dof_access>>
931  child(const unsigned int c) const;
932 
939  typename ::internal::DoFHandlerImplementation::
940  Iterators<1, spacedim, level_dof_access>::line_iterator
941  line(const unsigned int i) const;
942 
949  typename ::internal::DoFHandlerImplementation::
950  Iterators<1, spacedim, level_dof_access>::quad_iterator
951  quad(const unsigned int i) const;
952 
996  void
998  std::vector<types::global_dof_index> &dof_indices,
999  const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1000 
1007  void
1009  const int level,
1010  std::vector<types::global_dof_index> &dof_indices,
1011  const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1012 
1033  const unsigned int vertex,
1034  const unsigned int i,
1035  const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1036 
1055  dof_index(const unsigned int i,
1056  const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1057 
1076  unsigned int
1078 
1086  unsigned int
1087  nth_active_fe_index(const unsigned int n) const;
1088 
1097  bool
1098  fe_index_is_active(const unsigned int fe_index) const;
1099 
1105  const FiniteElement<1, spacedim> &
1106  get_fe(const unsigned int fe_index) const;
1107 
1118  "This accessor object has not been "
1119  "associated with any DoFHandler object.");
1152 
1153 protected:
1157  DoFHandler<1, spacedim> *dof_handler;
1158 
1162  template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1163  bool
1164  operator==(
1165  const DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1166 
1170  template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1171  bool
1172  operator!=(
1173  const DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1174 
1178  void
1179  set_dof_handler(DoFHandler<1, spacedim> *dh);
1180 
1199  void
1201  const unsigned int i,
1202  const types::global_dof_index index,
1203  const unsigned int fe_index = AccessorData::invalid_fe_index) const;
1204 
1205  // Iterator classes need to be friends because they need to access
1206  // operator== and operator!=.
1207  template <typename>
1208  friend class TriaRawIterator;
1209 
1210 
1211  // Make the DoFHandler class a friend so that it can call the set_xxx()
1212  // functions.
1213  template <int, int>
1214  friend class DoFHandler;
1215 
1216  friend struct ::internal::DoFHandlerImplementation::Policy::
1217  Implementation;
1218  friend struct ::internal::DoFHandlerImplementation::Implementation;
1219  friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
1220  friend struct ::internal::DoFCellAccessorImplementation::Implementation;
1221 };
1222 
1223 
1224 
1225 /* -------------------------------------------------------------------------- */
1226 
1227 
1248 template <int structdim, int dim, int spacedim = dim>
1249 class DoFInvalidAccessor : public InvalidAccessor<structdim, dim, spacedim>
1250 {
1251 public:
1257 
1266  const int level = -1,
1267  const int index = -1,
1268  const AccessorData * local_data = 0);
1269 
1278 
1283  template <typename OtherAccessor>
1284  DoFInvalidAccessor(const OtherAccessor &);
1285 
1292  dof_index(const unsigned int i,
1293  const unsigned int fe_index =
1295 
1301  void
1302  set_dof_index(const unsigned int i,
1303  const types::global_dof_index index,
1304  const unsigned int fe_index =
1306 };
1307 
1308 
1309 
1310 /* -------------------------------------------------------------------------- */
1311 
1312 
1324 template <int dimension_, int space_dimension_, bool level_dof_access>
1325 class DoFCellAccessor : public DoFAccessor<dimension_,
1326  dimension_,
1327  space_dimension_,
1328  level_dof_access>
1329 {
1330 public:
1334  static const unsigned int dim = dimension_;
1335 
1339  static const unsigned int spacedim = space_dimension_;
1340 
1341 
1346 
1351  using BaseClass =
1353 
1358 
1363  using face_iterator = TriaIterator<DoFAccessor<dimension_ - 1,
1364  dimension_,
1365  space_dimension_,
1366  level_dof_access>>;
1367 
1379  const int level,
1380  const int index,
1381  const AccessorData *local_data);
1382 
1395  template <int structdim2, int dim2, int spacedim2>
1397 
1402  template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1405 
1411  default;
1412 
1416  DoFCellAccessor( // NOLINT
1418  &&) = default; // NOLINT
1419 
1423  ~DoFCellAccessor() = default;
1424 
1437  delete;
1438 
1443  operator=( // NOLINT
1445  &&) = default; // NOLINT
1446 
1461  parent() const;
1462 
1476  neighbor(const unsigned int i) const;
1477 
1484  periodic_neighbor(const unsigned int i) const;
1485 
1492  neighbor_or_periodic_neighbor(const unsigned int i) const;
1493 
1500  child(const unsigned int i) const;
1501 
1505  boost::container::small_vector<
1506  TriaIterator<
1510 
1518  face(const unsigned int i) const;
1519 
1523  boost::container::small_vector<face_iterator,
1526 
1534  neighbor_child_on_subface(const unsigned int face_no,
1535  const unsigned int subface_no) const;
1536 
1544  periodic_neighbor_child_on_subface(const unsigned int face_no,
1545  const unsigned int subface_no) const;
1546 
1575  template <class InputVector, typename number>
1576  void
1577  get_dof_values(const InputVector &values, Vector<number> &local_values) const;
1578 
1596  template <class InputVector, typename ForwardIterator>
1597  void
1598  get_dof_values(const InputVector &values,
1599  ForwardIterator local_values_begin,
1600  ForwardIterator local_values_end) const;
1601 
1622  template <class InputVector, typename ForwardIterator>
1623  void
1626  const InputVector & values,
1627  ForwardIterator local_values_begin,
1628  ForwardIterator local_values_end) const;
1629 
1654  template <class OutputVector, typename number>
1655  void
1656  set_dof_values(const Vector<number> &local_values,
1657  OutputVector & values) const;
1658 
1690  template <class InputVector, typename number>
1691  void
1692  get_interpolated_dof_values(
1693  const InputVector &values,
1694  Vector<number> & interpolated_values,
1695  const unsigned int fe_index =
1697 
1759  template <class OutputVector, typename number>
1760  void
1761  set_dof_values_by_interpolation(
1762  const Vector<number> &local_values,
1763  OutputVector & values,
1764  const unsigned int fe_index =
1766  const bool perform_check = false) const;
1767 
1777  template <class OutputVector, typename number>
1778  void
1779  distribute_local_to_global_by_interpolation(
1780  const Vector<number> &local_values,
1781  OutputVector & values,
1782  const unsigned int fe_index =
1784 
1799  template <typename number, typename OutputVector>
1800  void
1802  OutputVector & global_destination) const;
1803 
1818  template <typename ForwardIterator, typename OutputVector>
1819  void
1820  distribute_local_to_global(ForwardIterator local_source_begin,
1821  ForwardIterator local_source_end,
1822  OutputVector & global_destination) const;
1823 
1837  template <typename ForwardIterator, typename OutputVector>
1838  void
1841  ForwardIterator local_source_begin,
1842  ForwardIterator local_source_end,
1843  OutputVector & global_destination) const;
1844 
1851  template <typename number, typename OutputMatrix>
1852  void
1854  OutputMatrix &global_destination) const;
1855 
1860  template <typename number, typename OutputMatrix, typename OutputVector>
1861  void
1863  const Vector<number> & local_vector,
1864  OutputMatrix & global_matrix,
1865  OutputVector & global_vector) const;
1866 
1891  void
1893  std::vector<types::global_dof_index> &dof_indices) const;
1894 
1926  void
1927  get_dof_indices(std::vector<types::global_dof_index> &dof_indices) const;
1928 
1933  void
1934  get_mg_dof_indices(std::vector<types::global_dof_index> &dof_indices) const;
1935 
1961  get_fe() const;
1962 
1988  unsigned int
1990 
2017  void
2018  set_active_fe_index(const unsigned int i) const;
2027  void
2028  set_dof_indices(const std::vector<types::global_dof_index> &dof_indices);
2029 
2033  void
2034  set_mg_dof_indices(const std::vector<types::global_dof_index> &dof_indices);
2035 
2039  void
2040  update_cell_dof_indices_cache() const;
2041 
2068  get_future_fe() const;
2069 
2088  unsigned int
2090 
2099  void
2100  set_future_fe_index(const unsigned int i) const;
2101 
2108  bool
2110 
2119  void
2125 private:
2126  // Make the DoFHandler class a friend so that it can call the
2127  // update_cell_dof_indices_cache() function
2128  template <int, int>
2129  friend class DoFHandler;
2130  friend struct ::internal::DoFCellAccessorImplementation::Implementation;
2131 };
2132 
2133 
2134 template <int structdim, int dim, int spacedim, bool level_dof_access>
2135 inline bool
2137 {
2138  return level_dof_access;
2139 }
2140 
2141 
2142 
2143 template <int structdim, int dim, int spacedim>
2144 template <typename OtherAccessor>
2146  const OtherAccessor &)
2147 {
2148  Assert(false,
2149  ExcMessage("You are attempting an illegal conversion between "
2150  "iterator/accessor types. The constructor you call "
2151  "only exists to make certain template constructs "
2152  "easier to write as dimension independent code but "
2153  "the conversion is not valid in the current context."));
2154 }
2155 
2156 
2157 
2159 
2160 // include more templates
2161 #include "dof_accessor.templates.h"
2162 
2163 
2164 #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:224
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=(DoFAccessor< structdim, dim, spacedim, level_dof_access > &&)
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)
const DoFHandler< dim, spacedim > & get_dof_handler() 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
std::set< unsigned int > get_active_fe_indices() const
static constexpr unsigned int dimension
Definition: dof_accessor.h:218
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
typename ::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::quad_iterator quad(const unsigned int i) const
DoFAccessor(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &)=default
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 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)
bool fe_index_is_active(const unsigned int 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 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
const FiniteElement< dim, spacedim > & get_fe(const unsigned int 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
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:231
DoFHandler< dim, spacedim > * dof_handler
Definition: dof_accessor.h:660
DoFAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
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< 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
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
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 > > child(const unsigned int i) const
boost::container::small_vector< face_iterator, GeometryInfo< dimension_ >::faces_per_cell > face_iterators() 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_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)
DoFCellAccessor(const DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &)=default
const FiniteElement< dimension_, space_dimension_ > & get_fe() const
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
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > parent() const
unsigned int future_fe_index() const
bool future_fe_index_set() const
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices) 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:44
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: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:4607
#define DeclException0(Exception0)
Definition: exceptions.h:464
static ::ExceptionBase & ExcVectorDoesNotMatch()
#define Assert(cond, exc)
Definition: exceptions.h:1473
static ::ExceptionBase & ExcNotActive()
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:487
static ::ExceptionBase & ExcInvalidObject()
static ::ExceptionBase & ExcMatrixDoesNotMatch()
static ::ExceptionBase & ExcCantCompareIterators()
static ::ExceptionBase & ExcVectorNotEmpty()
static ::ExceptionBase & ExcMessage(std::string arg1)
Definition: hp.h:118
Definition: types.h:32
unsigned int global_dof_index
Definition: types.h:76
const ::Triangulation< dim, spacedim > & tria