Reference documentation for deal.II version 8.5.1
dof_accessor.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2017 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 at
12 // the top level of the deal.II distribution.
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 #include <deal.II/grid/tria_accessor.h>
22 #include <deal.II/dofs/dof_handler.h>
23 #include <deal.II/hp/dof_handler.h>
24 
25 #include <vector>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 template <typename number> class FullMatrix;
30 template <typename number> class Vector;
31 class ConstraintMatrix;
32 
33 template <typename Accessor> class TriaRawIterator;
34 
35 template <int, int> class FiniteElement;
36 
37 
38 namespace internal
39 {
40  namespace DoFCellAccessor
41  {
42  struct Implementation;
43  }
44 
45  namespace DoFHandler
46  {
47  struct Implementation;
48  namespace Policy
49  {
50  struct Implementation;
51  }
52  }
53 
54  namespace hp
55  {
56  namespace DoFHandler
57  {
58  struct Implementation;
59  }
60  }
61 }
62 
63 // note: the file dof_accessor.templates.h is included at the end of
64 // this file. this includes a lot of templates and thus makes
65 // compilation slower, but at the same time allows for more aggressive
66 // inlining and thus faster code.
67 
68 
69 namespace internal
70 {
71  namespace DoFAccessor
72  {
90  template <int structdim, int dim, int spacedim>
91  struct Inheritance
92  {
97  typedef ::TriaAccessor<structdim,dim,spacedim> BaseClass;
98  };
99 
100 
106  template <int dim, int spacedim>
107  struct Inheritance<dim,dim,spacedim>
108  {
113  typedef ::CellAccessor<dim,spacedim> BaseClass;
114  };
115  }
116 }
117 
118 
119 /* -------------------------------------------------------------------------- */
120 
121 
122 
200 template <int structdim, typename DoFHandlerType, bool level_dof_access>
201 class DoFAccessor : public ::internal::DoFAccessor::Inheritance<structdim, DoFHandlerType::dimension, DoFHandlerType::space_dimension>::BaseClass
202 {
203 public:
204 
209  static const unsigned int dimension=DoFHandlerType::dimension;
210 
215  static const unsigned int space_dimension=DoFHandlerType::space_dimension;
216 
221  typedef
222  typename ::internal::DoFAccessor::Inheritance<structdim, dimension, space_dimension>::BaseClass
224 
228  typedef DoFHandlerType AccessorData;
229 
240  DoFAccessor ();
241 
259  const int level,
260  const int index,
261  const DoFHandlerType *dof_handler);
262 
275  template <int structdim2, int dim2, int spacedim2>
277 
282  template <int dim2, class DoFHandlerType2, bool level_dof_access2>
284 
288  template <bool level_dof_access2>
297  const DoFHandlerType &
298  get_dof_handler () const;
299 
303  template <bool level_dof_access2>
305 
311 
316  static bool is_level_cell();
317 
329  child (const unsigned int c) const;
330 
336  typename ::internal::DoFHandler::Iterators<DoFHandlerType, level_dof_access>::line_iterator
337  line (const unsigned int i) const;
338 
344  typename ::internal::DoFHandler::Iterators<DoFHandlerType, level_dof_access>::quad_iterator
345  quad (const unsigned int i) const;
346 
392  void get_dof_indices (std::vector<types::global_dof_index> &dof_indices,
393  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
394 
401  void get_mg_dof_indices (const int level,
402  std::vector<types::global_dof_index> &dof_indices,
403  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
404 
408  void set_mg_dof_indices (const int level,
409  const std::vector<types::global_dof_index> &dof_indices,
410  const unsigned int fe_index = DoFHandlerType::default_fe_index);
411 
430  (const unsigned int vertex,
431  const unsigned int i,
432  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
433 
440  (const int level,
441  const unsigned int vertex,
442  const unsigned int i,
443  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
444 
473  (const unsigned int i,
474  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
475 
479  types::global_dof_index mg_dof_index (const int level, const unsigned int i) const;
480 
502  unsigned int
503  n_active_fe_indices () const;
504 
512  unsigned int
513  nth_active_fe_index (const unsigned int n) const;
514 
523  bool
524  fe_index_is_active (const unsigned int fe_index) const;
525 
532  get_fe (const unsigned int fe_index) const;
533 
576 
577 protected:
578 
582  DoFHandlerType *dof_handler;
583 public:
597  template <int dim2, class DoFHandlerType2, bool level_dof_access2>
599 
603  template <int dim2, class DoFHandlerType2, bool level_dof_access2>
605 protected:
609  void set_dof_handler (DoFHandlerType *dh);
610 
628  void set_dof_index
629  (const unsigned int i,
631  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
632 
633  void set_mg_dof_index (const int level, const unsigned int i, const types::global_dof_index index) const;
634 
653  (const unsigned int vertex,
654  const unsigned int i,
656  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
657 
658  void set_mg_vertex_dof_index
659  (const int level,
660  const unsigned int vertex,
661  const unsigned int i,
663  const unsigned int fe_index = DoFHandlerType::default_fe_index) const;
664 
669  template <typename> friend class TriaRawIterator;
670  template <int, class, bool> friend class DoFAccessor;
671 
672 private:
685 
690  template <int dim, int spacedim> friend class DoFHandler;
691  template <int dim, int spacedim> friend class hp::DoFHandler;
692 
693  friend struct ::internal::DoFHandler::Policy::Implementation;
694  friend struct ::internal::DoFHandler::Implementation;
695  friend struct ::internal::hp::DoFHandler::Implementation;
696  friend struct ::internal::DoFCellAccessor::Implementation;
697  friend struct ::internal::DoFAccessor::Implementation;
698 };
699 
700 
701 
711 template <template <int, int> class DoFHandlerType, int spacedim, bool level_dof_access>
712 class DoFAccessor<0,DoFHandlerType<1,spacedim>, level_dof_access> : public TriaAccessor<0,1,spacedim>
713 {
714 public:
715 
720  static const unsigned int dimension=1;
721 
726  static const unsigned int space_dimension=spacedim;
727 
733 
737  typedef DoFHandlerType<1,spacedim> AccessorData;
738 
749  DoFAccessor ();
750 
768  const typename TriaAccessor<0,1,spacedim>::VertexKind vertex_kind,
769  const unsigned int vertex_index,
770  const DoFHandlerType<1,spacedim> *dof_handler);
771 
778  const int = 0,
779  const int = 0,
780  const DoFHandlerType<1,spacedim> *dof_handler = 0);
781 
794  template <int structdim2, int dim2, int spacedim2>
796 
801  template <int dim2, class DoFHandlerType2, bool level_dof_access2>
803 
811  const DoFHandlerType<1,spacedim> &
812  get_dof_handler () const;
813 
817  DoFAccessor<0,DoFHandlerType<1,spacedim>, level_dof_access> &
818  operator = (const DoFAccessor<0,DoFHandlerType<1,spacedim>, level_dof_access> &da);
819 
823  template <bool level_dof_access2>
824  void copy_from (const DoFAccessor<0, DoFHandlerType<1,spacedim>, level_dof_access2> &a);
825 
831 
845  child (const unsigned int c) const;
846 
853  typename ::internal::DoFHandler::Iterators<DoFHandlerType<1,spacedim>, level_dof_access>::line_iterator
854  line (const unsigned int i) const;
855 
862  typename ::internal::DoFHandler::Iterators<DoFHandlerType<1,spacedim>, level_dof_access>::quad_iterator
863  quad (const unsigned int i) const;
864 
907  void get_dof_indices (std::vector<types::global_dof_index> &dof_indices,
908  const unsigned int fe_index = AccessorData::default_fe_index) const;
909 
916  void get_mg_dof_indices (const int level,
917  std::vector<types::global_dof_index> &dof_indices,
918  const unsigned int fe_index = AccessorData::default_fe_index) const;
919 
937  types::global_dof_index vertex_dof_index (const unsigned int vertex,
938  const unsigned int i,
939  const unsigned int fe_index = AccessorData::default_fe_index) const;
940 
955  types::global_dof_index dof_index (const unsigned int i,
956  const unsigned int fe_index = AccessorData::default_fe_index) const;
957 
976  unsigned int
977  n_active_fe_indices () const;
978 
986  unsigned int
987  nth_active_fe_index (const unsigned int n) const;
988 
997  bool
998  fe_index_is_active (const unsigned int fe_index) const;
999 
1005  const FiniteElement<DoFHandlerType<1,spacedim>::dimension,DoFHandlerType<1,spacedim>::space_dimension> &
1006  get_fe (const unsigned int fe_index) const;
1007 
1050 
1051 protected:
1052 
1056  DoFHandlerType<1,spacedim> *dof_handler;
1057 
1061  template <int dim2, class DoFHandlerType2, bool level_dof_access2>
1063 
1067  template <int dim2, class DoFHandlerType2, bool level_dof_access2>
1069 
1073  void set_dof_handler (DoFHandlerType<1,spacedim> *dh);
1074 
1092  void set_dof_index (const unsigned int i,
1094  const unsigned int fe_index = AccessorData::default_fe_index) const;
1095 
1113  void set_vertex_dof_index (const unsigned int vertex,
1114  const unsigned int i,
1116  const unsigned int fe_index = AccessorData::default_fe_index) const;
1117 
1122  template <typename> friend class TriaRawIterator;
1123 
1124 
1129  template <int, int> friend class DoFHandler;
1130  template <int, int> friend class hp::DoFHandler;
1131 
1132  friend struct ::internal::DoFHandler::Policy::Implementation;
1133  friend struct ::internal::DoFHandler::Implementation;
1134  friend struct ::internal::hp::DoFHandler::Implementation;
1135  friend struct ::internal::DoFCellAccessor::Implementation;
1136 };
1137 
1138 
1139 /* -------------------------------------------------------------------------- */
1140 
1141 
1154 template <typename DoFHandlerType, bool level_dof_access>
1155 class DoFCellAccessor : public DoFAccessor<DoFHandlerType::dimension,DoFHandlerType, level_dof_access>
1156 {
1157 public:
1161  static const unsigned int dim = DoFHandlerType::dimension;
1162 
1166  static const unsigned int spacedim = DoFHandlerType::space_dimension;
1167 
1168 
1172  typedef DoFHandlerType AccessorData;
1173 
1179 
1183  typedef DoFHandlerType Container;
1184 
1189  typedef
1190  TriaIterator<DoFAccessor<DoFHandlerType::dimension-1, DoFHandlerType, level_dof_access> >
1192 
1204  const int level,
1205  const int index,
1206  const AccessorData *local_data);
1207 
1220  template <int structdim2, int dim2, int spacedim2>
1222 
1227  template <int dim2, class DoFHandlerType2, bool level_dof_access2>
1228  explicit
1230 
1245  parent () const;
1246 
1260  neighbor (const unsigned int i) const;
1261 
1268  periodic_neighbor (const unsigned int i) const;
1269 
1276  neighbor_or_periodic_neighbor (const unsigned int i) const;
1277 
1284  child (const unsigned int i) const;
1285 
1293  face (const unsigned int i) const;
1294 
1302  neighbor_child_on_subface (const unsigned int face_no,
1303  const unsigned int subface_no) const;
1304 
1312  periodic_neighbor_child_on_subface (const unsigned int face_no,
1313  const unsigned int subface_no) const;
1314 
1340  template <class InputVector, typename number>
1341  void get_dof_values (const InputVector &values,
1342  Vector<number> &local_values) const;
1343 
1358  template <class InputVector, typename ForwardIterator>
1359  void get_dof_values (const InputVector &values,
1360  ForwardIterator local_values_begin,
1361  ForwardIterator local_values_end) const;
1362 
1379  template <class InputVector, typename ForwardIterator>
1380  void get_dof_values (const ConstraintMatrix &constraints,
1381  const InputVector &values,
1382  ForwardIterator local_values_begin,
1383  ForwardIterator local_values_end) const;
1384 
1406  template <class OutputVector, typename number>
1407  void set_dof_values (const Vector<number> &local_values,
1408  OutputVector &values) const;
1409 
1441  template <class InputVector, typename number>
1442  void get_interpolated_dof_values (const InputVector &values,
1443  Vector<number> &interpolated_values,
1444  const unsigned int fe_index
1445  = DoFHandlerType::default_fe_index) const;
1446 
1499  template <class OutputVector, typename number>
1500  void set_dof_values_by_interpolation (const Vector<number> &local_values,
1501  OutputVector &values,
1502  const unsigned int fe_index
1503  = DoFHandlerType::default_fe_index) const;
1504 
1513  template <typename number, typename OutputVector>
1514  void
1515  distribute_local_to_global (const Vector<number> &local_source,
1516  OutputVector &global_destination) const;
1517 
1526  template <typename ForwardIterator, typename OutputVector>
1527  void
1528  distribute_local_to_global (ForwardIterator local_source_begin,
1529  ForwardIterator local_source_end,
1530  OutputVector &global_destination) const;
1531 
1542  template <typename ForwardIterator, typename OutputVector>
1543  void
1544  distribute_local_to_global (const ConstraintMatrix &constraints,
1545  ForwardIterator local_source_begin,
1546  ForwardIterator local_source_end,
1547  OutputVector &global_destination) const;
1548 
1555  template <typename number, typename OutputMatrix>
1556  void
1557  distribute_local_to_global (const FullMatrix<number> &local_source,
1558  OutputMatrix &global_destination) const;
1559 
1564  template <typename number, typename OutputMatrix, typename OutputVector>
1565  void
1566  distribute_local_to_global (const FullMatrix<number> &local_matrix,
1567  const Vector<number> &local_vector,
1568  OutputMatrix &global_matrix,
1569  OutputVector &global_vector) const;
1570 
1595  void get_active_or_mg_dof_indices (std::vector<types::global_dof_index> &dof_indices) const;
1596 
1632  void get_dof_indices (std::vector<types::global_dof_index> &dof_indices) const;
1633 
1641  void get_mg_dof_indices (std::vector<types::global_dof_index> &dof_indices) const;
1642 
1668  get_fe () const;
1669 
1682  unsigned int active_fe_index () const;
1683 
1697  void set_active_fe_index (const unsigned int i);
1706  void set_dof_indices (const std::vector<types::global_dof_index> &dof_indices);
1707 
1711  void set_mg_dof_indices (const std::vector<types::global_dof_index> &dof_indices);
1712 
1716  void update_cell_dof_indices_cache () const;
1717 
1718 private:
1731 
1736  template <int dim, int spacedim> friend class DoFHandler;
1737  friend struct ::internal::DoFCellAccessor::Implementation;
1738 };
1739 
1740 
1741 template <int sd, typename DoFHandlerType, bool level_dof_access>
1742 inline
1743 bool
1745 {
1746  return level_dof_access;
1747 }
1748 
1749 
1750 
1751 DEAL_II_NAMESPACE_CLOSE
1752 
1753 // include more templates
1754 #include "dof_accessor.templates.h"
1755 
1756 
1757 #endif
DoFHandlerType * dof_handler
Definition: dof_accessor.h:582
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
::TriaAccessor< structdim, dim, spacedim > BaseClass
Definition: dof_accessor.h:97
static bool is_level_cell()
typename ::internal::DoFHandler::Iterators< DoFHandlerType, level_dof_access >::quad_iterator quad(const unsigned int i) const
face_iterator face(const unsigned int i) const
static const unsigned int spacedim
DoFHandlerType Container
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > parent() const
void copy_from(const DoFAccessor< structdim, DoFHandlerType, level_dof_access2 > &a)
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > neighbor(const unsigned int i) const
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > periodic_neighbor(const unsigned int i) const
void set_dof_values_by_interpolation(const Vector< number > &local_values, OutputVector &values, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
bool operator==(const DoFAccessor< dim2, DoFHandlerType2, level_dof_access2 > &) const
TriaIterator< DoFAccessor< DoFHandlerType::dimension-1, DoFHandlerType, level_dof_access > > face_iterator
void set_mg_dof_indices(const int level, const std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandlerType::default_fe_index)
void set_dof_handler(DoFHandlerType *dh)
int level() const
void get_mg_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
void get_interpolated_dof_values(const InputVector &values, Vector< number > &interpolated_values, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
unsigned int vertex_index(const unsigned int i) const
types::global_dof_index dof_index(const unsigned int i, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
DoFHandlerType AccessorData
Definition: dof_accessor.h:228
void update_cell_dof_indices_cache() const
Definition: dof_accessor.cc:46
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
unsigned int active_fe_index() const
void get_mg_dof_indices(const int level, std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
void set_active_fe_index(const unsigned int i)
static ::ExceptionBase & ExcCantCompareIterators()
typename ::internal::DoFAccessor::Inheritance< structdim, dimension, space_dimension >::BaseClass BaseClass
Definition: dof_accessor.h:223
typename ::internal::DoFHandler::Iterators< DoFHandlerType, level_dof_access >::line_iterator line(const unsigned int i) const
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > child(const unsigned int i) const
unsigned int global_dof_index
Definition: types.h:88
static ::ExceptionBase & ExcMatrixDoesNotMatch()
Point< spacedim > & vertex(const unsigned int i) const
void set_mg_dof_indices(const std::vector< types::global_dof_index > &dof_indices)
static const unsigned int dimension
Definition: dof_accessor.h:209
const Triangulation< dim, spacedim > * tria
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
Definition: dof_accessor.cc:77
#define DeclException0(Exception0)
Definition: exceptions.h:541
void set_dof_indices(const std::vector< types::global_dof_index > &dof_indices)
Definition: dof_accessor.cc:61
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > neighbor_or_periodic_neighbor(const unsigned int i) const
static ::ExceptionBase & ExcVectorNotEmpty()
int index() const
DoFCellAccessor< DoFHandlerType, level_dof_access > & operator=(const DoFCellAccessor< DoFHandlerType, level_dof_access > &da)
void get_dof_values(const InputVector &values, Vector< number > &local_values) const
Definition: hp.h:102
types::global_dof_index mg_vertex_dof_index(const int level, const unsigned int vertex, const unsigned int i, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
bool operator!=(const DoFAccessor< dim2, DoFHandlerType2, level_dof_access2 > &) const
void set_dof_values(const Vector< number > &local_values, OutputVector &values) const
unsigned int n_active_fe_indices() const
static const unsigned int dim
void distribute_local_to_global(const Vector< number > &local_source, OutputVector &global_destination) const
void set_dof_index(const unsigned int i, const types::global_dof_index index, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
DoFHandlerType AccessorData
static ::ExceptionBase & ExcNotActive()
types::global_dof_index vertex_dof_index(const unsigned int vertex, const unsigned int i, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
static ::ExceptionBase & ExcInvalidObject()
DoFCellAccessor(const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > *tria, const int level, const int index, const AccessorData *local_data)
static const unsigned int space_dimension
Definition: dof_accessor.h:215
unsigned int nth_active_fe_index(const unsigned int n) const
void get_active_or_mg_dof_indices(std::vector< types::global_dof_index > &dof_indices) 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=DoFHandlerType::default_fe_index) const
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > periodic_neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
Definition: dof_accessor.cc:89
DoFAccessor< structdim, DoFHandlerType, level_dof_access > & operator=(const DoFAccessor< structdim, DoFHandlerType, level_dof_access > &da)
const FiniteElement< DoFHandlerType::dimension, DoFHandlerType::space_dimension > & get_fe() const
const DoFHandlerType & get_dof_handler() const
const FiniteElement< DoFHandlerType::dimension, DoFHandlerType::space_dimension > & get_fe(const unsigned int fe_index) const
types::global_dof_index mg_dof_index(const int level, const unsigned int i) const
bool fe_index_is_active(const unsigned int fe_index) const
static ::ExceptionBase & ExcVectorDoesNotMatch()
TriaIterator< DoFAccessor< structdim, DoFHandlerType, level_dof_access > > child(const unsigned int c) const
DoFAccessor< DoFHandlerType::dimension, DoFHandlerType, level_dof_access > BaseClass