Reference documentation for deal.II version 9.1.0-pre
dof_accessor.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 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 
22 #include <deal.II/dofs/dof_handler.h>
23 
24 #include <deal.II/grid/tria_accessor.h>
25 
26 #include <deal.II/hp/dof_handler.h>
27 
28 #include <vector>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 template <typename number>
33 class FullMatrix;
34 template <typename number>
35 class Vector;
36 template <typename number>
38 
39 template <typename Accessor>
41 
42 template <int, int>
44 
45 
46 namespace internal
47 {
48  namespace DoFCellAccessorImplementation
49  {
50  struct Implementation;
51  }
52 
53  namespace DoFHandlerImplementation
54  {
55  struct Implementation;
56  namespace Policy
57  {
58  struct Implementation;
59  }
60  } // namespace DoFHandlerImplementation
61 
62  namespace hp
63  {
64  namespace DoFHandlerImplementation
65  {
66  struct Implementation;
67  }
68  } // namespace hp
69 } // namespace internal
70 
71 // note: the file dof_accessor.templates.h is included at the end of
72 // this file. this includes a lot of templates and thus makes
73 // compilation slower, but at the same time allows for more aggressive
74 // inlining and thus faster code.
75 
76 
77 namespace internal
78 {
79  namespace DoFAccessorImplementation
80  {
98  template <int structdim, int dim, int spacedim>
99  struct Inheritance
100  {
106  };
107 
108 
114  template <int dim, int spacedim>
115  struct Inheritance<dim, dim, spacedim>
116  {
122  };
123  } // namespace DoFAccessorImplementation
124 } // namespace internal
125 
126 
127 /* -------------------------------------------------------------------------- */
128 
129 
130 
208 template <int structdim, typename DoFHandlerType, bool level_dof_access>
211  structdim,
212  DoFHandlerType::dimension,
213  DoFHandlerType::space_dimension>::BaseClass
214 {
215 public:
220  static const unsigned int dimension = DoFHandlerType::dimension;
221 
226  static const unsigned int space_dimension = DoFHandlerType::space_dimension;
227 
232  using BaseClass = typename ::internal::DoFAccessorImplementation::
233  Inheritance<structdim, dimension, space_dimension>::BaseClass;
234 
238  using AccessorData = DoFHandlerType;
239 
250  DoFAccessor();
251 
268  DoFAccessor(const Triangulation<DoFHandlerType::dimension,
269  DoFHandlerType::space_dimension> *tria,
270  const int level,
271  const int index,
272  const DoFHandlerType *dof_handler);
273 
286  template <int structdim2, int dim2, int spacedim2>
288 
293  template <int dim2, class DoFHandlerType2, bool level_dof_access2>
295 
299  template <bool level_dof_access2>
300  DoFAccessor(
302 
314  &da) = delete;
315 
323  const DoFHandlerType &
324  get_dof_handler() const;
325 
329  template <bool level_dof_access2>
330  void
332 
337  void
338  copy_from(const TriaAccessorBase<structdim,
339  DoFHandlerType::dimension,
340  DoFHandlerType::space_dimension> &da);
341 
346  static bool
347  is_level_cell();
348 
360  child(const unsigned int c) const;
361 
367  typename ::internal::DoFHandlerImplementation::
368  Iterators<DoFHandlerType, level_dof_access>::line_iterator
369  line(const unsigned int i) const;
370 
376  typename ::internal::DoFHandlerImplementation::
377  Iterators<DoFHandlerType, level_dof_access>::quad_iterator
378  quad(const unsigned int i) const;
379 
425  void
426  get_dof_indices(
427  std::vector<types::global_dof_index> &dof_indices,
428  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
429 
436  void
437  get_mg_dof_indices(
438  const int level,
439  std::vector<types::global_dof_index> &dof_indices,
440  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
441 
445  void
446  set_mg_dof_indices(
447  const int level,
448  const std::vector<types::global_dof_index> &dof_indices,
449  const unsigned int fe_index = DoFHandlerType::default_fe_index);
450 
469  vertex_dof_index(
470  const unsigned int vertex,
471  const unsigned int i,
472  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
473 
480  mg_vertex_dof_index(
481  const int level,
482  const unsigned int vertex,
483  const unsigned int i,
484  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
485 
514  dof_index(
515  const unsigned int i,
516  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
517 
522  mg_dof_index(const int level, const unsigned int i) const;
523 
545  unsigned int
546  n_active_fe_indices() const;
547 
555  unsigned int
556  nth_active_fe_index(const unsigned int n) const;
557 
566  bool
567  fe_index_is_active(const unsigned int fe_index) const;
568 
574  const FiniteElement<DoFHandlerType::dimension,
575  DoFHandlerType::space_dimension> &
576  get_fe(const unsigned int fe_index) const;
577 
587  DeclExceptionMsg(ExcInvalidObject,
588  "This accessor object has not been "
589  "associated with any DoFHandler object.");
595  DeclException0(ExcVectorNotEmpty);
601  DeclException0(ExcVectorDoesNotMatch);
607  DeclException0(ExcMatrixDoesNotMatch);
615  DeclException0(ExcNotActive);
622 
623 protected:
627  DoFHandlerType *dof_handler;
628 
629 public:
643  template <int dim2, class DoFHandlerType2, bool level_dof_access2>
644  bool
645  operator==(
647 
651  template <int dim2, class DoFHandlerType2, bool level_dof_access2>
652  bool
653  operator!=(
655 
656 protected:
660  void
661  set_dof_handler(DoFHandlerType *dh);
662 
680  void
681  set_dof_index(
682  const unsigned int i,
683  const types::global_dof_index index,
684  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
685 
686  void
687  set_mg_dof_index(const int level,
688  const unsigned int i,
689  const types::global_dof_index index) const;
690 
708  void
709  set_vertex_dof_index(
710  const unsigned int vertex,
711  const unsigned int i,
712  const types::global_dof_index index,
713  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
714 
715  void
716  set_mg_vertex_dof_index(
717  const int level,
718  const unsigned int vertex,
719  const unsigned int i,
720  const types::global_dof_index index,
721  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
722 
727  template <typename>
728  friend class TriaRawIterator;
729  template <int, class, bool>
730  friend class DoFAccessor;
731 
732 private:
737  template <int dim, int spacedim>
738  friend class DoFHandler;
739  template <int dim, int spacedim>
740  friend class hp::DoFHandler;
741 
742  friend struct ::internal::DoFHandlerImplementation::Policy::
743  Implementation;
744  friend struct ::internal::DoFHandlerImplementation::Implementation;
745  friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
746  friend struct ::internal::DoFCellAccessorImplementation::Implementation;
747  friend struct ::internal::DoFAccessorImplementation::Implementation;
748 };
749 
750 
751 
761 template <template <int, int> class DoFHandlerType,
762  int spacedim,
763  bool level_dof_access>
764 class DoFAccessor<0, DoFHandlerType<1, spacedim>, level_dof_access>
765  : public TriaAccessor<0, 1, spacedim>
766 {
767 public:
772  static const unsigned int dimension = 1;
773 
778  static const unsigned int space_dimension = spacedim;
779 
785 
789  using AccessorData = DoFHandlerType<1, spacedim>;
790 
801  DoFAccessor();
802 
819  DoFAccessor(
820  const Triangulation<1, spacedim> * tria,
821  const typename TriaAccessor<0, 1, spacedim>::VertexKind vertex_kind,
822  const unsigned int vertex_index,
823  const DoFHandlerType<1, spacedim> * dof_handler);
824 
831  const int = 0,
832  const int = 0,
833  const DoFHandlerType<1, spacedim> *dof_handler = 0);
834 
847  template <int structdim2, int dim2, int spacedim2>
849 
854  template <int dim2, class DoFHandlerType2, bool level_dof_access2>
856 
866  DoFAccessor<0, DoFHandlerType<1, spacedim>, level_dof_access> &
867  operator=(const DoFAccessor<0, DoFHandlerType<1, spacedim>, level_dof_access>
868  &da) = delete;
869 
877  const DoFHandlerType<1, spacedim> &
878  get_dof_handler() const;
879 
883  template <bool level_dof_access2>
884  void
885  copy_from(
886  const DoFAccessor<0, DoFHandlerType<1, spacedim>, level_dof_access2> &a);
887 
892  void
893  copy_from(const TriaAccessorBase<0, 1, spacedim> &da);
894 
908  child(const unsigned int c) const;
909 
916  typename ::internal::DoFHandlerImplementation::
917  Iterators<DoFHandlerType<1, spacedim>, level_dof_access>::line_iterator
918  line(const unsigned int i) const;
919 
926  typename ::internal::DoFHandlerImplementation::
927  Iterators<DoFHandlerType<1, spacedim>, level_dof_access>::quad_iterator
928  quad(const unsigned int i) const;
929 
972  void
973  get_dof_indices(
974  std::vector<types::global_dof_index> &dof_indices,
975  const unsigned int fe_index = AccessorData::default_fe_index) const;
976 
983  void
984  get_mg_dof_indices(
985  const int level,
986  std::vector<types::global_dof_index> &dof_indices,
987  const unsigned int fe_index = AccessorData::default_fe_index) const;
988 
1007  vertex_dof_index(
1008  const unsigned int vertex,
1009  const unsigned int i,
1010  const unsigned int fe_index = AccessorData::default_fe_index) const;
1011 
1027  dof_index(const unsigned int i,
1028  const unsigned int fe_index = AccessorData::default_fe_index) const;
1029 
1048  unsigned int
1049  n_active_fe_indices() const;
1050 
1058  unsigned int
1059  nth_active_fe_index(const unsigned int n) const;
1060 
1069  bool
1070  fe_index_is_active(const unsigned int fe_index) const;
1071 
1078  DoFHandlerType<1, spacedim>::space_dimension> &
1079  get_fe(const unsigned int fe_index) const;
1080 
1090  DeclExceptionMsg(ExcInvalidObject,
1091  "This accessor object has not been "
1092  "associated with any DoFHandler object.");
1098  DeclException0(ExcVectorNotEmpty);
1104  DeclException0(ExcVectorDoesNotMatch);
1110  DeclException0(ExcMatrixDoesNotMatch);
1118  DeclException0(ExcNotActive);
1125 
1126 protected:
1130  DoFHandlerType<1, spacedim> *dof_handler;
1131 
1135  template <int dim2, class DoFHandlerType2, bool level_dof_access2>
1136  bool
1137  operator==(
1139 
1143  template <int dim2, class DoFHandlerType2, bool level_dof_access2>
1144  bool
1145  operator!=(
1147 
1151  void set_dof_handler(DoFHandlerType<1, spacedim> *dh);
1152 
1170  void
1171  set_dof_index(
1172  const unsigned int i,
1173  const types::global_dof_index index,
1174  const unsigned int fe_index = AccessorData::default_fe_index) const;
1175 
1193  void
1194  set_vertex_dof_index(
1195  const unsigned int vertex,
1196  const unsigned int i,
1197  const types::global_dof_index index,
1198  const unsigned int fe_index = AccessorData::default_fe_index) const;
1199 
1204  template <typename>
1205  friend class TriaRawIterator;
1206 
1207 
1212  template <int, int>
1213  friend class DoFHandler;
1214  template <int, int>
1215  friend class hp::DoFHandler;
1216 
1217  friend struct ::internal::DoFHandlerImplementation::Policy::
1218  Implementation;
1219  friend struct ::internal::DoFHandlerImplementation::Implementation;
1220  friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
1221  friend struct ::internal::DoFCellAccessorImplementation::Implementation;
1222 };
1223 
1224 
1225 
1226 /* -------------------------------------------------------------------------- */
1227 
1228 
1250 template <int structdim, int dim, int spacedim = dim>
1251 class DoFInvalidAccessor : public InvalidAccessor<structdim, dim, spacedim>
1252 {
1253 public:
1257  using AccessorData =
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 
1293  void
1294  set_dof_index(const unsigned int i,
1295  const types::global_dof_index index,
1296  const unsigned int fe_index =
1298 };
1299 
1300 
1301 
1302 /* -------------------------------------------------------------------------- */
1303 
1304 
1317 template <typename DoFHandlerType, bool level_dof_access>
1318 class DoFCellAccessor : public DoFAccessor<DoFHandlerType::dimension,
1319  DoFHandlerType,
1320  level_dof_access>
1321 {
1322 public:
1326  static const unsigned int dim = DoFHandlerType::dimension;
1327 
1331  static const unsigned int spacedim = DoFHandlerType::space_dimension;
1332 
1333 
1337  using AccessorData = DoFHandlerType;
1338 
1343  using BaseClass =
1345 
1349  using Container = DoFHandlerType;
1350 
1355  using face_iterator = TriaIterator<DoFAccessor<DoFHandlerType::dimension - 1,
1356  DoFHandlerType,
1357  level_dof_access>>;
1358 
1369  DoFCellAccessor(const Triangulation<DoFHandlerType::dimension,
1370  DoFHandlerType::space_dimension> *tria,
1371  const int level,
1372  const int index,
1373  const AccessorData *local_data);
1374 
1387  template <int structdim2, int dim2, int spacedim2>
1389 
1394  template <int dim2, class DoFHandlerType2, bool level_dof_access2>
1395  explicit DoFCellAccessor(
1397 
1409  delete;
1410 
1425  parent() const;
1426 
1440  neighbor(const unsigned int i) const;
1441 
1448  periodic_neighbor(const unsigned int i) const;
1449 
1456  neighbor_or_periodic_neighbor(const unsigned int i) const;
1457 
1464  child(const unsigned int i) const;
1465 
1473  face(const unsigned int i) const;
1474 
1482  neighbor_child_on_subface(const unsigned int face_no,
1483  const unsigned int subface_no) const;
1484 
1492  periodic_neighbor_child_on_subface(const unsigned int face_no,
1493  const unsigned int subface_no) const;
1494 
1520  template <class InputVector, typename number>
1521  void
1522  get_dof_values(const InputVector &values, Vector<number> &local_values) const;
1523 
1538  template <class InputVector, typename ForwardIterator>
1539  void
1540  get_dof_values(const InputVector &values,
1541  ForwardIterator local_values_begin,
1542  ForwardIterator local_values_end) const;
1543 
1561  template <class InputVector, typename ForwardIterator>
1562  void
1563  get_dof_values(
1565  const InputVector & values,
1566  ForwardIterator local_values_begin,
1567  ForwardIterator local_values_end) const;
1568 
1590  template <class OutputVector, typename number>
1591  void
1592  set_dof_values(const Vector<number> &local_values,
1593  OutputVector & values) const;
1594 
1626  template <class InputVector, typename number>
1627  void
1628  get_interpolated_dof_values(
1629  const InputVector &values,
1630  Vector<number> & interpolated_values,
1631  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
1632 
1685  template <class OutputVector, typename number>
1686  void
1687  set_dof_values_by_interpolation(
1688  const Vector<number> &local_values,
1689  OutputVector & values,
1690  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
1691 
1700  template <typename number, typename OutputVector>
1701  void
1702  distribute_local_to_global(const Vector<number> &local_source,
1703  OutputVector & global_destination) const;
1704 
1713  template <typename ForwardIterator, typename OutputVector>
1714  void
1715  distribute_local_to_global(ForwardIterator local_source_begin,
1716  ForwardIterator local_source_end,
1717  OutputVector & global_destination) const;
1718 
1729  template <typename ForwardIterator, typename OutputVector>
1730  void
1731  distribute_local_to_global(
1733  ForwardIterator local_source_begin,
1734  ForwardIterator local_source_end,
1735  OutputVector & global_destination) const;
1736 
1743  template <typename number, typename OutputMatrix>
1744  void
1745  distribute_local_to_global(const FullMatrix<number> &local_source,
1746  OutputMatrix &global_destination) const;
1747 
1752  template <typename number, typename OutputMatrix, typename OutputVector>
1753  void
1754  distribute_local_to_global(const FullMatrix<number> &local_matrix,
1755  const Vector<number> & local_vector,
1756  OutputMatrix & global_matrix,
1757  OutputVector & global_vector) const;
1758 
1783  void
1784  get_active_or_mg_dof_indices(
1785  std::vector<types::global_dof_index> &dof_indices) const;
1786 
1822  void
1823  get_dof_indices(std::vector<types::global_dof_index> &dof_indices) const;
1824 
1832  void
1833  get_mg_dof_indices(std::vector<types::global_dof_index> &dof_indices) const;
1834 
1859  const FiniteElement<DoFHandlerType::dimension,
1860  DoFHandlerType::space_dimension> &
1861  get_fe() const;
1862 
1884  unsigned int
1885  active_fe_index() const;
1886 
1914  void
1915  set_active_fe_index(const unsigned int i) const;
1924  void
1925  set_dof_indices(const std::vector<types::global_dof_index> &dof_indices);
1926 
1930  void
1931  set_mg_dof_indices(const std::vector<types::global_dof_index> &dof_indices);
1932 
1936  void
1937  update_cell_dof_indices_cache() const;
1938 
1939 private:
1944  template <int dim, int spacedim>
1945  friend class DoFHandler;
1946  friend struct ::internal::DoFCellAccessorImplementation::Implementation;
1947 };
1948 
1949 
1950 template <int sd, typename DoFHandlerType, bool level_dof_access>
1951 inline bool
1953 {
1954  return level_dof_access;
1955 }
1956 
1957 
1958 
1959 template <int structdim, int dim, int spacedim>
1960 template <typename OtherAccessor>
1962  const OtherAccessor &)
1963 {
1964  Assert(false,
1965  ExcMessage("You are attempting an illegal conversion between "
1966  "iterator/accessor types. The constructor you call "
1967  "only exists to make certain template constructs "
1968  "easier to write as dimension independent code but "
1969  "the conversion is not valid in the current context."));
1970 }
1971 
1972 
1973 
1974 DEAL_II_NAMESPACE_CLOSE
1975 
1976 // include more templates
1977 #include "dof_accessor.templates.h"
1978 
1979 
1980 #endif
DoFHandlerType * dof_handler
Definition: dof_accessor.h:627
static bool is_level_cell()
typename::internal::DoFAccessorImplementation::Inheritance< structdim, dimension, space_dimension >::BaseClass BaseClass
Definition: dof_accessor.h:233
typename InvalidAccessor< structdim, dim, spacedim >::AccessorData AccessorData
typename TriaAccessorBase< structdim, dim, spacedim >::AccessorData AccessorData
unsigned long long int global_dof_index
Definition: types.h:72
static::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1227
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:397
#define DeclException0(Exception0)
Definition: exceptions.h:385
Definition: hp.h:102
DoFHandlerType Container
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
DoFHandlerType AccessorData
static::ExceptionBase & ExcCantCompareIterators()