Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
dof_accessor.cc
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 #include <deal.II/dofs/dof_accessor.h>
17 #include <deal.II/dofs/dof_handler.h>
18 #include <deal.II/dofs/dof_levels.h>
19 
20 #include <deal.II/fe/fe.h>
21 
22 #include <deal.II/grid/tria_iterator.h>
23 #include <deal.II/grid/tria_iterator.templates.h>
24 
25 #include <deal.II/hp/dof_handler.h>
26 
27 #include <vector>
28 
29 DEAL_II_NAMESPACE_OPEN
30 
31 /* --------------------- Static variables: DoFAccessor --------------------- */
32 
33 template <int structdim, typename DoFHandlerType, bool level_dof_access>
34 const unsigned int
36 
37 template <int structdim, typename DoFHandlerType, bool level_dof_access>
38 const unsigned int
40 
41 
42 
43 /* ---------------------- Functions: DoFInvalidAccessor -------------------- */
44 
45 template <int structdim, int dim, int spacedim>
48  const int,
49  const int,
50  const AccessorData *)
51 {
52  Assert(false,
53  ExcMessage("You are attempting an illegal conversion between "
54  "iterator/accessor types. The constructor you call "
55  "only exists to make certain template constructs "
56  "easier to write as dimension independent code but "
57  "the conversion is not valid in the current context."));
58 }
59 
60 
61 
62 template <int structdim, int dim, int spacedim>
64  const DoFInvalidAccessor &i)
65  : InvalidAccessor<structdim, dim, spacedim>(
66  static_cast<const InvalidAccessor<structdim, dim, spacedim> &>(i))
67 {
68  Assert(false,
69  ExcMessage("You are attempting an illegal conversion between "
70  "iterator/accessor types. The constructor you call "
71  "only exists to make certain template constructs "
72  "easier to write as dimension independent code but "
73  "the conversion is not valid in the current context."));
74 }
75 
76 
77 
78 template <int structdim, int dim, int spacedim>
79 void
81  const unsigned int,
83  const unsigned int) const
84 {
85  Assert(false, ExcInternalError());
86 }
87 
88 
89 
90 /*------------------------- Functions: DoFCellAccessor -----------------------*/
91 
92 
93 
94 template <typename DoFHandlerType, bool lda>
95 void
97 {
98  Assert(static_cast<unsigned int>(this->present_level) <
99  this->dof_handler->levels.size(),
100  ExcMessage("DoFHandler not initialized"));
101 
102  Assert(this->dof_handler != nullptr, typename BaseClass::ExcInvalidObject());
103 
104  internal::DoFCellAccessorImplementation::Implementation::
105  update_cell_dof_indices_cache(*this);
106 }
107 
108 
109 
110 template <typename DoFHandlerType, bool lda>
111 void
113  const std::vector<types::global_dof_index> &local_dof_indices)
114 {
115  Assert(static_cast<unsigned int>(this->present_level) <
116  this->dof_handler->levels.size(),
117  ExcMessage("DoFHandler not initialized"));
118 
119  Assert(this->dof_handler != nullptr, typename BaseClass::ExcInvalidObject());
120 
121  internal::DoFCellAccessorImplementation::Implementation::set_dof_indices(
122  *this, local_dof_indices);
123 }
124 
125 
126 
127 template <typename DoFHandlerType, bool lda>
130  const unsigned int face,
131  const unsigned int subface) const
132 {
136  this->dof_handler);
137 }
138 
139 
140 
141 template <typename DoFHandlerType, bool lda>
144  const unsigned int face,
145  const unsigned int subface) const
146 {
149  subface);
151  this->dof_handler);
152 }
153 
154 
155 
156 template <typename DoFHandlerType, bool lda>
159  const unsigned int face) const
160 {
164  this->dof_handler);
165 }
166 
167 
168 
169 template <typename DoFHandlerType, bool lda>
172  const unsigned int face) const
173 {
177  this->dof_handler);
178 }
179 
180 
181 
182 // --------------------------------------------------------------------------
183 // explicit instantiations
184 #include "dof_accessor.inst"
185 
186 DEAL_II_NAMESPACE_CLOSE
void set_dof_index(const unsigned int i, const types::global_dof_index index, const unsigned int fe_index=DoFHandler< dim, spacedim >::default_fe_index) const
Definition: dof_accessor.cc:80
TriaIterator< CellAccessor< dim, spacedim > > neighbor_or_periodic_neighbor(const unsigned int i) const
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > periodic_neighbor(const unsigned int i) const
typename InvalidAccessor< structdim, dim, spacedim >::AccessorData AccessorData
TriaIterator< CellAccessor< dim, spacedim > > neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
void update_cell_dof_indices_cache() const
Definition: dof_accessor.cc:96
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1407
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
void set_dof_indices(const std::vector< types::global_dof_index > &dof_indices)
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > neighbor_or_periodic_neighbor(const unsigned int i) const
TriaIterator< CellAccessor< dim, spacedim > > periodic_neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
unsigned int global_dof_index
Definition: types.h:89
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
TriaIterator< CellAccessor< dim, spacedim > > periodic_neighbor(const unsigned int i) const
TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > periodic_neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
static ::ExceptionBase & ExcInternalError()