deal.II version 9.7.0
\(\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\}}\)
Loading...
Searching...
No Matches
dof_accessor.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1998 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_dof_accessor_h
16#define dealii_dof_accessor_h
17
18
19#include <deal.II/base/config.h>
20
21#include <deal.II/base/types.h>
22
26
28
32
33#include <boost/container/small_vector.hpp>
34
35#include <limits>
36#include <set>
37#include <type_traits>
38#include <vector>
39
40
42
43// Forward declarations
44#ifndef DOXYGEN
45template <typename number>
46class FullMatrix;
47template <typename number>
48class Vector;
49template <typename number>
51
52template <typename Accessor>
53class TriaRawIterator;
54
55template <int, int>
56class FiniteElement;
57
58namespace internal
59{
61 {
62 struct Implementation;
63 }
64
66 {
67 struct Implementation;
68 namespace Policy
69 {
70 struct Implementation;
71 }
72 } // namespace DoFHandlerImplementation
73
74 namespace hp
75 {
77 {
78 struct Implementation;
79 }
80 } // namespace hp
81} // namespace internal
82#endif
83
84
85namespace internal
86{
88 {
105 template <int structdim, int dim, int spacedim>
114
115
121 template <int dim, int spacedim>
122 struct Inheritance<dim, dim, spacedim>
123 {
129 };
130
131 struct Implementation;
132 } // namespace DoFAccessorImplementation
133} // namespace internal
134
135
136/* -------------------------------------------------------------------------- */
137
138
139
212template <int structdim, int dim, int spacedim, bool level_dof_access>
214 Inheritance<structdim, dim, spacedim>::BaseClass
215{
216public:
221 static constexpr unsigned int dimension = dim;
222
227 static constexpr unsigned int space_dimension = spacedim;
228
233 using BaseClass = typename ::internal::DoFAccessorImplementation::
234 Inheritance<structdim, dimension, space_dimension>::BaseClass;
235
240
247
252
269 const int level,
270 const int index,
272
277 default;
278
282 DoFAccessor( // NOLINT
284 default; // NOLINT
285
289 ~DoFAccessor() = default;
290
303 template <int structdim2, int dim2, int spacedim2>
305
310 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
313
317 template <bool level_dof_access2>
319
331 delete;
332
337 operator=( // NOLINT
339 default; // NOLINT
340
344
350
354 template <bool level_dof_access2>
355 void
357
362 void
364
369 static bool
371
378
383 child(const unsigned int c) const;
384
390 typename ::internal::DoFHandlerImplementation::
391 Iterators<dim, spacedim, level_dof_access>::line_iterator
392 line(const unsigned int i) const;
393
399 typename ::internal::DoFHandlerImplementation::
400 Iterators<dim, spacedim, level_dof_access>::quad_iterator
401 quad(const unsigned int i) const;
402
406
413
448 void
450 std::vector<types::global_dof_index> &dof_indices,
451 const types::fe_index fe_index = numbers::invalid_fe_index) const;
452
459 void
461 const int level,
462 std::vector<types::global_dof_index> &dof_indices,
463 const types::fe_index fe_index = numbers::invalid_fe_index) const;
464
468 void
470 const int level,
471 const std::vector<types::global_dof_index> &dof_indices,
473
497 const unsigned int vertex,
498 const unsigned int i,
499 const types::fe_index fe_index = numbers::invalid_fe_index) const;
500
508 const int level,
509 const unsigned int vertex,
510 const unsigned int i,
511 const types::fe_index fe_index = numbers::invalid_fe_index) const;
512
541 dof_index(const unsigned int i,
542 const types::fe_index fe_index = numbers::invalid_fe_index) const;
543
548 mg_dof_index(const int level, const unsigned int i) const;
549
553
560
571 unsigned int
573
582 nth_active_fe_index(const unsigned int n) const;
583
590 std::set<types::fe_index>
592
602 bool
603 fe_index_is_active(const types::fe_index fe_index) const;
604
611 get_fe(const types::fe_index fe_index) const;
612
616
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
718 const unsigned int i,
719 const types::global_dof_index index,
720 const types::fe_index fe_index = numbers::invalid_fe_index) const;
721
722 void
723 set_mg_dof_index(const int level,
724 const unsigned int i,
725 const types::global_dof_index index) const;
726
727 void
729 const int level,
730 const unsigned int vertex,
731 const unsigned int i,
732 const types::global_dof_index index,
733 const types::fe_index fe_index = numbers::invalid_fe_index) const;
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 friend class DoFHandler<dim, spacedim>;
746
747 friend struct ::internal::DoFHandlerImplementation::Policy::
748 Implementation;
749 friend struct ::internal::DoFHandlerImplementation::Implementation;
750 friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
751 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
752 friend struct ::internal::DoFAccessorImplementation::Implementation;
753};
754
755
756
764template <int spacedim, bool level_dof_access>
765class DoFAccessor<0, 1, spacedim, level_dof_access>
766 : public TriaAccessor<0, 1, spacedim>
767{
768public:
773 static constexpr unsigned int dimension = 1;
774
779 static constexpr unsigned int space_dimension = spacedim;
780
786
791
798
802 DoFAccessor();
803
822 const typename TriaAccessor<0, 1, spacedim>::VertexKind vertex_kind,
823 const unsigned int vertex_index,
825
833 const int = 0,
834 const int = 0,
836
849 template <int structdim2, int dim2, int spacedim2>
851
856 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
859
864
868 // NOLINTNEXTLINE OSX does not compile with noexcept
870
874 ~DoFAccessor() = default;
875
887
893 default;
894
898
902 const DoFHandler<1, spacedim> &
903 get_dof_handler() const;
904
908 template <bool level_dof_access2>
909 void
910 copy_from(const DoFAccessor<0, 1, spacedim, level_dof_access2> &a);
911
916 void
917 copy_from(const TriaAccessorBase<0, 1, spacedim> &da);
918
925
931 TriaIterator<DoFAccessor<0, 1, spacedim, level_dof_access>>
932 child(const unsigned int c) const;
933
940 typename ::internal::DoFHandlerImplementation::
941 Iterators<1, spacedim, level_dof_access>::line_iterator
942 line(const unsigned int i) const;
943
950 typename ::internal::DoFHandlerImplementation::
951 Iterators<1, spacedim, level_dof_access>::quad_iterator
952 quad(const unsigned int i) const;
953
957
964
997 void
999 std::vector<types::global_dof_index> &dof_indices,
1000 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1001
1008 void
1010 const int level,
1011 std::vector<types::global_dof_index> &dof_indices,
1012 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1013
1032 types::global_dof_index
1034 const unsigned int vertex,
1035 const unsigned int i,
1036 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1037
1055 types::global_dof_index
1056 dof_index(const unsigned int i,
1057 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1058
1062
1069
1077 unsigned int
1078 n_active_fe_indices() const;
1079
1087 types::fe_index
1088 nth_active_fe_index(const unsigned int n) const;
1089
1098 bool
1099 fe_index_is_active(const types::fe_index fe_index) const;
1100
1106 const FiniteElement<1, spacedim> &
1107 get_fe(const types::fe_index fe_index) const;
1108
1112
1119 "This accessor object has not been "
1120 "associated with any DoFHandler object.");
1153
1154protected:
1159
1163 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1164 bool
1165 operator==(
1166 const DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1167
1171 template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
1172 bool
1173 operator!=(
1174 const DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1175
1179 void
1180 set_dof_handler(DoFHandler<1, spacedim> *dh);
1181
1200 void
1202 const unsigned int i,
1203 const types::global_dof_index index,
1204 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1205
1206 // Iterator classes need to be friends because they need to access
1207 // operator== and operator!=.
1208 template <typename>
1209 friend class TriaRawIterator;
1210
1211
1212 // Make the DoFHandler class a friend so that it can call the set_xxx()
1213 // functions.
1214 friend class DoFHandler<1, spacedim>;
1215
1216 friend struct ::internal::DoFHandlerImplementation::Policy::
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
1248template <int structdim, int dim, int spacedim = dim>
1249class DoFInvalidAccessor : public InvalidAccessor<structdim, dim, spacedim>
1250{
1251public:
1257
1265 DoFInvalidAccessor(const void *parent = nullptr,
1266 const int level = -1,
1267 const int index = -1,
1268 const AccessorData *local_data = nullptr);
1269
1278
1283 template <typename OtherAccessor>
1284 DoFInvalidAccessor(const OtherAccessor &);
1285
1292 dof_index(const unsigned int i,
1293 const types::fe_index fe_index =
1295
1301 void
1303 const unsigned int i,
1305 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1306};
1307
1308
1309
1310/* -------------------------------------------------------------------------- */
1311
1312
1324template <int dimension_, int space_dimension_, bool level_dof_access>
1325class DoFCellAccessor : public DoFAccessor<dimension_,
1326 dimension_,
1327 space_dimension_,
1328 level_dof_access>
1329{
1330public:
1334 static const unsigned int dim = dimension_;
1335
1339 static const unsigned int spacedim = space_dimension_;
1340
1341
1346
1353
1358
1364 dimension_,
1365 space_dimension_,
1366 level_dof_access>>;
1367
1374
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
1418 &&) = default; // NOLINT
1419
1423 ~DoFCellAccessor() = default;
1424
1437 delete;
1438
1443 operator=( // NOLINT
1445 &&) = default; // NOLINT
1446
1450
1461 parent() const;
1462
1469
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<
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
1550
1557
1575 template <class InputVector, typename number>
1576 void
1577 get_dof_values(const InputVector &values, Vector<number> &local_values) const;
1578
1596 template <typename Number, typename ForwardIterator>
1597 void
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 <typename Number>
1691 void
1693 const ReadVector<Number> &values,
1694 Vector<Number> &interpolated_values,
1695 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1696
1758 template <class OutputVector, typename number>
1759 void
1761 const Vector<number> &local_values,
1762 OutputVector &values,
1764 const bool perform_check = false) const;
1765
1775 template <class OutputVector, typename number>
1776 void
1778 const Vector<number> &local_values,
1779 OutputVector &values,
1780 const types::fe_index fe_index = numbers::invalid_fe_index) const;
1781
1796 template <typename number, typename OutputVector>
1797 void
1799 OutputVector &global_destination) const;
1800
1815 template <typename ForwardIterator, typename OutputVector>
1816 void
1817 distribute_local_to_global(ForwardIterator local_source_begin,
1818 ForwardIterator local_source_end,
1819 OutputVector &global_destination) const;
1820
1834 template <typename ForwardIterator, typename OutputVector>
1835 void
1838 ForwardIterator local_source_begin,
1839 ForwardIterator local_source_end,
1840 OutputVector &global_destination) const;
1841
1848 template <typename number, typename OutputMatrix>
1849 void
1851 OutputMatrix &global_destination) const;
1852
1857 template <typename number, typename OutputMatrix, typename OutputVector>
1858 void
1860 const Vector<number> &local_vector,
1861 OutputMatrix &global_matrix,
1862 OutputVector &global_vector) const;
1863
1867
1871
1875
1888 void
1890 std::vector<types::global_dof_index> &dof_indices) const;
1891
1923 void
1924 get_dof_indices(std::vector<types::global_dof_index> &dof_indices) const;
1925
1930 void
1931 get_mg_dof_indices(std::vector<types::global_dof_index> &dof_indices) const;
1932
1936
1943
1958 get_fe() const;
1959
1987
2014 void
2019
2024 void
2025 set_dof_indices(const std::vector<types::global_dof_index> &dof_indices);
2026
2030 void
2031 set_mg_dof_indices(const std::vector<types::global_dof_index> &dof_indices);
2032
2039
2060
2081
2090 void
2092
2099 bool
2101
2110 void
2115
2116private:
2117 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
2118};
2119
2120
2121template <int structdim, int dim, int spacedim, bool level_dof_access>
2122inline bool
2127
2128
2129
2130template <int structdim, int dim, int spacedim>
2131template <typename OtherAccessor>
2133 const OtherAccessor &)
2134{
2135 Assert(false,
2136 ExcMessage("You are attempting an illegal conversion between "
2137 "iterator/accessor types. The constructor you call "
2138 "only exists to make certain template constructs "
2139 "easier to write as dimension independent code but "
2140 "the conversion is not valid in the current context."));
2141}
2142
2143
2144
2145/*------------------------- Functions: DoFAccessor ---------------------------*/
2146
2147
2148template <int structdim, int dim, int spacedim, bool level_dof_access>
2153
2154
2155
2156template <int structdim, int dim, int spacedim, bool level_dof_access>
2158 const Triangulation<dim, spacedim> *tria,
2159 const int level,
2160 const int index,
2162 : ::internal::DoFAccessorImplementation::
2163 Inheritance<structdim, dim, spacedim>::BaseClass(tria, level, index)
2164 , dof_handler(const_cast<DoFHandler<dim, spacedim> *>(dof_handler))
2165{
2166 Assert(
2167 tria == nullptr || &dof_handler->get_triangulation() == tria,
2168 ExcMessage(
2169 "You can't create a DoF accessor in which the DoFHandler object "
2170 "uses a different triangulation than the one you pass as argument."));
2171}
2172
2173
2174
2175template <int structdim, int dim, int spacedim, bool level_dof_access>
2176template <int structdim2, int dim2, int spacedim2>
2182
2183
2184
2185template <int structdim, int dim, int spacedim, bool level_dof_access>
2186template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
2189 : BaseClass(other)
2190 , dof_handler(nullptr)
2191{
2192 Assert(false,
2193 ExcMessage(
2194 "You are trying to assign iterators that are incompatible. "
2195 "The reason for incompatibility is that they refer to objects of "
2196 "different dimensionality (e.g., assigning a line iterator "
2197 "to a quad iterator)."));
2198}
2199
2200
2201
2202template <int structdim, int dim, int spacedim, bool level_dof_access>
2203template <bool level_dof_access2>
2209
2210
2211
2212template <int structdim, int dim, int spacedim, bool level_dof_access>
2213inline void
2220
2221
2222
2223template <int structdim, int dim, int spacedim, bool level_dof_access>
2224inline const DoFHandler<dim, spacedim> &
2230
2231
2232
2233template <int structdim, int dim, int spacedim, bool level_dof_access>
2234inline void
2241
2242
2243
2244template <int structdim, int dim, int spacedim, bool level_dof_access>
2245template <bool level_dof_access2>
2246inline void
2253
2254
2255
2256template <int structdim, int dim, int spacedim, bool level_dof_access>
2257template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
2258inline bool
2266
2267
2268
2269template <int structdim, int dim, int spacedim, bool level_dof_access>
2270template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
2271inline bool
2279
2280
2281
2282template <int structdim, int dim, int spacedim, bool level_dof_access>
2285 const unsigned int i) const
2286{
2287 Assert(static_cast<unsigned int>(this->present_level) <
2288 this->dof_handler->object_dof_indices.size(),
2289 ExcMessage("DoFHandler not initialized"));
2290
2293
2295 *t, this->dof_handler);
2296 return q;
2297}
2298
2299
2300namespace internal
2301{
2302 namespace DoFAccessorImplementation
2303 {
2308 template <int structdim, int dim, int spacedim, bool level_dof_access>
2312 const types::fe_index fe_index)
2313 {
2314 if (cell.get_dof_handler().has_hp_capabilities() == false)
2315 {
2316 // No hp enabled, and the argument is at its default value -> we
2317 // can translate to the default active fe index
2318 Assert(
2319 (fe_index == numbers::invalid_fe_index) ||
2321 ExcMessage(
2322 "It is not possible to specify a FE index if no hp support is used!"));
2323
2325 }
2326 else
2327 {
2328 // Otherwise: If anything other than the default is provided by
2329 // the caller, then we should take just that. As an exception, if
2330 // we are on a cell (rather than a face/edge/vertex), then we know
2331 // that there is only one active fe index on this cell and we can
2332 // use that:
2333 if ((dim == structdim) && (fe_index == numbers::invalid_fe_index))
2334 {
2336
2337 return cell.nth_active_fe_index(0);
2338 }
2339
2340 Assert((fe_index != numbers::invalid_fe_index),
2341 ExcMessage(
2342 "You need to specify a FE index if hp support is used!"));
2343
2344 return fe_index;
2345 }
2346 }
2347
2353 {
2363 boost::container::small_vector<::types::global_dof_index, 27>;
2364
2375 template <int dim,
2376 int spacedim,
2377 int structdim,
2378 typename GlobalIndexType,
2379 typename DoFPProcessor>
2380 static void
2382 const unsigned int obj_level,
2383 const unsigned int obj_index,
2384 const types::fe_index fe_index,
2385 const unsigned int local_index,
2386 const std::integral_constant<int, structdim> &,
2387 GlobalIndexType &global_index,
2388 const DoFPProcessor &process)
2389 {
2390 Assert(structdim == dim || obj_level == 0, ExcNotImplemented());
2391
2392 // 1) no hp used -> fe_index == 0
2393 if (dof_handler.hp_capability_enabled == false)
2394 {
2395 AssertDimension(fe_index,
2397
2398 process(
2399 dof_handler.object_dof_indices
2400 [obj_level][structdim]
2401 [dof_handler.object_dof_ptr[obj_level][structdim][obj_index] +
2402 local_index],
2403 global_index);
2404
2405 return;
2406 }
2407
2408 // 2) cell and hp is used -> there is only one fe_index
2409 if (structdim == dim)
2410 {
2411 process(
2412 dof_handler.object_dof_indices
2413 [obj_level][structdim]
2414 [dof_handler.object_dof_ptr[obj_level][structdim][obj_index] +
2415 local_index],
2416 global_index);
2417 return;
2418 }
2419
2420 // 3) general entity and hp is used
2421 AssertIndexRange(obj_level, dof_handler.object_dof_indices.size());
2422 AssertIndexRange(structdim,
2423 dof_handler.object_dof_indices[obj_level].size());
2424
2426
2427 AssertIndexRange(structdim, dof_handler.hp_object_fe_ptr.size());
2428 AssertIndexRange(obj_index,
2429 dof_handler.hp_object_fe_ptr[structdim].size());
2430
2431 const auto ptr =
2432 std::find(dof_handler.hp_object_fe_indices[structdim].begin() +
2433 dof_handler.hp_object_fe_ptr[structdim][obj_index],
2434 dof_handler.hp_object_fe_indices[structdim].begin() +
2435 dof_handler.hp_object_fe_ptr[structdim][obj_index + 1],
2436 fe_index);
2437
2438 Assert(ptr != dof_handler.hp_object_fe_indices[structdim].begin() +
2439 dof_handler.hp_object_fe_ptr[structdim][obj_index + 1],
2440 ExcMessage(
2441 "You are requesting an active FE index that is not assigned "
2442 "to any of the cells connected to this entity."));
2443
2444 const types::fe_index fe_index_ =
2445 std::distance(dof_handler.hp_object_fe_indices[structdim].begin() +
2446 dof_handler.hp_object_fe_ptr[structdim][obj_index],
2447 ptr);
2448
2450 dof_handler.hp_capability_enabled ?
2451 (dof_handler.hp_object_fe_ptr[structdim][obj_index] + fe_index_) :
2452 obj_index,
2453 dof_handler.object_dof_ptr[obj_level][structdim].size());
2454
2456 dof_handler.object_dof_ptr
2457 [obj_level][structdim]
2458 [dof_handler.hp_capability_enabled ?
2459 (dof_handler.hp_object_fe_ptr[structdim][obj_index] +
2460 fe_index_) :
2461 obj_index] +
2462 local_index,
2463 dof_handler.object_dof_indices[obj_level][structdim].size());
2464
2465 process(dof_handler.object_dof_indices
2466 [obj_level][structdim]
2467 [dof_handler.object_dof_ptr
2468 [obj_level][structdim]
2469 [dof_handler.hp_capability_enabled ?
2470 (dof_handler.hp_object_fe_ptr[structdim][obj_index] +
2471 fe_index_) :
2472 obj_index] +
2473 local_index],
2474 global_index);
2475 }
2476
2486 template <int dim, int spacedim, int structdim>
2487 static std::pair<unsigned int, unsigned int>
2489 const unsigned int obj_level,
2490 const unsigned int obj_index,
2491 const types::fe_index fe_index,
2492 const std::integral_constant<int, structdim> &)
2493 {
2494 Assert(structdim == dim || obj_level == 0, ExcNotImplemented());
2495
2496 // determine range of dofs in global data structure
2497 // 1) cell
2498 if (structdim == dim)
2499 {
2500 const unsigned int ptr_0 =
2501 dof_handler.object_dof_ptr[obj_level][structdim][obj_index];
2502 const unsigned int length =
2503 dof_handler.get_fe(fe_index).template n_dofs_per_object<dim>(0);
2504
2505 return {ptr_0, length};
2506 }
2507
2508 // 2) hp is not used -> fe_index == 0
2509 if (dof_handler.hp_capability_enabled == false)
2510 {
2511 AssertDimension(fe_index,
2513
2514 const unsigned int ptr_0 =
2515 dof_handler.object_dof_ptr[obj_level][structdim][obj_index];
2516 const unsigned int length =
2517 dof_handler.object_dof_ptr[obj_level][structdim][obj_index + 1] -
2518 ptr_0;
2519
2520 return {ptr_0, length};
2521 }
2522
2523 // 3) hp is used
2524 AssertIndexRange(obj_level, dof_handler.object_dof_indices.size());
2525 AssertIndexRange(structdim,
2526 dof_handler.object_dof_indices[obj_level].size());
2527
2528 AssertIndexRange(structdim, dof_handler.hp_object_fe_ptr.size());
2529 AssertIndexRange(obj_index,
2530 dof_handler.hp_object_fe_ptr[structdim].size());
2531
2532 const auto fe_index_local_ptr =
2533 std::find(dof_handler.hp_object_fe_indices[structdim].begin() +
2534 dof_handler.hp_object_fe_ptr[structdim][obj_index],
2535 dof_handler.hp_object_fe_indices[structdim].begin() +
2536 dof_handler.hp_object_fe_ptr[structdim][obj_index + 1],
2537 fe_index);
2538
2539 Assert(fe_index_local_ptr !=
2540 dof_handler.hp_object_fe_indices[structdim].begin() +
2541 dof_handler.hp_object_fe_ptr[structdim][obj_index + 1],
2542 ExcMessage(
2543 "You tried to call a function accessing DoF indices, but "
2544 "they appear not be available (yet) or inconsistent. "
2545 "Did you call distribute_dofs() first? Alternatively, if "
2546 "you are using different elements on different cells (i.e., "
2547 "you are using the hp capabilities of deal.II), did you "
2548 "change the active_fe_index of a cell since you last "
2549 "called distribute_dofs()?"));
2550
2551 const types::fe_index fe_index_local =
2552 std::distance(dof_handler.hp_object_fe_indices[structdim].begin() +
2553 dof_handler.hp_object_fe_ptr[structdim][obj_index],
2554 fe_index_local_ptr);
2555
2557 dof_handler.hp_object_fe_ptr[structdim][obj_index] + fe_index_local,
2558 dof_handler.object_dof_ptr[obj_level][structdim].size());
2559
2560 const unsigned int ptr_0 =
2561 dof_handler
2562 .object_dof_ptr[obj_level][structdim]
2563 [dof_handler.hp_object_fe_ptr[structdim][obj_index] +
2564 fe_index_local];
2565 const unsigned int ptr_1 =
2566 dof_handler
2567 .object_dof_ptr[obj_level][structdim]
2568 [dof_handler.hp_object_fe_ptr[structdim][obj_index] +
2569 fe_index_local + 1];
2570
2571 return {ptr_0, ptr_1 - ptr_0};
2572 }
2573
2574 template <int dim, int spacedim, int structdim, bool level_dof_access>
2575 static std::pair<unsigned int, unsigned int>
2577 const ::DoFAccessor<structdim, dim, spacedim, level_dof_access>
2578 accessor,
2579 const types::fe_index fe_index)
2580 {
2581 return process_object_range(accessor.get_dof_handler(),
2582 accessor.level(),
2583 accessor.index(),
2584 fe_index,
2585 std::integral_constant<int, structdim>());
2586 }
2587
2588 template <int dim, int spacedim, int structdim>
2589 static std::pair<unsigned int, unsigned int>
2591 const unsigned int)
2592 {
2594
2595 return {0, 0};
2596 }
2597
2598
2599
2608 template <int dim,
2609 int spacedim,
2610 int structdim,
2611 typename DoFProcessor,
2612 typename DoFMapping>
2613 static DEAL_II_ALWAYS_INLINE void
2615 const unsigned int obj_level,
2616 const unsigned int obj_index,
2617 const types::fe_index fe_index,
2618 const DoFMapping &mapping,
2619 const std::integral_constant<int, structdim> &dd,
2620 types::global_dof_index *&dof_indices_ptr,
2621 const DoFProcessor &process)
2622 {
2623 Assert(structdim == dim || obj_level == 0, ExcNotImplemented());
2624
2625 // determine range of dofs in global data structure
2626 const auto range =
2627 process_object_range(dof_handler, obj_level, obj_index, fe_index, dd);
2628 if (range.second == 0)
2629 return;
2630
2631 std::vector<types::global_dof_index> &object_dof_indices =
2632 dof_handler
2633 .object_dof_indices[structdim < dim ? 0 : obj_level][structdim];
2634 AssertIndexRange(range.first, object_dof_indices.size());
2636 object_dof_indices.data() + range.first;
2637
2638 // process dofs
2639 for (unsigned int i = 0; i < range.second; ++i)
2640 {
2641 process(
2642 stored_indices[(structdim == 0 || structdim == dim) ? i :
2643 mapping(i)],
2644 dof_indices_ptr);
2645 if (dof_indices_ptr != nullptr)
2646 ++dof_indices_ptr;
2647 }
2648 }
2649
2650
2651
2662 template <int dim, int spacedim, int structdim>
2663 static void
2665 const unsigned int obj_level,
2666 const unsigned int obj_index,
2667 const types::fe_index fe_index,
2668 const unsigned int local_index,
2669 const std::integral_constant<int, structdim> &dd,
2670 const types::global_dof_index global_index)
2671 {
2672 process_dof_index(dof_handler,
2673 obj_level,
2674 obj_index,
2675 fe_index,
2676 local_index,
2677 dd,
2678 global_index,
2679 [](auto &ptr, const auto &value) { ptr = value; });
2680 }
2681
2682
2693 template <int dim, int spacedim, int structdim>
2696 const unsigned int obj_level,
2697 const unsigned int obj_index,
2698 const types::fe_index fe_index,
2699 const unsigned int local_index,
2700 const std::integral_constant<int, structdim> &dd)
2701 {
2702 types::global_dof_index global_index;
2703 process_dof_index(dof_handler,
2704 obj_level,
2705 obj_index,
2706 fe_index,
2707 local_index,
2708 dd,
2709 global_index,
2710 [](const auto &ptr, auto &value) { value = ptr; });
2711 return global_index;
2712 }
2713
2714
2715 template <int dim, int spacedim>
2718 const int level,
2719 const unsigned int vertex_index,
2720 const unsigned int i)
2721 {
2722 Assert(dof_handler.hp_capability_enabled == false,
2723 ExcMessage(
2724 "DoFHandler in hp-mode does not implement multilevel DoFs."));
2725
2726 return dof_handler.mg_vertex_dofs[vertex_index].access_index(
2727 level, i, dof_handler.get_fe().n_dofs_per_vertex());
2728 }
2729
2730
2731
2741 template <int dim, int spacedim, int structdim>
2742 static unsigned int
2744 const unsigned int obj_level,
2745 const unsigned int obj_index,
2746 const std::integral_constant<int, structdim> &)
2747 {
2748 (void)obj_level;
2749
2750 Assert(structdim == dim || obj_level == 0, ExcNotImplemented());
2751
2752 // 1) no hp used -> fe_index == 0
2753 if (dof_handler.hp_capability_enabled == false)
2754 return 1;
2755
2756 // 2) cell and hp is used -> there is only one fe_index
2757 if (structdim == dim)
2758 return 1;
2759
2760 // 3) general entity and hp is used
2761 AssertIndexRange(structdim, dof_handler.hp_object_fe_ptr.size());
2762 AssertIndexRange(obj_index + 1,
2763 dof_handler.hp_object_fe_ptr[structdim].size());
2764
2765 return dof_handler.hp_object_fe_ptr[structdim][obj_index + 1] -
2766 dof_handler.hp_object_fe_ptr[structdim][obj_index];
2767 }
2768
2769
2770
2780 template <int dim, int spacedim, int structdim>
2781 static types::fe_index
2783 const unsigned int obj_level,
2784 const unsigned int obj_index,
2785 const unsigned int local_index,
2786 const std::integral_constant<int, structdim> &)
2787 {
2788 Assert(structdim == dim || obj_level == 0, ExcNotImplemented());
2789
2790 // for cells only one active FE index available
2791 Assert(((structdim == dim) &&
2793 false,
2795
2796 // 1) no hp used -> fe_index == 0
2797 if (dof_handler.hp_capability_enabled == false)
2799
2800 // 2) cell and hp is used -> there is only one fe_index
2801 if (structdim == dim)
2802 return dof_handler.hp_cell_active_fe_indices[obj_level][obj_index];
2803
2804 // 3) general entity and hp is used
2805 AssertIndexRange(structdim, dof_handler.hp_object_fe_indices.size());
2806 AssertIndexRange(structdim, dof_handler.hp_object_fe_ptr.size());
2807 AssertIndexRange(obj_index,
2808 dof_handler.hp_object_fe_ptr[structdim].size());
2809 AssertIndexRange(dof_handler.hp_object_fe_ptr[structdim][obj_index] +
2810 local_index,
2811 dof_handler.hp_object_fe_indices[structdim].size());
2812
2813 return dof_handler.hp_object_fe_indices
2814 [structdim]
2815 [dof_handler.hp_object_fe_ptr[structdim][obj_index] + local_index];
2816 }
2817
2818
2819
2832 template <int dim, int spacedim, int structdim>
2833 static std::set<types::fe_index>
2835 const unsigned int obj_level,
2836 const unsigned int obj_index,
2837 const std::integral_constant<int, structdim> &t)
2838 {
2839 Assert(structdim == dim || obj_level == 0, ExcNotImplemented());
2840
2841 // 1) no hp used -> fe_index == 0
2842 if (dof_handler.hp_capability_enabled == false)
2844
2845 // 2) cell and hp is used -> there is only one fe_index
2846 if (structdim == dim)
2847 return {dof_handler.hp_cell_active_fe_indices[obj_level][obj_index]};
2848
2849 // 3) general entity and hp is used
2850 std::set<types::fe_index> active_fe_indices;
2851 for (unsigned int i = 0;
2852 i < n_active_fe_indices(dof_handler, obj_level, obj_index, t);
2853 ++i)
2854 active_fe_indices.insert(
2855 nth_active_fe_index(dof_handler, obj_level, obj_index, i, t));
2856 return active_fe_indices;
2857 }
2858
2859
2860
2861 template <int dim, int spacedim, int structdim>
2862 static bool
2864 const unsigned int obj_level,
2865 const unsigned int obj_index,
2866 const types::fe_index fe_index,
2867 const std::integral_constant<int, structdim> &)
2868 {
2869 Assert(structdim == dim || obj_level == 0, ExcNotImplemented());
2870
2871 // 1) no hp used -> fe_index == 0
2872 if (dof_handler.hp_capability_enabled == false)
2874
2875 // 2) cell and hp is used -> there is only one fe_index
2876 if (structdim == dim)
2877 return dof_handler.hp_cell_active_fe_indices[obj_level][obj_index] ==
2878 fe_index;
2879
2880 // 3) general entity and hp is used
2881 return std::find(
2882 dof_handler.hp_object_fe_indices[structdim].begin() +
2883 dof_handler.hp_object_fe_ptr[structdim][obj_index],
2884 dof_handler.hp_object_fe_indices[structdim].begin() +
2885 dof_handler.hp_object_fe_ptr[structdim][obj_index + 1],
2886 fe_index) !=
2887 (dof_handler.hp_object_fe_indices[structdim].begin() +
2888 dof_handler.hp_object_fe_ptr[structdim][obj_index + 1]);
2889 }
2890
2891
2892
2893 template <typename InputVector, typename ForwardIterator>
2894 static void
2895 extract_subvector_to(const InputVector &values,
2896 const types::global_dof_index *cache,
2897 const types::global_dof_index *cache_end,
2898 ForwardIterator local_values_begin)
2899 {
2900 values.extract_subvector_to(cache, cache_end, local_values_begin);
2901 }
2902
2903
2904
2905#ifdef DEAL_II_WITH_TRILINOS
2906 static std::vector<unsigned int>
2908 const types::global_dof_index *v_end)
2909 {
2910 // initialize original index locations
2911 std::vector<unsigned int> idx(v_end - v_begin);
2912 std::iota(idx.begin(), idx.end(), 0u);
2913
2914 // sort indices based on comparing values in v
2915 std::sort(idx.begin(),
2916 idx.end(),
2917 [&v_begin](unsigned int i1, unsigned int i2) {
2918 return *(v_begin + i1) < *(v_begin + i2);
2919 });
2920
2921 return idx;
2922 }
2923
2924
2925
2926# ifdef DEAL_II_TRILINOS_WITH_TPETRA
2927 template <typename ForwardIterator, typename Number, typename MemorySpace>
2928 static void
2931 &values,
2932 const types::global_dof_index *cache_begin,
2933 const types::global_dof_index *cache_end,
2934 ForwardIterator local_values_begin)
2935 {
2936 std::vector<unsigned int> sorted_indices_pos =
2937 sort_indices(cache_begin, cache_end);
2938 const unsigned int cache_size = cache_end - cache_begin;
2939 std::vector<types::global_dof_index> cache_indices(cache_size);
2940 for (unsigned int i = 0; i < cache_size; ++i)
2941 cache_indices[i] = *(cache_begin + sorted_indices_pos[i]);
2942
2943 IndexSet index_set(cache_indices.back() + 1);
2944 index_set.add_indices(cache_indices.begin(), cache_indices.end());
2945 index_set.compress();
2946 LinearAlgebra::ReadWriteVector<Number> read_write_vector(index_set);
2947 read_write_vector.import_elements(values, VectorOperation::insert);
2948
2949 // Copy the elements from read_write_vector and reorder them.
2950 for (unsigned int i = 0; i < cache_size; ++i, ++local_values_begin)
2951 *local_values_begin = read_write_vector[sorted_indices_pos[i]];
2952 }
2953# endif
2954
2955
2956
2957 template <typename ForwardIterator>
2958 static void
2960 const types::global_dof_index *cache_begin,
2961 const types::global_dof_index *cache_end,
2962 ForwardIterator local_values_begin)
2963 {
2964 std::vector<unsigned int> sorted_indices_pos =
2965 sort_indices(cache_begin, cache_end);
2966 const unsigned int cache_size = cache_end - cache_begin;
2967 std::vector<types::global_dof_index> cache_indices(cache_size);
2968 for (unsigned int i = 0; i < cache_size; ++i)
2969 cache_indices[i] = *(cache_begin + sorted_indices_pos[i]);
2970
2971 IndexSet index_set(cache_indices.back() + 1);
2972 index_set.add_indices(cache_indices.begin(), cache_indices.end());
2973 index_set.compress();
2974 LinearAlgebra::ReadWriteVector<double> read_write_vector(index_set);
2975 read_write_vector.import_elements(values, VectorOperation::insert);
2976
2977 // Copy the elements from read_write_vector and reorder them.
2978 for (unsigned int i = 0; i < cache_size; ++i, ++local_values_begin)
2979 *local_values_begin = read_write_vector[sorted_indices_pos[i]];
2980 }
2981#endif
2982
2987 template <int dim, int spacedim, bool level_dof_access, int structdim>
2988 static unsigned int
2990 const ::DoFAccessor<structdim, dim, spacedim, level_dof_access>
2991 &accessor,
2992 const types::fe_index fe_index_,
2993 const bool count_level_dofs)
2994 {
2995 // note: we cannot rely on the template parameter level_dof_access here,
2996 // since the function get_mg_dof_indices()/set_mg_dof_indices() can be
2997 // called even if level_dof_access==false.
2998 if (count_level_dofs)
2999 {
3000 const auto &fe = accessor.get_fe(fe_index_);
3001
3002 const unsigned int //
3003 dofs_per_vertex = fe.n_dofs_per_vertex(), //
3004 dofs_per_line = fe.n_dofs_per_line(), //
3005 dofs_per_quad = fe.n_dofs_per_quad(0 /*dummy*/), //
3006 dofs_per_hex = fe.n_dofs_per_hex(); //
3007
3008 unsigned int index = 0;
3009
3010 // 1) VERTEX dofs
3011 index += dofs_per_vertex * accessor.n_vertices();
3012
3013 // 2) LINE dofs
3014 if (structdim == 2 || structdim == 3)
3015 index += dofs_per_line * accessor.n_lines();
3016
3017 // 3) FACE dofs
3018 if (structdim == 3)
3019 index += dofs_per_quad * accessor.n_faces();
3020
3021 // 4) INNER dofs
3022 const unsigned int interior_dofs =
3023 structdim == 1 ? dofs_per_line :
3024 (structdim == 2 ? dofs_per_quad : dofs_per_hex);
3025
3026 index += interior_dofs;
3027
3028 return index;
3029 }
3030 else
3031 {
3032 const auto fe_index =
3034 accessor, fe_index_);
3035
3036 unsigned int index = 0;
3037
3038 // 1) VERTEX dofs
3039 for (const auto vertex : accessor.vertex_indices())
3040 index += process_object_range(accessor.get_dof_handler(),
3041 0,
3042 accessor.vertex_index(vertex),
3043 fe_index,
3044 std::integral_constant<int, 0>())
3045 .second;
3046
3047 // 2) LINE dofs
3048 if (structdim == 2 || structdim == 3)
3049 for (const auto line : accessor.line_indices())
3050 index +=
3051 process_object_range(*accessor.line(line), fe_index).second;
3052
3053 // 3) FACE dofs
3054 if (structdim == 3)
3055 for (const auto face : accessor.face_indices())
3056 index +=
3057 process_object_range(*accessor.quad(face), fe_index).second;
3058
3059 // 4) INNER dofs
3060 index += process_object_range(accessor, fe_index).second;
3061
3062 return index;
3063 }
3064 }
3065
3066
3067
3068 // The next few internal helper functions are needed to support various
3069 // DoFIndicesType kinds, e.g. actual vectors of DoFIndices or empty
3070 // types that we use when we only want to work on the internally stored
3071 // DoFs and never extract any number.
3072 template <typename ArrayType>
3073 static unsigned int
3074 get_array_length(const ArrayType &array)
3075 {
3076 return array.size();
3077 }
3078
3079 static unsigned int
3080 get_array_length(const std::tuple<> &)
3081 {
3082 return 0;
3083 }
3084
3085 template <typename ArrayType>
3087 get_array_ptr(const ArrayType &array)
3088 {
3089 return const_cast<types::global_dof_index *>(array.data());
3090 }
3091
3093 get_array_ptr(const std::tuple<> &)
3094 {
3095 return nullptr;
3096 }
3097
3098
3099
3105 template <int dim,
3106 int spacedim,
3107 bool level_dof_access,
3108 int structdim,
3109 typename DoFIndicesType,
3110 typename DoFOperation,
3111 typename DoFProcessor>
3112 static void
3114 const ::DoFAccessor<structdim, dim, spacedim, level_dof_access>
3115 &accessor,
3116 const DoFIndicesType &const_dof_indices,
3117 const types::fe_index fe_index_,
3118 const DoFOperation &dof_operation,
3119 const DoFProcessor &dof_processor,
3120 const bool count_level_dofs)
3121 {
3122 const types::fe_index fe_index =
3124 accessor, fe_index_);
3125
3126 // we cannot rely on the template parameter level_dof_access here, since
3127 // the function get_mg_dof_indices()/set_mg_dof_indices() can be called
3128 // even if level_dof_access==false.
3129 (void)count_level_dofs;
3130
3131 const auto &fe = accessor.get_fe(fe_index);
3132
3133 // we want to pass in rvalue 'std::tuple<>' types as `DoFIndicesType`,
3134 // but we need non-const references for std::vector<> types, so get in
3135 // a const reference here and immediately cast the constness away -
3136 // note that any use of the dereferenced invalid type will result in a
3137 // segfault
3138 types::global_dof_index *dof_indices_ptr =
3139 get_array_ptr(const_dof_indices);
3140 types::global_dof_index *end_dof_indices =
3141 dof_indices_ptr + get_array_length(const_dof_indices);
3142
3143 // 1) VERTEX dofs, only step into the functions if we actually have
3144 // DoFs on them
3145 if (fe.n_dofs_per_vertex() > 0)
3146 for (const auto vertex : accessor.vertex_indices())
3147 dof_operation.process_vertex_dofs(*accessor.dof_handler,
3148 accessor.vertex_index(vertex),
3149 fe_index,
3150 dof_indices_ptr,
3151 dof_processor);
3152
3153 // 2) copy dof numbers from the LINE, accounting for the possibility of
3154 // reversed line orientations.
3155 if (structdim > 1 && fe.n_dofs_per_line() > 0)
3156 {
3157 const auto line_indices = internal::TriaAccessorImplementation::
3158 Implementation::get_line_indices_of_cell(accessor);
3159 const auto line_orientations =
3160 internal::TriaAccessorImplementation::Implementation::
3161 get_line_orientations_of_cell(accessor);
3162
3163 for (const auto line : accessor.line_indices())
3164 {
3165 const auto line_orientation = line_orientations[line];
3166 if (line_orientation == numbers::default_geometric_orientation)
3167 dof_operation.process_dofs(
3168 accessor.get_dof_handler(),
3169 0,
3170 line_indices[line],
3171 fe_index,
3172 [](const auto d) { return d; },
3173 std::integral_constant<int, 1>(),
3174 dof_indices_ptr,
3175 dof_processor);
3176 else
3177 {
3178 Assert(line_orientation ==
3181 dof_operation.process_dofs(
3182 accessor.get_dof_handler(),
3183 0,
3184 line_indices[line],
3185 fe_index,
3186 [&fe, line_orientation](const auto d) {
3187 return fe.adjust_line_dof_index_for_line_orientation(
3188 d, line_orientation);
3189 },
3190 std::integral_constant<int, 1>(),
3191 dof_indices_ptr,
3192 dof_processor);
3193 }
3194 }
3195 }
3196
3197 // 3) copy dof numbers from the QUAD (i.e., 3d faces). Like lines we
3198 // only adjust dof indices (which has some cost) if we are not in the
3199 // default orientation.
3200 if (structdim == 3 && fe.max_dofs_per_quad() > 0)
3201 for (const auto face_no : accessor.face_indices())
3202 {
3203 const auto combined_orientation =
3204 accessor.combined_face_orientation(face_no);
3205 const unsigned int quad_index = accessor.quad_index(face_no);
3206 if (combined_orientation ==
3208 dof_operation.process_dofs(
3209 accessor.get_dof_handler(),
3210 0,
3211 quad_index,
3212 fe_index,
3213 [](const auto d) { return d; },
3214 std::integral_constant<int, 2>(),
3215 dof_indices_ptr,
3216 dof_processor);
3217 else
3218 dof_operation.process_dofs(
3219 accessor.get_dof_handler(),
3220 0,
3221 quad_index,
3222 fe_index,
3223 [&](const auto d) {
3224 return fe.adjust_quad_dof_index_for_face_orientation(
3225 d, face_no, combined_orientation);
3226 },
3227 std::integral_constant<int, 2>(),
3228 dof_indices_ptr,
3229 dof_processor);
3230 }
3231
3232 // 4) INNER dofs (i.e., line dofs in 1d, quad dofs in 2d, or hex dofs in
3233 // 3d) - here we need to make sure that the shortcut to not run the
3234 // function does not miss the faces of wedge and pyramid elements where
3235 // n_dofs_per_object might not return the largest possible value
3236 if (((dim == 3 && structdim == 2) ?
3237 fe.max_dofs_per_quad() :
3238 fe.template n_dofs_per_object<structdim>()) > 0)
3239 dof_operation.process_dofs(
3240 accessor.get_dof_handler(),
3241 accessor.level(),
3242 accessor.index(),
3243 fe_index,
3244 [&](const auto d) { return d; },
3245 std::integral_constant<int, structdim>(),
3246 dof_indices_ptr,
3247 dof_processor);
3248
3249 if (dof_indices_ptr != nullptr)
3250 {
3251 AssertDimension(n_dof_indices(accessor, fe_index, count_level_dofs),
3252 dof_indices_ptr - get_array_ptr(const_dof_indices));
3253 }
3254
3255 // PM: This is a part that should not be reached since it indicates that
3256 // an object (and/or its subobjects) is not active. However,
3257 // unfortunately this function is called by
3258 // DoFTools::set_periodicity_constraints() indirectly by
3259 // get_dof_indices() also for artificial faces to determine if a face
3260 // is artificial.
3262 for (; dof_indices_ptr < end_dof_indices; ++dof_indices_ptr)
3263 dof_processor(invalid_index, dof_indices_ptr);
3264 }
3265
3266
3267
3272 template <int dim, int spacedim>
3274 {
3278 template <typename DoFProcessor>
3281 const unsigned int vertex_index,
3282 const types::fe_index fe_index,
3283 types::global_dof_index *&dof_indices_ptr,
3284 const DoFProcessor &dof_processor) const
3285 {
3287 dof_handler,
3288 0,
3289 vertex_index,
3290 fe_index,
3291 [](const auto d) {
3293 return d;
3294 },
3295 std::integral_constant<int, 0>(),
3296 dof_indices_ptr,
3297 dof_processor);
3298 }
3299
3303 template <int structdim, typename DoFMapping, typename DoFProcessor>
3306 const unsigned int obj_level,
3307 const unsigned int obj_index,
3308 const types::fe_index fe_index,
3309 const DoFMapping &mapping,
3310 const std::integral_constant<int, structdim>,
3311 types::global_dof_index *&dof_indices_ptr,
3312 const DoFProcessor &dof_processor) const
3313 {
3315 dof_handler,
3316 obj_level,
3317 obj_index,
3318 fe_index,
3319 mapping,
3320 std::integral_constant<int, std::min(structdim, dim)>(),
3321 dof_indices_ptr,
3322 dof_processor);
3323 }
3324 };
3325
3326
3327
3332 template <int dim, int spacedim>
3334 {
3338 MGDoFIndexProcessor(const unsigned int level)
3339 : level(level)
3340 {}
3341
3345 template <typename DoFProcessor>
3348 const unsigned int vertex_index,
3349 const types::fe_index,
3350 types::global_dof_index *&dof_indices_ptr,
3351 const DoFProcessor &dof_processor) const
3352 {
3353 const unsigned int n_indices =
3354 dof_handler.get_fe(0).template n_dofs_per_object<0>();
3355 types::global_dof_index *stored_indices =
3356 &dof_handler.mg_vertex_dofs[vertex_index].access_index(level,
3357 0,
3358 n_indices);
3359 for (unsigned int d = 0; d < n_indices; ++d, ++dof_indices_ptr)
3360 dof_processor(stored_indices[d], dof_indices_ptr);
3361 }
3362
3366 template <int structdim, typename DoFMapping, typename DoFProcessor>
3369 const unsigned int,
3370 const unsigned int obj_index,
3371 const types::fe_index fe_index,
3372 const DoFMapping &mapping,
3373 const std::integral_constant<int, structdim>,
3374 types::global_dof_index *&dof_indices_ptr,
3375 const DoFProcessor &dof_processor) const
3376 {
3377 const unsigned int n_indices =
3378 dof_handler.get_fe(0).template n_dofs_per_object<structdim>();
3379 types::global_dof_index *stored_indices = &get_mg_dof_index(
3380 dof_handler,
3381 dof_handler.mg_levels[level],
3382 dof_handler.mg_faces,
3383 obj_index,
3384 fe_index,
3385 0,
3386 std::integral_constant<int, std::min(structdim, dim)>());
3387 for (unsigned int d = 0; d < n_indices; ++d, ++dof_indices_ptr)
3388 dof_processor(stored_indices[structdim < dim ? mapping(d) : d],
3389 dof_indices_ptr);
3390 }
3391
3392 private:
3393 const unsigned int level;
3394 };
3395
3396
3397
3398 template <int dim, int spacedim, bool level_dof_access, int structdim>
3399 static void
3401 const ::DoFAccessor<structdim, dim, spacedim, level_dof_access>
3402 &accessor,
3403 std::vector<types::global_dof_index> &dof_indices,
3404 const types::fe_index fe_index)
3405 {
3407 accessor,
3408 dof_indices,
3409 fe_index,
3411 [](auto stored_index, auto dof_ptr) { *dof_ptr = stored_index; },
3412 false);
3413 }
3414
3415
3416
3417 template <int dim, int spacedim, bool level_dof_access, int structdim>
3418 static void
3420 const ::DoFAccessor<structdim, dim, spacedim, level_dof_access>
3421 &accessor,
3422 const std::vector<types::global_dof_index> &dof_indices,
3423 const types::fe_index fe_index)
3424 {
3425 // Note: this function is as general as `get_dof_indices()`. This
3426 // assert is placed here since it is currently only used by the
3427 // function DoFCellAccessor::set_dof_indices(), which is called by
3428 // internal::DoFHandlerImplementation::Policy::Implementation::distribute_dofs().
3429 // In the case of new use cases, this assert can be removed.
3430 Assert(
3431 dim == structdim,
3432 ExcMessage(
3433 "This function is intended to be used for DoFCellAccessor, i.e., "
3434 "dimension == structdim."));
3435
3437 accessor,
3438 dof_indices,
3439 fe_index,
3441 [](auto &stored_index, auto dof_ptr) { stored_index = *dof_ptr; },
3442 false);
3443 }
3444
3445
3446
3447 template <int dim, int spacedim, bool level_dof_access, int structdim>
3448 static void
3450 const ::DoFAccessor<structdim, dim, spacedim, level_dof_access>
3451 &accessor,
3452 const int level,
3453 std::vector<types::global_dof_index> &dof_indices,
3454 const types::fe_index fe_index)
3455 {
3457 ExcMessage("MG DoF indices cannot be queried in hp case"));
3459 accessor,
3460 dof_indices,
3461 fe_index,
3463 [](auto stored_index, auto dof_ptr) { *dof_ptr = stored_index; },
3464 true);
3465 }
3466
3467
3468
3469 template <int dim, int spacedim, bool level_dof_access, int structdim>
3470 static void
3472 const ::DoFAccessor<structdim, dim, spacedim, level_dof_access>
3473 &accessor,
3474 const int level,
3475 const std::vector<types::global_dof_index> &dof_indices,
3476 const types::fe_index fe_index)
3477 {
3479 ExcMessage("MG DoF indices cannot be queried in hp case"));
3480
3481 // Note: this function is as general as `get_mg_dof_indices()`. This
3482 // assert is placed here since it is currently only used by the
3483 // function DoFCellAccessor::set_mg_dof_indices(), which is called by
3484 // internal::DoFHandlerImplementation::Policy::Implementation::distribute_mg_dofs().
3485 // In the case of new use cases, this assert can be removed.
3486 Assert(dim == structdim,
3487 ExcMessage("This function is intended to be used for "
3488 "DoFCellAccessor, i.e., dimension == structdim."));
3489
3491 accessor,
3492 dof_indices,
3493 fe_index,
3495 [](auto &stored_index, auto dof_ptr) { stored_index = *dof_ptr; },
3496 true);
3497 }
3498
3499
3500
3501 template <int dim, int spacedim>
3504 const DoFHandler<dim, spacedim> &dof_handler,
3506 &mg_level,
3508 &,
3509 const unsigned int obj_index,
3510 const types::fe_index fe_index,
3511 const unsigned int local_index,
3512 const std::integral_constant<int, dim>)
3513 {
3514 Assert(dof_handler.hp_capability_enabled == false,
3516
3517 return mg_level->dof_object.access_dof_index(
3518 static_cast<const DoFHandler<dim, spacedim> &>(dof_handler),
3519 obj_index,
3520 fe_index,
3521 local_index);
3522 }
3523
3524
3525
3526 template <int dim, int spacedim, std::enable_if_t<(dim > 1), int> = 0>
3529 const DoFHandler<dim, spacedim> &dof_handler,
3531 &,
3533 &mg_faces,
3534 const unsigned int obj_index,
3535 const types::fe_index fe_index,
3536 const unsigned int local_index,
3537 const std::integral_constant<int, 1>)
3538 {
3539 return mg_faces->lines.access_dof_index(
3540 static_cast<const DoFHandler<dim, spacedim> &>(dof_handler),
3541 obj_index,
3542 fe_index,
3543 local_index);
3544 }
3545
3546
3547
3548 template <int spacedim>
3551 const DoFHandler<3, spacedim> &dof_handler,
3553 &,
3555 &mg_faces,
3556 const unsigned int obj_index,
3557 const types::fe_index fe_index,
3558 const unsigned int local_index,
3559 const std::integral_constant<int, 2>)
3560 {
3561 Assert(dof_handler.hp_capability_enabled == false,
3563 return mg_faces->quads.access_dof_index(
3564 static_cast<const DoFHandler<3, spacedim> &>(dof_handler),
3565 obj_index,
3566 fe_index,
3567 local_index);
3568 }
3569 };
3570
3571
3572
3573 template <int dim, int spacedim, bool level_dof_access>
3574 void
3576 const ::DoFCellAccessor<dim, spacedim, level_dof_access> &accessor,
3578 const unsigned int fe_index);
3579 } // namespace DoFAccessorImplementation
3580} // namespace internal
3581
3582
3583
3584template <int structdim, int dim, int spacedim, bool level_dof_access>
3587 const unsigned int i,
3588 const types::fe_index fe_index_) const
3589{
3590 const auto fe_index =
3592 fe_index_);
3593
3594 // access the respective DoF
3595 return ::internal::DoFAccessorImplementation::Implementation::
3596 get_dof_index(*this->dof_handler,
3597 this->level(),
3598 this->index(),
3599 fe_index,
3600 i,
3601 std::integral_constant<int, structdim>());
3602}
3603
3604
3605template <int structdim, int dim, int spacedim, bool level_dof_access>
3608 const int level,
3609 const unsigned int i) const
3610{
3612 *this->dof_handler,
3613 this->dof_handler->mg_levels[level],
3614 this->dof_handler->mg_faces,
3615 this->index(),
3616 0,
3617 i,
3618 std::integral_constant<int, structdim>());
3619}
3620
3621
3622template <int structdim, int dim, int spacedim, bool level_dof_access>
3623inline void
3625 const unsigned int i,
3626 const types::global_dof_index index,
3627 const types::fe_index fe_index_) const
3628{
3629 const auto fe_index =
3631 fe_index_);
3632
3633 // access the respective DoF
3635 *this->dof_handler,
3636 this->level(),
3637 this->index(),
3638 fe_index,
3639 i,
3640 std::integral_constant<int, structdim>(),
3641 index);
3642}
3643
3644
3645
3646template <int structdim, int dim, int spacedim, bool level_dof_access>
3647inline void
3649 const int level,
3650 const unsigned int i,
3651 const types::global_dof_index index) const
3652{
3654 *this->dof_handler,
3655 this->dof_handler->mg_levels[level],
3656 this->dof_handler->mg_faces,
3657 this->index(),
3658 0,
3659 i,
3660 std::integral_constant<int, structdim>()) = index;
3661}
3662
3663
3664
3665template <int structdim, int dim, int spacedim, bool level_dof_access>
3666inline unsigned int
3668 const
3669{
3670 // access the respective DoF
3671 return ::internal::DoFAccessorImplementation::Implementation::
3672 n_active_fe_indices(*this->dof_handler,
3673 this->level(),
3674 this->index(),
3675 std::integral_constant<int, structdim>());
3676}
3677
3678
3679
3680template <int structdim, int dim, int spacedim, bool level_dof_access>
3681inline types::fe_index
3683 const unsigned int n) const
3684{
3685 // access the respective DoF
3686 return ::internal::DoFAccessorImplementation::Implementation::
3687 nth_active_fe_index(*this->dof_handler,
3688 this->level(),
3689 this->index(),
3690 n,
3691 std::integral_constant<int, structdim>());
3692}
3693
3694
3695
3696template <int structdim, int dim, int spacedim, bool level_dof_access>
3697inline std::set<types::fe_index>
3699 const
3700{
3701 std::set<types::fe_index> active_fe_indices;
3702 for (unsigned int i = 0; i < n_active_fe_indices(); ++i)
3703 active_fe_indices.insert(nth_active_fe_index(i));
3704 return active_fe_indices;
3705}
3706
3707
3708
3709template <int structdim, int dim, int spacedim, bool level_dof_access>
3710inline bool
3712 const types::fe_index fe_index) const
3713{
3714 // access the respective DoF
3715 return ::internal::DoFAccessorImplementation::Implementation::
3716 fe_index_is_active(*this->dof_handler,
3717 this->level(),
3718 this->index(),
3719 fe_index,
3720 std::integral_constant<int, structdim>());
3721}
3722
3723
3724
3725template <int structdim, int dim, int spacedim, bool level_dof_access>
3728 const unsigned int vertex,
3729 const unsigned int i,
3730 const types::fe_index fe_index_) const
3731{
3732 const types::fe_index fe_index =
3733 (((this->dof_handler->hp_capability_enabled == false) &&
3734 (fe_index_ == numbers::invalid_fe_index)) ?
3735 // No hp enabled, and the argument is at its default value -> we
3736 // can translate to the default active fe index
3738 // Otherwise: If anything other than the default is provided by
3739 // the caller, then we should take just that. As an exception, if
3740 // we are on a cell (rather than a face/edge/vertex), then we know
3741 // that there is only one active fe index on this cell and we can
3742 // use that:
3743 ((dim == structdim) && (fe_index_ == numbers::invalid_fe_index) ?
3744 this->nth_active_fe_index(0) :
3745 fe_index_));
3746
3747 return ::internal::DoFAccessorImplementation::Implementation::
3748 get_dof_index(*this->dof_handler,
3749 0,
3750 this->vertex_index(vertex),
3751 fe_index,
3752 i,
3753 std::integral_constant<int, 0>());
3754}
3755
3756
3757template <int structdim, int dim, int spacedim, bool level_dof_access>
3760 const int level,
3761 const unsigned int vertex,
3762 const unsigned int i,
3763 const types::fe_index fe_index_) const
3764{
3765 const auto fe_index =
3767 fe_index_);
3768 (void)fe_index;
3769 Assert(this->dof_handler != nullptr, ExcInvalidObject());
3770 Assert(this->dof_handler->mg_vertex_dofs.size() > 0,
3771 ExcMessage("Multigrid DoF indices can only be accessed after "
3772 "DoFHandler::distribute_mg_dofs() has been called!"));
3773 AssertIndexRange(vertex, this->n_vertices());
3774 AssertIndexRange(i, this->dof_handler->get_fe(fe_index).n_dofs_per_vertex());
3775
3776 Assert(dof_handler->hp_capability_enabled == false,
3777 ExcMessage(
3778 "DoFHandler in hp-mode does not implement multilevel DoFs."));
3779
3780 return this->dof_handler->mg_vertex_dofs[this->vertex_index(vertex)]
3781 .access_index(level, i, this->dof_handler->get_fe().n_dofs_per_vertex());
3782}
3783
3784
3785
3786template <int structdim, int dim, int spacedim, bool level_dof_access>
3787inline void
3789 set_mg_vertex_dof_index(const int level,
3790 const unsigned int vertex,
3791 const unsigned int i,
3792 const types::global_dof_index index,
3793 const types::fe_index fe_index_) const
3794{
3795 const auto fe_index =
3797 fe_index_);
3798 (void)fe_index;
3799 Assert(this->dof_handler != nullptr, ExcInvalidObject());
3800 AssertIndexRange(vertex, this->n_vertices());
3801 AssertIndexRange(i, this->dof_handler->get_fe(fe_index).n_dofs_per_vertex());
3802
3803 Assert(dof_handler->hp_capability_enabled == false,
3804 ExcMessage(
3805 "DoFHandler in hp-mode does not implement multilevel DoFs."));
3806
3807 this->dof_handler->mg_vertex_dofs[this->vertex_index(vertex)].access_index(
3808 level, i, this->dof_handler->get_fe().n_dofs_per_vertex()) = index;
3809}
3810
3811
3812
3813template <int structdim, int dim, int spacedim, bool level_dof_access>
3814inline const FiniteElement<dim, spacedim> &
3816 const types::fe_index fe_index) const
3817{
3818 Assert(fe_index_is_active(fe_index) == true,
3819 ExcMessage("This function can only be called for active FE indices"));
3820
3821 return this->dof_handler->get_fe(fe_index);
3822}
3823
3824
3825
3826template <int structdim, int dim, int spacedim, bool level_dof_access>
3827inline typename ::internal::DoFHandlerImplementation::
3828 Iterators<dim, spacedim, level_dof_access>::line_iterator
3830 const unsigned int i) const
3831{
3832 // if we are asking for a particular line and this object refers to
3833 // a line, then the only valid index is i==0 and we should return
3834 // *this
3835 if (structdim == 1)
3836 {
3837 Assert(i == 0,
3838 ExcMessage("You can only ask for line zero if the "
3839 "current object is a line itself."));
3840 return typename ::internal::DoFHandlerImplementation::
3841 Iterators<dim, spacedim, level_dof_access>::cell_iterator(
3842 &this->get_triangulation(),
3843 this->level(),
3844 this->index(),
3845 &this->get_dof_handler());
3846 }
3847
3848 // otherwise we need to be in structdim>=2
3849 Assert(structdim > 1, ExcImpossibleInDim(structdim));
3850 Assert(dim > 1, ExcImpossibleInDim(dim));
3851
3852 // checking of 'i' happens in line_index(i)
3853 return typename ::internal::DoFHandlerImplementation::
3854 Iterators<dim, spacedim, level_dof_access>::line_iterator(
3855 this->tria,
3856 0, // only sub-objects are allowed, which have no level
3857 this->line_index(i),
3858 this->dof_handler);
3859}
3860
3861
3862template <int structdim, int dim, int spacedim, bool level_dof_access>
3863inline typename ::internal::DoFHandlerImplementation::
3864 Iterators<dim, spacedim, level_dof_access>::quad_iterator
3866 const unsigned int i) const
3867{
3868 // if we are asking for a
3869 // particular quad and this object
3870 // refers to a quad, then the only
3871 // valid index is i==0 and we
3872 // should return *this
3873 if (structdim == 2)
3874 {
3875 Assert(i == 0,
3876 ExcMessage("You can only ask for quad zero if the "
3877 "current object is a quad itself."));
3878 return typename ::internal::DoFHandlerImplementation::
3879 Iterators<dim, spacedim>::cell_iterator(&this->get_triangulation(),
3880 this->level(),
3881 this->index(),
3882 &this->get_dof_handler());
3883 }
3884
3885 // otherwise we need to be in structdim>=3
3886 Assert(structdim > 2, ExcImpossibleInDim(structdim));
3887 Assert(dim > 2, ExcImpossibleInDim(dim));
3888
3889 // checking of 'i' happens in quad_index(i)
3890 return typename ::internal::DoFHandlerImplementation::
3891 Iterators<dim, spacedim, level_dof_access>::quad_iterator(
3892 this->tria,
3893 0, // only sub-objects are allowed, which have no level
3894 this->quad_index(i),
3895 this->dof_handler);
3896}
3897
3898
3899/*----------------- Functions: DoFAccessor<0,1,spacedim> --------------------*/
3900
3901
3902template <int spacedim, bool level_dof_access>
3907
3908
3909
3910template <int spacedim, bool level_dof_access>
3913 const typename TriaAccessor<0, 1, spacedim>::VertexKind vertex_kind,
3914 const unsigned int vertex_index,
3916 : BaseClass(tria, vertex_kind, vertex_index)
3917 , dof_handler(const_cast<DoFHandler<1, spacedim> *>(dof_handler))
3918{}
3919
3920
3921
3922template <int spacedim, bool level_dof_access>
3925 const int level,
3926 const int index,
3928 // This is the constructor signature for "ordinary" (non-vertex)
3929 // accessors and we shouldn't be calling it altogether. But it is also
3930 // the constructor that the default-constructor of TriaRawIterator
3931 // calls when default-constructing an iterator object. If so, this
3932 // happens with level==-2 and index==-2, and this is the only case we
3933 // would like to support. We do this by just forwarding to the
3934 // other constructor of this class, and then asserting the condition
3935 // on level and index.
3936 : DoFAccessor<0, 1, spacedim, level_dof_access>(
3937 tria,
3938 TriaAccessor<0, 1, spacedim>::interior_vertex,
3939 0U,
3941{
3942 (void)level;
3943 (void)index;
3944 Assert((tria == nullptr) && (level == -2) && (index == -2) &&
3945 (dof_handler == nullptr),
3946 ExcMessage(
3947 "This constructor can not be called for face iterators in 1d, "
3948 "except to default-construct iterator objects."));
3949}
3950
3951
3952
3953template <int spacedim, bool level_dof_access>
3954template <int structdim2, int dim2, int spacedim2>
3960
3961
3962
3963template <int spacedim, bool level_dof_access>
3964template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
3970
3971
3972
3973template <int spacedim, bool level_dof_access>
3974inline void
3981
3982
3983
3984template <int spacedim, bool level_dof_access>
3985inline void
3987 const unsigned int /*i*/,
3988 const types::global_dof_index /*index*/,
3989 const types::fe_index /*fe_index*/) const
3990{
3992}
3993
3994
3995
3996template <int spacedim, bool level_dof_access>
3997inline const DoFHandler<1, spacedim> &
4002
4003
4004
4005template <int spacedim, bool level_dof_access>
4006inline void
4008 std::vector<types::global_dof_index> &dof_indices,
4009 const types::fe_index fe_index_) const
4010{
4011 const auto fe_index =
4013 fe_index_);
4014
4015 for (unsigned int i = 0; i < dof_indices.size(); ++i)
4016 dof_indices[i] = ::internal::DoFAccessorImplementation::
4017 Implementation::get_dof_index(*dof_handler,
4018 0,
4019 this->global_vertex_index,
4020 fe_index,
4021 i,
4022 std::integral_constant<int, 0>());
4023}
4024
4025
4026
4027template <int spacedim, bool level_dof_access>
4028inline void
4030 const int level,
4031 std::vector<types::global_dof_index> &dof_indices,
4032 const types::fe_index fe_index_) const
4033{
4034 const auto fe_index =
4036 fe_index_);
4037 (void)fe_index;
4039
4040 for (unsigned int i = 0; i < dof_indices.size(); ++i)
4041 dof_indices[i] =
4043 mg_vertex_dof_index(*dof_handler, level, this->global_vertex_index, i);
4044}
4045
4046
4047
4048template <int spacedim, bool level_dof_access>
4051 const unsigned int vertex,
4052 const unsigned int i,
4053 const types::fe_index fe_index_) const
4054{
4055 const auto fe_index =
4057 fe_index_);
4058
4059 (void)vertex;
4061 return ::internal::DoFAccessorImplementation::Implementation::
4062 get_dof_index(*dof_handler,
4063 0,
4064 this->global_vertex_index,
4065 fe_index,
4066 i,
4067 std::integral_constant<int, 0>());
4068}
4069
4070
4071
4072template <int spacedim, bool level_dof_access>
4075 const unsigned int i,
4076 const types::fe_index fe_index_) const
4077{
4078 const auto fe_index =
4080 fe_index_);
4081
4082 return ::internal::DoFAccessorImplementation::Implementation::
4083 get_dof_index(*this->dof_handler,
4084 0,
4085 this->vertex_index(0),
4086 fe_index,
4087 i,
4088 std::integral_constant<int, 0>());
4089}
4090
4091
4092
4093template <int spacedim, bool level_dof_access>
4094inline unsigned int
4099
4100
4101
4102template <int spacedim, bool level_dof_access>
4103inline types::fe_index
4105 const unsigned int /*n*/) const
4106{
4107 return 0;
4108}
4109
4110
4111
4112template <int spacedim, bool level_dof_access>
4113inline bool
4115 const types::fe_index /*fe_index*/) const
4116{
4118 return false;
4119}
4120
4121
4122
4123template <int spacedim, bool level_dof_access>
4124inline const FiniteElement<1, spacedim> &
4126 const types::fe_index fe_index) const
4127{
4128 Assert(this->dof_handler != nullptr, ExcInvalidObject());
4129 return dof_handler->get_fe(fe_index);
4130}
4131
4132
4133
4134template <int spacedim, bool level_dof_access>
4135inline void
4138{
4139 Assert(this->dof_handler != nullptr, ExcInvalidObject());
4140 BaseClass::copy_from(da);
4141}
4142
4143
4144
4145template <int spacedim, bool level_dof_access>
4146template <bool level_dof_access2>
4147inline void
4154
4155
4156
4157template <int spacedim, bool level_dof_access>
4164
4165
4166
4167template <int spacedim, bool level_dof_access>
4168inline typename ::internal::DoFHandlerImplementation::
4169 Iterators<1, spacedim, level_dof_access>::line_iterator
4171 const unsigned int /*c*/) const
4172{
4174 return typename ::internal::DoFHandlerImplementation::
4175 Iterators<1, spacedim, level_dof_access>::line_iterator();
4176}
4177
4178
4179
4180template <int spacedim, bool level_dof_access>
4181inline typename ::internal::DoFHandlerImplementation::
4182 Iterators<1, spacedim, level_dof_access>::quad_iterator
4184 const unsigned int /*c*/) const
4185{
4187 return typename ::internal::DoFHandlerImplementation::
4188 Iterators<1, spacedim, level_dof_access>::quad_iterator();
4189}
4190
4191
4192
4193template <int spacedim, bool level_dof_access>
4194template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
4195inline bool
4203
4204
4205
4206template <int spacedim, bool level_dof_access>
4207template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
4208inline bool
4216
4217
4218
4219/*------------------------- Functions: DoFCellAccessor -----------------------*/
4220
4221
4222namespace internal
4223{
4225 {
4231 {
4236 template <int dim, int spacedim, bool level_dof_access>
4237 static types::fe_index
4240 {
4241 if (accessor.dof_handler->hp_capability_enabled == false)
4243
4244 Assert(accessor.dof_handler != nullptr,
4245 (typename std::decay_t<decltype(accessor)>::ExcInvalidObject()));
4246 Assert(static_cast<unsigned int>(accessor.level()) <
4247 accessor.dof_handler->hp_cell_future_fe_indices.size(),
4248 ExcMessage("DoFHandler not initialized"));
4249
4250 return accessor.dof_handler
4251 ->hp_cell_active_fe_indices[accessor.level()][accessor.present_index];
4252 }
4253
4254
4255
4260 template <int dim, int spacedim, bool level_dof_access>
4261 static void
4264 const types::fe_index i)
4265 {
4266 if (accessor.dof_handler->hp_capability_enabled == false)
4267 {
4269 return;
4270 }
4271
4272 Assert(accessor.dof_handler != nullptr,
4273 (typename std::decay_t<decltype(accessor)>::ExcInvalidObject()));
4274 Assert(static_cast<unsigned int>(accessor.level()) <
4275 accessor.dof_handler->hp_cell_future_fe_indices.size(),
4276 ExcMessage("DoFHandler not initialized"));
4278 ExcMessage("Invalid finite element index."));
4279
4280 accessor.dof_handler
4281 ->hp_cell_active_fe_indices[accessor.level()]
4282 [accessor.present_index] = i;
4283 }
4284
4285
4286
4291 template <int dim, int spacedim, bool level_dof_access>
4292 static types::fe_index
4295 {
4296 if (accessor.dof_handler->hp_capability_enabled == false)
4298
4299 Assert(accessor.dof_handler != nullptr,
4300 (typename std::decay_t<decltype(accessor)>::ExcInvalidObject()));
4301 Assert(static_cast<unsigned int>(accessor.level()) <
4302 accessor.dof_handler->hp_cell_future_fe_indices.size(),
4303 ExcMessage("DoFHandler not initialized"));
4304
4305 if (future_fe_index_set(accessor))
4306 return accessor.dof_handler
4307 ->hp_cell_future_fe_indices[accessor.level()]
4308 [accessor.present_index];
4309 else
4310 return accessor.dof_handler
4311 ->hp_cell_active_fe_indices[accessor.level()]
4312 [accessor.present_index];
4313 }
4314
4315
4320 template <int dim, int spacedim, bool level_dof_access>
4321 static void
4324 const types::fe_index i)
4325 {
4326 if (accessor.dof_handler->hp_capability_enabled == false)
4327 {
4329 return;
4330 }
4331
4332 Assert(accessor.dof_handler != nullptr,
4333 (typename std::decay_t<decltype(accessor)>::ExcInvalidObject()));
4334 Assert(static_cast<unsigned int>(accessor.level()) <
4335 accessor.dof_handler->hp_cell_future_fe_indices.size(),
4336 ExcMessage("DoFHandler not initialized"));
4338 ExcMessage("Invalid finite element index."));
4339
4340 accessor.dof_handler
4341 ->hp_cell_future_fe_indices[accessor.level()]
4342 [accessor.present_index] = i;
4343 }
4344
4345
4346
4351 template <int dim, int spacedim, bool level_dof_access>
4352 static bool
4355 {
4356 if (accessor.dof_handler->hp_capability_enabled == false)
4357 return false;
4358
4359 Assert(accessor.dof_handler != nullptr,
4360 (typename std::decay_t<decltype(accessor)>::ExcInvalidObject()));
4361 Assert(static_cast<unsigned int>(accessor.level()) <
4362 accessor.dof_handler->hp_cell_future_fe_indices.size(),
4363 ExcMessage("DoFHandler not initialized"));
4364
4365 return accessor.dof_handler
4366 ->hp_cell_future_fe_indices[accessor.level()]
4367 [accessor.present_index] !=
4369 }
4370
4371
4372
4377 template <int dim, int spacedim, bool level_dof_access>
4378 static void
4381 {
4382 if (accessor.dof_handler->hp_capability_enabled == false)
4383 return;
4384
4385 Assert(accessor.dof_handler != nullptr,
4386 (typename std::decay_t<decltype(accessor)>::ExcInvalidObject()));
4387 Assert(static_cast<unsigned int>(accessor.level()) <
4388 accessor.dof_handler->hp_cell_future_fe_indices.size(),
4389 ExcMessage("DoFHandler not initialized"));
4390
4391 accessor.dof_handler
4392 ->hp_cell_future_fe_indices[accessor.level()]
4393 [accessor.present_index] =
4395 }
4396 };
4397 } // namespace DoFCellAccessorImplementation
4398} // namespace internal
4399
4400
4401
4402template <int dimension_, int space_dimension_, bool level_dof_access>
4405 const int level,
4406 const int index,
4407 const AccessorData *local_data)
4408 : DoFAccessor<dimension_, dimension_, space_dimension_, level_dof_access>(
4409 tria,
4410 level,
4411 index,
4412 local_data)
4413{}
4414
4415
4416
4417template <int dimension_, int space_dimension_, bool level_dof_access>
4418template <int structdim2, int dim2, int spacedim2>
4424
4425
4426
4427template <int dimension_, int space_dimension_, bool level_dof_access>
4428template <int structdim2, int dim2, int spacedim2, bool level_dof_access2>
4434
4435
4436
4437template <int dimension_, int space_dimension_, bool level_dof_access>
4438inline TriaIterator<
4441 const unsigned int i) const
4442{
4444 q(this->tria,
4445 this->neighbor_level(i),
4446 this->neighbor_index(i),
4447 this->dof_handler);
4448
4449 if constexpr (running_in_debug_mode())
4450 {
4452 Assert(q->used(), ExcInternalError());
4453 }
4454 return q;
4455}
4456
4457
4458
4459template <int dimension_, int space_dimension_, bool level_dof_access>
4460inline TriaIterator<
4463 const unsigned int i) const
4464{
4466 q(this->tria, this->level() + 1, this->child_index(i), this->dof_handler);
4467
4468 if constexpr (running_in_debug_mode())
4469 {
4471 Assert(q->used(), ExcInternalError());
4472 }
4473 return q;
4474}
4475
4476
4477
4478template <int dimension_, int space_dimension_, bool level_dof_access>
4479inline boost::container::small_vector<
4483 child_iterators() const
4484{
4485 boost::container::small_vector<
4489 child_iterators(this->n_children());
4490
4491 for (unsigned int i = 0; i < this->n_children(); ++i)
4492 child_iterators[i] = this->child(i);
4493
4494 return child_iterators;
4495}
4496
4497
4498
4499template <int dimension_, int space_dimension_, bool level_dof_access>
4500inline TriaIterator<
4503{
4505 q(this->tria, this->level() - 1, this->parent_index(), this->dof_handler);
4506
4507 return q;
4508}
4509
4510
4511
4512namespace internal
4513{
4514 namespace DoFCellAccessorImplementation
4515 {
4516 template <int dim, int spacedim, bool level_dof_access>
4517 inline TriaIterator<
4518 ::DoFAccessor<dim - 1, dim, spacedim, level_dof_access>>
4520 const ::DoFCellAccessor<dim, spacedim, level_dof_access> &cell,
4521 const unsigned int i,
4522 const std::integral_constant<int, 1>)
4523 {
4525 &cell.get_triangulation(),
4526 ((i == 0) && cell.at_boundary(0) ?
4528 ((i == 1) && cell.at_boundary(1) ?
4531 cell.vertex_index(i),
4532 &cell.get_dof_handler());
4533 return ::TriaIterator<
4535 }
4536
4537
4538 template <int dim, int spacedim, bool level_dof_access>
4539 inline TriaIterator<
4540 ::DoFAccessor<dim - 1, dim, spacedim, level_dof_access>>
4542 const ::DoFCellAccessor<dim, spacedim, level_dof_access> &cell,
4543 const unsigned int i,
4544 const std::integral_constant<int, 2>)
4545 {
4546 return cell.line(i);
4547 }
4548
4549
4550 template <int dim, int spacedim, bool level_dof_access>
4551 inline TriaIterator<
4552 ::DoFAccessor<dim - 1, dim, spacedim, level_dof_access>>
4554 const ::DoFCellAccessor<dim, spacedim, level_dof_access> &cell,
4555 const unsigned int i,
4556 const std::integral_constant<int, 3>)
4557 {
4558 return cell.quad(i);
4559 }
4560 } // namespace DoFCellAccessorImplementation
4561} // namespace internal
4562
4563
4564
4565template <int dimension_, int space_dimension_, bool level_dof_access>
4566inline typename DoFCellAccessor<dimension_,
4567 space_dimension_,
4568 level_dof_access>::face_iterator
4570 const unsigned int i) const
4571{
4572 AssertIndexRange(i, this->n_faces());
4573
4574 return ::internal::DoFCellAccessorImplementation::get_face(
4575 *this, i, std::integral_constant<int, dimension_>());
4576}
4577
4578
4579
4580template <int dimension_, int space_dimension_, bool level_dof_access>
4581inline boost::container::small_vector<
4583 face_iterator,
4586 face_iterators() const
4587{
4588 boost::container::small_vector<
4590 face_iterator,
4592 face_iterators(this->n_faces());
4593
4594 for (const unsigned int i : this->face_indices())
4595 face_iterators[i] =
4597 *this, i, std::integral_constant<int, dimension_>());
4598
4599 return face_iterators;
4600}
4601
4602
4603
4604template <int dimension_, int space_dimension_, bool level_dof_access>
4605inline void
4608 std::vector<types::global_dof_index> &dof_indices) const
4609{
4610 if (level_dof_access)
4611 get_mg_dof_indices(dof_indices);
4612 else
4613 get_dof_indices(dof_indices);
4614}
4615
4616
4617
4618template <int dimension_, int space_dimension_, bool level_dof_access>
4619template <class InputVector, typename number>
4620inline void
4622 const InputVector &values,
4623 Vector<number> &local_values) const
4624{
4625 get_dof_values(values, local_values.begin(), local_values.end());
4626}
4627
4628
4629
4630template <int dimension_, int space_dimension_, bool level_dof_access>
4631template <typename Number, typename ForwardIterator>
4632inline void
4634 const ReadVector<Number> &values,
4635 ForwardIterator local_values_begin,
4636 ForwardIterator local_values_end) const
4637{
4638 (void)local_values_end;
4639 Assert(this->is_artificial() == false,
4640 ExcMessage("Can't ask for DoF indices on artificial cells."));
4641 Assert(this->is_active(), ExcMessage("Cell must be active."));
4642 Assert(this->dof_handler != nullptr, typename BaseClass::ExcInvalidObject());
4643
4644 Assert(static_cast<unsigned int>(local_values_end - local_values_begin) ==
4645 this->get_fe().n_dofs_per_cell(),
4647 Assert(values.size() == this->get_dof_handler().n_dofs(),
4649
4651 dof_indices(this->get_fe().n_dofs_per_cell());
4653 *this, dof_indices, this->active_fe_index());
4654
4655 boost::container::small_vector<Number, 27> values_temp(local_values_end -
4656 local_values_begin);
4657 auto view = make_array_view(values_temp.begin(), values_temp.end());
4658 values.extract_subvector_to(make_array_view(dof_indices.begin(),
4659 dof_indices.end()),
4660 view);
4661 using view_type = std::remove_reference_t<decltype(*local_values_begin)>;
4662 ArrayView<view_type> values_view2(&*local_values_begin,
4663 local_values_end - local_values_begin);
4664 std::copy(values_temp.begin(), values_temp.end(), values_view2.begin());
4665}
4666
4667
4668
4669template <int dimension_, int space_dimension_, bool level_dof_access>
4670template <class InputVector, typename ForwardIterator>
4671inline void
4674 const InputVector &values,
4675 ForwardIterator local_values_begin,
4676 ForwardIterator local_values_end) const
4677{
4678 Assert(this->is_artificial() == false,
4679 ExcMessage("Can't ask for DoF indices on artificial cells."));
4680 Assert(this->is_active(), ExcMessage("Cell must be active."));
4681
4682 Assert(static_cast<unsigned int>(local_values_end - local_values_begin) ==
4683 this->get_fe().n_dofs_per_cell(),
4685 Assert(values.size() == this->get_dof_handler().n_dofs(),
4687
4688
4690 dof_indices(this->get_fe().n_dofs_per_cell());
4692 *this, dof_indices, this->active_fe_index());
4693
4694 constraints.get_dof_values(values,
4695 dof_indices.data(),
4696 local_values_begin,
4697 local_values_end);
4698}
4699
4700
4701
4702template <int dimension_, int space_dimension_, bool level_dof_access>
4703template <class OutputVector, typename number>
4704inline void
4706 const Vector<number> &local_values,
4707 OutputVector &values) const
4708{
4709 Assert(this->is_artificial() == false,
4710 ExcMessage("Can't ask for DoF indices on artificial cells."));
4711 Assert(this->is_active(), ExcMessage("Cell must be active."));
4712
4713 Assert(static_cast<unsigned int>(local_values.size()) ==
4714 this->get_fe().n_dofs_per_cell(),
4716 Assert(values.size() == this->get_dof_handler().n_dofs(),
4718
4719
4720 Assert(this->dof_handler != nullptr, typename BaseClass::ExcInvalidObject());
4722 dof_indices(this->get_fe().n_dofs_per_cell());
4724 *this, dof_indices, this->active_fe_index());
4725
4726 for (unsigned int i = 0; i < this->get_fe().n_dofs_per_cell(); ++i)
4728 dof_indices[i],
4729 values);
4730}
4731
4732
4733
4734template <int dimension_, int space_dimension_, bool level_dof_access>
4737{
4738 Assert(this->dof_handler != nullptr, typename BaseClass::ExcInvalidObject());
4739 Assert((this->dof_handler->hp_capability_enabled == false) ||
4740 this->is_active(),
4741 ExcMessage(
4742 "For DoFHandler objects in hp-mode, finite elements are only "
4743 "associated with active cells. Consequently, you can not ask "
4744 "for the active finite element on cells with children."));
4745
4746 const auto &fe = this->dof_handler->get_fe(active_fe_index());
4747
4748 Assert(this->reference_cell() == fe.reference_cell(),
4749 internal::ExcNonMatchingReferenceCellTypes(this->reference_cell(),
4750 fe.reference_cell()));
4751
4752 return fe;
4753}
4754
4755
4756
4757template <int dimension_, int space_dimension_, bool level_dof_access>
4758inline types::fe_index
4760 active_fe_index() const
4761{
4762 Assert((this->dof_handler->hp_capability_enabled == false) ||
4763 this->is_active(),
4764 ExcMessage(
4765 "You can not ask for the active FE index on a cell that has "
4766 "children because no degrees of freedom are assigned "
4767 "to this cell and, consequently, no finite element "
4768 "is associated with it."));
4769 Assert((this->dof_handler->hp_capability_enabled == false) ||
4770 (this->is_locally_owned() || this->is_ghost()),
4771 ExcMessage("You can only query active FE index information on cells "
4772 "that are either locally owned or (after distributing "
4773 "degrees of freedom) are ghost cells."));
4774
4775 return ::internal::DoFCellAccessorImplementation::Implementation::
4776 active_fe_index(*this);
4777}
4778
4779
4780
4781template <int dimension_, int space_dimension_, bool level_dof_access>
4782inline void
4785{
4786 Assert((this->dof_handler->hp_capability_enabled == false) ||
4787 this->is_active(),
4788 ExcMessage("You can not set the active FE index on a cell that has "
4789 "children because no degrees of freedom will be assigned "
4790 "to this cell."));
4791
4792 Assert((this->dof_handler->hp_capability_enabled == false) ||
4793 this->is_locally_owned(),
4794 ExcMessage("You can only set active FE index information on cells "
4795 "that are locally owned. On ghost cells, this information "
4796 "will automatically be propagated from the owning process "
4797 "of that cell, and there is no information at all on "
4798 "artificial cells."));
4799
4801 set_active_fe_index(*this, i);
4802}
4803
4804
4805
4806template <int dimension_, int space_dimension_, bool level_dof_access>
4809 const
4810{
4811 Assert(this->dof_handler != nullptr, typename BaseClass::ExcInvalidObject());
4812 Assert((this->dof_handler->hp_capability_enabled == false) ||
4813 this->is_active(),
4814 ExcMessage(
4815 "For DoFHandler objects in hp-mode, finite elements are only "
4816 "associated with active cells. Consequently, you can not ask "
4817 "for the future finite element on cells with children."));
4818
4819 return this->dof_handler->get_fe(future_fe_index());
4820}
4821
4822
4823
4824template <int dimension_, int space_dimension_, bool level_dof_access>
4825inline types::fe_index
4827 future_fe_index() const
4828{
4829 Assert((this->dof_handler->hp_capability_enabled == false) ||
4830 (this->has_children() == false),
4831 ExcMessage(
4832 "You can not ask for the future FE index on a cell that has "
4833 "children because no degrees of freedom are assigned "
4834 "to this cell and, consequently, no finite element "
4835 "is associated with it."));
4836 Assert((this->dof_handler->hp_capability_enabled == false) ||
4837 (this->is_locally_owned()),
4838 ExcMessage("You can only query future FE index information on cells "
4839 "that are locally owned."));
4840
4841 return ::internal::DoFCellAccessorImplementation::Implementation::
4842 future_fe_index(*this);
4843}
4844
4845
4846
4847template <int dimension_, int space_dimension_, bool level_dof_access>
4848inline void
4851{
4852 Assert((this->dof_handler->hp_capability_enabled == false) ||
4853 (this->has_children() == false),
4854 ExcMessage("You can not set the future FE index on a cell that has "
4855 "children because no degrees of freedom will be assigned "
4856 "to this cell."));
4857
4858 Assert((this->dof_handler->hp_capability_enabled == false) ||
4859 this->is_locally_owned(),
4860 ExcMessage("You can only set future FE index information on cells "
4861 "that are locally owned."));
4862
4864 set_future_fe_index(*this, i);
4865}
4866
4867
4868
4869template <int dimension_, int space_dimension_, bool level_dof_access>
4870inline bool
4873{
4874 Assert((this->dof_handler->hp_capability_enabled == false) ||
4875 (this->has_children() == false),
4876 ExcMessage(
4877 "You can not ask for the future FE index on a cell that has "
4878 "children because no degrees of freedom are assigned "
4879 "to this cell and, consequently, no finite element "
4880 "is associated with it."));
4881 Assert((this->dof_handler->hp_capability_enabled == false) ||
4882 (this->is_locally_owned()),
4883 ExcMessage("You can only query future FE index information on cells "
4884 "that are locally owned."));
4885
4886 return ::internal::DoFCellAccessorImplementation::Implementation::
4887 future_fe_index_set(*this);
4888}
4889
4890
4891
4892template <int dimension_, int space_dimension_, bool level_dof_access>
4893inline void
4896{
4897 Assert((this->dof_handler->hp_capability_enabled == false) ||
4898 (this->has_children() == false),
4899 ExcMessage(
4900 "You can not ask for the future FE index on a cell that has "
4901 "children because no degrees of freedom are assigned "
4902 "to this cell and, consequently, no finite element "
4903 "is associated with it."));
4904 Assert((this->dof_handler->hp_capability_enabled == false) ||
4905 (this->is_locally_owned()),
4906 ExcMessage("You can only query future FE index information on cells "
4907 "that are locally owned."));
4908
4911}
4912
4913
4914
4915template <int dimension_, int space_dimension_, bool level_dof_access>
4916template <typename number, typename OutputVector>
4917inline void
4920 OutputVector &global_destination) const
4921{
4922 this->distribute_local_to_global(local_source.begin(),
4923 local_source.end(),
4924 global_destination);
4925}
4926
4927
4928
4929template <int dimension_, int space_dimension_, bool level_dof_access>
4930template <typename ForwardIterator, typename OutputVector>
4931inline void
4933 distribute_local_to_global(ForwardIterator local_source_begin,
4934 ForwardIterator local_source_end,
4935 OutputVector &global_destination) const
4936{
4937 Assert(this->dof_handler != nullptr,
4938 (typename std::decay_t<decltype(*this)>::ExcInvalidObject()));
4939 Assert(static_cast<unsigned int>(local_source_end - local_source_begin) ==
4940 this->get_fe().n_dofs_per_cell(),
4941 (typename std::decay_t<decltype(*this)>::ExcVectorDoesNotMatch()));
4942 Assert(this->dof_handler->n_dofs() == global_destination.size(),
4943 (typename std::decay_t<decltype(*this)>::ExcVectorDoesNotMatch()));
4944
4945 Assert(!this->has_children(), ExcMessage("Cell must be active"));
4946
4947 Assert(
4948 internal::ArrayViewHelper::is_contiguous(local_source_begin,
4949 local_source_end),
4950 ExcMessage(
4951 "This function can not be called with iterator types that do not point to contiguous memory."));
4952
4953 const unsigned int n_dofs = local_source_end - local_source_begin;
4954
4956 dof_indices(n_dofs);
4958 *this, dof_indices, this->active_fe_index());
4959
4960 // distribute cell vector
4961 global_destination.add(n_dofs, dof_indices.data(), &(*local_source_begin));
4962}
4963
4964
4965
4966template <int dimension_, int space_dimension_, bool level_dof_access>
4967template <typename ForwardIterator, typename OutputVector>
4968inline void
4972 ForwardIterator local_source_begin,
4973 ForwardIterator local_source_end,
4974 OutputVector &global_destination) const
4975{
4976 Assert(this->dof_handler != nullptr,
4977 (typename std::decay_t<decltype(*this)>::ExcInvalidObject()));
4978 Assert(local_source_end - local_source_begin ==
4979 this->get_fe().n_dofs_per_cell(),
4980 (typename std::decay_t<decltype(*this)>::ExcVectorDoesNotMatch()));
4981 Assert(this->dof_handler->n_dofs() == global_destination.size(),
4982 (typename std::decay_t<decltype(*this)>::ExcVectorDoesNotMatch()));
4983
4984 Assert(!this->has_children(), ExcMessage("Cell must be active."));
4985
4987 dof_indices(this->get_fe().n_dofs_per_cell());
4989 *this, dof_indices, this->active_fe_index());
4990
4991 // distribute cell vector
4992 constraints.distribute_local_to_global(local_source_begin,
4993 local_source_end,
4994 dof_indices.data(),
4995 global_destination);
4996}
4997
4998
4999
5000template <int dimension_, int space_dimension_, bool level_dof_access>
5001template <typename number, typename OutputMatrix>
5002inline void
5005 OutputMatrix &global_destination) const
5006{
5007 Assert(this->dof_handler != nullptr,
5008 (typename std::decay_t<decltype(*this)>::ExcInvalidObject()));
5009 Assert(local_source.m() == this->get_fe().n_dofs_per_cell(),
5010 (typename std::decay_t<decltype(*this)>::ExcMatrixDoesNotMatch()));
5011 Assert(local_source.n() == this->get_fe().n_dofs_per_cell(),
5012 (typename std::decay_t<decltype(*this)>::ExcMatrixDoesNotMatch()));
5013 Assert(this->dof_handler->n_dofs() == global_destination.m(),
5014 (typename std::decay_t<decltype(*this)>::ExcMatrixDoesNotMatch()));
5015 Assert(this->dof_handler->n_dofs() == global_destination.n(),
5016 (typename std::decay_t<decltype(*this)>::ExcMatrixDoesNotMatch()));
5017
5018 Assert(!this->has_children(), ExcMessage("Cell must be active."));
5019
5020 const unsigned int n_dofs = local_source.m();
5021
5023 dof_indices(n_dofs);
5025 *this, dof_indices, this->active_fe_index());
5026
5027 // distribute cell matrix
5028 for (unsigned int i = 0; i < n_dofs; ++i)
5029 global_destination.add(dof_indices[i],
5030 n_dofs,
5031 dof_indices.data(),
5032 &local_source(i, 0));
5033}
5034
5035
5036
5037template <int dimension_, int space_dimension_, bool level_dof_access>
5038template <typename number, typename OutputMatrix, typename OutputVector>
5039inline void
5042 const Vector<number> &local_vector,
5043 OutputMatrix &global_matrix,
5044 OutputVector &global_vector) const
5045{
5046 Assert(this->dof_handler != nullptr,
5047 (typename std::decay_t<decltype(*this)>::ExcInvalidObject()));
5048 Assert(local_matrix.m() == this->get_fe().n_dofs_per_cell(),
5049 (typename std::decay_t<decltype(*this)>::ExcMatrixDoesNotMatch()));
5050 Assert(local_matrix.n() == this->get_fe().n_dofs_per_cell(),
5051 (typename std::decay_t<decltype(*this)>::ExcVectorDoesNotMatch()));
5052 Assert(this->dof_handler->n_dofs() == global_matrix.m(),
5053 (typename std::decay_t<decltype(*this)>::ExcMatrixDoesNotMatch()));
5054 Assert(this->dof_handler->n_dofs() == global_matrix.n(),
5055 (typename std::decay_t<decltype(*this)>::ExcMatrixDoesNotMatch()));
5056 Assert(local_vector.size() == this->get_fe().n_dofs_per_cell(),
5057 (typename std::decay_t<decltype(*this)>::ExcVectorDoesNotMatch()));
5058 Assert(this->dof_handler->n_dofs() == global_vector.size(),
5059 (typename std::decay_t<decltype(*this)>::ExcVectorDoesNotMatch()));
5060
5061 Assert(!this->has_children(), ExcMessage("Cell must be active."));
5062
5063 const unsigned int n_dofs = this->get_fe().n_dofs_per_cell();
5065 dof_indices(n_dofs);
5067 *this, dof_indices, this->active_fe_index());
5068
5069 // distribute cell matrices
5070 for (unsigned int i = 0; i < n_dofs; ++i)
5071 {
5072 global_matrix.add(dof_indices[i],
5073 n_dofs,
5074 dof_indices.data(),
5075 &local_matrix(i, 0));
5076 global_vector(dof_indices[i]) += local_vector(i);
5077 }
5078}
5079
5080
5081
5083
5084
5085#endif
ArrayView< std::remove_reference_t< typename std::iterator_traits< Iterator >::reference >, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition array_view.h:954
void get_dof_values(const VectorType &global_vector, ForwardIteratorInd local_indices_begin, ForwardIteratorVec local_vector_begin, ForwardIteratorVec local_vector_end) const
void distribute_local_to_global(const InVector &local_vector, const std::vector< size_type > &local_dof_indices, OutVector &global_vector) const
iterator begin() const
Definition array_view.h:707
friend struct ::internal::DoFHandlerImplementation::Policy::Implementation
void get_mg_dof_indices(const int level, std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index) const
DoFAccessor< 0, 1, spacedim, level_dof_access > & operator=(const DoFAccessor< 0, 1, spacedim, level_dof_access > &da)=delete
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index) const
typename::internal::DoFHandlerImplementation::Iterators< 1, spacedim, level_dof_access >::line_iterator line(const unsigned int i) const
DoFAccessor< 0, 1, spacedim, level_dof_access > & operator=(DoFAccessor< 0, 1, spacedim, level_dof_access > &&) noexcept=default
types::global_dof_index vertex_dof_index(const unsigned int vertex, const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
TriaAccessor< 0, 1, spacedim > BaseClass
DoFAccessor(DoFAccessor< 0, 1, spacedim, level_dof_access > &&)=default
const DoFHandler< 1, spacedim > & get_dof_handler() const
types::global_dof_index dof_index(const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
typename::internal::DoFHandlerImplementation::Iterators< 1, spacedim, level_dof_access >::quad_iterator quad(const unsigned int i) const
DoFAccessor(const DoFAccessor< 0, 1, spacedim, level_dof_access > &)=default
types::fe_index nth_active_fe_index(const unsigned int n) const
const FiniteElement< 1, spacedim > & get_fe(const types::fe_index fe_index) const
void set_dof_index(const unsigned int i, const types::global_dof_index index, const types::fe_index fe_index=numbers::invalid_fe_index) const
static constexpr unsigned int space_dimension
void set_dof_handler(DoFHandler< 1, spacedim > *dh)
TriaIterator< DoFAccessor< 0, 1, spacedim, level_dof_access > > child(const unsigned int c) const
bool fe_index_is_active(const types::fe_index fe_index) const
void copy_from(const DoFAccessor< 0, 1, spacedim, level_dof_access2 > &a)
typename::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::line_iterator line(const unsigned int i) const
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
friend class DoFAccessor
DoFAccessor< structdim, dim, spacedim, level_dof_access > & operator=(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &da)=delete
types::fe_index nth_active_fe_index(const unsigned int n) const
DoFAccessor(DoFAccessor< structdim, dim, spacedim, level_dof_access > &&)
DoFAccessor(const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &)
TriaIterator< DoFAccessor< structdim, dim, spacedim, level_dof_access > > child(const unsigned int c) 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 types::fe_index fe_index=numbers::invalid_fe_index) const
~DoFAccessor()=default
void set_dof_handler(DoFHandler< dim, spacedim > *dh)
DoFAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
DoFAccessor< structdim, dim, spacedim, level_dof_access > & operator=(DoFAccessor< structdim, dim, spacedim, level_dof_access > &&)
std::set< types::fe_index > get_active_fe_indices() const
DoFAccessor(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &)=default
void set_mg_dof_indices(const int level, const std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index)
void get_mg_dof_indices(const int level, std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index) const
types::global_dof_index mg_vertex_dof_index(const int level, const unsigned int vertex, const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
typename ::internal::DoFAccessorImplementation:: Inheritance< structdim, dimension, space_dimension >::BaseClass BaseClass
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 copy_from(const TriaAccessorBase< structdim, dim, spacedim > &da)
types::global_dof_index dof_index(const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
typename::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::quad_iterator quad(const unsigned int i) const
DoFHandler< dimension, space_dimension > AccessorData
void copy_from(const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &a)
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index fe_index) const
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::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
void set_dof_index(const unsigned int i, const types::global_dof_index index, const types::fe_index fe_index=numbers::invalid_fe_index) const
unsigned int n_active_fe_indices() const
static bool is_level_cell()
types::global_dof_index vertex_dof_index(const unsigned int vertex, const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
DoFAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
bool fe_index_is_active(const types::fe_index fe_index) const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > periodic_neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) 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
DoFHandler< dimension_, space_dimension_ > AccessorData
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
boost::container::small_vector< face_iterator, GeometryInfo< dimension_ >::faces_per_cell > face_iterators() const
types::fe_index future_fe_index() 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
DoFAccessor< dimension_, dimension_, space_dimension_, level_dof_access > BaseClass
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > child(const unsigned int i) const
TriaIterator< DoFAccessor< dimension_ - 1, dimension_, space_dimension_, level_dof_access > > face_iterator
DoFCellAccessor(const Triangulation< dimension_, space_dimension_ > *tria, const int level, const int index, const AccessorData *local_data)
void set_future_fe_index(const types::fe_index i) const
void distribute_local_to_global(const Vector< number > &local_source, OutputVector &global_destination) const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > neighbor_or_periodic_neighbor(const unsigned int i) 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
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
static const unsigned int dim
void set_active_fe_index(const types::fe_index i) const
void distribute_local_to_global_by_interpolation(const Vector< number > &local_values, OutputVector &values, const types::fe_index fe_index=numbers::invalid_fe_index) const
DoFCellAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > periodic_neighbor(const unsigned int i) const
DoFCellAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
static const unsigned int spacedim
void get_interpolated_dof_values(const ReadVector< Number > &values, Vector< Number > &interpolated_values, const types::fe_index fe_index=numbers::invalid_fe_index) const
face_iterator face(const unsigned int i) const
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 set_dof_indices(const std::vector< types::global_dof_index > &dof_indices)
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_dof_values(const Vector< number > &local_values, OutputVector &values) const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
void set_dof_values_by_interpolation(const Vector< number > &local_values, OutputVector &values, const types::fe_index fe_index=numbers::invalid_fe_index, const bool perform_check=false) 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
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
void get_dof_values(const ReadVector< Number > &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) const
~DoFCellAccessor()=default
DoFHandler< dimension_, space_dimension_ > Container
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
types::fe_index active_fe_index() const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > parent() const
std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel< dim > > > mg_levels
static const types::fe_index default_fe_index
std::vector< std::vector< types::fe_index > > hp_cell_future_fe_indices
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
std::vector< std::vector< types::fe_index > > hp_cell_active_fe_indices
std::vector< MGVertexDoFs > mg_vertex_dofs
std::vector< std::array< std::vector< types::global_dof_index >, dim+1 > > object_dof_indices
std::vector< std::array< std::vector< offset_type >, dim+1 > > object_dof_ptr
std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces< dim > > mg_faces
bool hp_capability_enabled
bool has_hp_capabilities() const
std::array< std::vector< offset_type >, dim+1 > hp_object_fe_ptr
std::array< std::vector< types::fe_index >, dim+1 > hp_object_fe_indices
typename InvalidAccessor< structdim, dim, spacedim >::AccessorData AccessorData
DoFInvalidAccessor(const void *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
unsigned int n_dofs_per_vertex() const
size_type n() const
size_type m() const
void compress() const
Definition index_set.h:1799
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Definition index_set.h:1846
InvalidAccessor(const void *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
void import_elements(const ::Vector< Number > &vec, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
TriaAccessorBase(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *=nullptr)
const Triangulation< dim, spacedim > * tria
TriaIterator< TriaAccessor< structdim, dim, spacedim > > child(const unsigned int i) const
unsigned int vertex_index(const unsigned int i) const
Point< spacedim > & vertex(const unsigned int i) const
TriaAccessor(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
IteratorState::IteratorStates state() const
virtual size_type size() const override
iterator end()
iterator begin()
#define DEAL_II_ALWAYS_INLINE
Definition config.h:166
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_RESTRICT
Definition config.h:167
constexpr bool running_in_debug_mode()
Definition config.h:78
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
#define DeclException0(Exception0)
static ::ExceptionBase & ExcVectorDoesNotMatch()
static ::ExceptionBase & ExcNonMatchingReferenceCellTypes(ReferenceCell arg1, ReferenceCell arg2)
static ::ExceptionBase & ExcVectorNotEmpty()
static ::ExceptionBase & ExcMatrixDoesNotMatch()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
static ::ExceptionBase & ExcNotImplementedWithHP()
static ::ExceptionBase & ExcCantCompareIterators()
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
@ past_the_end
Iterator reached end of container.
Definition hp.h:117
types::fe_index get_fe_index_or_default(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &cell, const types::fe_index fe_index)
void get_cell_dof_indices(const ::DoFCellAccessor< dim, spacedim, level_dof_access > &accessor, Implementation::dof_index_vector_type &dof_indices, const unsigned int fe_index)
TriaIterator< ::DoFAccessor< dim - 1, dim, spacedim, level_dof_access > > get_face(const ::DoFCellAccessor< dim, spacedim, level_dof_access > &cell, const unsigned int i, const std::integral_constant< int, 1 >)
constexpr types::global_dof_index invalid_dof_index
Definition types.h:269
constexpr types::geometric_orientation reverse_line_orientation
Definition types.h:365
constexpr types::geometric_orientation default_geometric_orientation
Definition types.h:352
constexpr types::fe_index invalid_fe_index
Definition types.h:260
STL namespace.
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
Definition types.h:32
unsigned short int fe_index
Definition types.h:72
unsigned int global_dof_index
Definition types.h:94
static constexpr unsigned int faces_per_cell
static constexpr unsigned int max_children_per_cell
void process_vertex_dofs(DoFHandler< dim, spacedim > &dof_handler, const unsigned int vertex_index, const types::fe_index fe_index, types::global_dof_index *&dof_indices_ptr, const DoFProcessor &dof_processor) const
void process_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const types::fe_index fe_index, const DoFMapping &mapping, const std::integral_constant< int, structdim >, types::global_dof_index *&dof_indices_ptr, const DoFProcessor &dof_processor) const
void process_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int, const unsigned int obj_index, const types::fe_index fe_index, const DoFMapping &mapping, const std::integral_constant< int, structdim >, types::global_dof_index *&dof_indices_ptr, const DoFProcessor &dof_processor) const
void process_vertex_dofs(DoFHandler< dim, spacedim > &dof_handler, const unsigned int vertex_index, const types::fe_index, types::global_dof_index *&dof_indices_ptr, const DoFProcessor &dof_processor) const
static types::global_dof_index * get_array_ptr(const std::tuple<> &)
static void process_dof_indices(const ::DoFAccessor< structdim, dim, spacedim, level_dof_access > &accessor, const DoFIndicesType &const_dof_indices, const types::fe_index fe_index_, const DoFOperation &dof_operation, const DoFProcessor &dof_processor, const bool count_level_dofs)
static void get_dof_indices(const ::DoFAccessor< structdim, dim, spacedim, level_dof_access > &accessor, std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index)
static void extract_subvector_to(const LinearAlgebra::TpetraWrappers::Vector< Number, MemorySpace > &values, const types::global_dof_index *cache_begin, const types::global_dof_index *cache_end, ForwardIterator local_values_begin)
static std::set< types::fe_index > get_active_fe_indices(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const std::integral_constant< int, structdim > &t)
static types::global_dof_index * get_array_ptr(const ArrayType &array)
static unsigned int get_array_length(const std::tuple<> &)
static unsigned int n_active_fe_indices(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const std::integral_constant< int, structdim > &)
static void set_dof_index(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const types::fe_index fe_index, const unsigned int local_index, const std::integral_constant< int, structdim > &dd, const types::global_dof_index global_index)
static void extract_subvector_to(const LinearAlgebra::EpetraWrappers::Vector &values, const types::global_dof_index *cache_begin, const types::global_dof_index *cache_end, ForwardIterator local_values_begin)
static void get_mg_dof_indices(const ::DoFAccessor< structdim, dim, spacedim, level_dof_access > &accessor, const int level, std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index)
static bool fe_index_is_active(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const types::fe_index fe_index, const std::integral_constant< int, structdim > &)
static types::global_dof_index & get_mg_dof_index(const DoFHandler< dim, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< dim > > &, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< dim > > &mg_faces, const unsigned int obj_index, const types::fe_index fe_index, const unsigned int local_index, const std::integral_constant< int, 1 >)
static types::global_dof_index & get_mg_dof_index(const DoFHandler< dim, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< dim > > &mg_level, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< dim > > &, const unsigned int obj_index, const types::fe_index fe_index, const unsigned int local_index, const std::integral_constant< int, dim >)
static std::pair< unsigned int, unsigned int > process_object_range(const ::DoFAccessor< structdim, dim, spacedim, level_dof_access > accessor, const types::fe_index fe_index)
static types::global_dof_index & mg_vertex_dof_index(DoFHandler< dim, spacedim > &dof_handler, const int level, const unsigned int vertex_index, const unsigned int i)
boost::container::small_vector<::types::global_dof_index, 27 > dof_index_vector_type
static void set_mg_dof_indices(const ::DoFAccessor< structdim, dim, spacedim, level_dof_access > &accessor, const int level, const std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index)
static void set_dof_indices(const ::DoFAccessor< structdim, dim, spacedim, level_dof_access > &accessor, const std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index)
static void process_dof_index(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const types::fe_index fe_index, const unsigned int local_index, const std::integral_constant< int, structdim > &, GlobalIndexType &global_index, const DoFPProcessor &process)
static types::global_dof_index & get_mg_dof_index(const DoFHandler< 3, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< 3 > > &, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< 3 > > &mg_faces, const unsigned int obj_index, const types::fe_index fe_index, const unsigned int local_index, const std::integral_constant< int, 2 >)
static std::pair< unsigned int, unsigned int > process_object_range(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const types::fe_index fe_index, const std::integral_constant< int, structdim > &)
static void extract_subvector_to(const InputVector &values, const types::global_dof_index *cache, const types::global_dof_index *cache_end, ForwardIterator local_values_begin)
static std::vector< unsigned int > sort_indices(const types::global_dof_index *v_begin, const types::global_dof_index *v_end)
static types::fe_index nth_active_fe_index(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const unsigned int local_index, const std::integral_constant< int, structdim > &)
static void process_object(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const types::fe_index fe_index, const DoFMapping &mapping, const std::integral_constant< int, structdim > &dd, types::global_dof_index *&dof_indices_ptr, const DoFProcessor &process)
static types::global_dof_index get_dof_index(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const types::fe_index fe_index, const unsigned int local_index, const std::integral_constant< int, structdim > &dd)
static unsigned int n_dof_indices(const ::DoFAccessor< structdim, dim, spacedim, level_dof_access > &accessor, const types::fe_index fe_index_, const bool count_level_dofs)
static std::pair< unsigned int, unsigned int > process_object_range(::DoFInvalidAccessor< structdim, dim, spacedim >, const unsigned int)
static unsigned int get_array_length(const ArrayType &array)
::TriaAccessor< structdim, dim, spacedim > BaseClass
static types::fe_index future_fe_index(const DoFCellAccessor< dim, spacedim, level_dof_access > &accessor)
static void clear_future_fe_index(const DoFCellAccessor< dim, spacedim, level_dof_access > &accessor)
static bool future_fe_index_set(const DoFCellAccessor< dim, spacedim, level_dof_access > &accessor)
static void set_active_fe_index(const DoFCellAccessor< dim, spacedim, level_dof_access > &accessor, const types::fe_index i)
static void set_future_fe_index(const DoFCellAccessor< dim, spacedim, level_dof_access > &accessor, const types::fe_index i)
static types::fe_index active_fe_index(const DoFCellAccessor< dim, spacedim, level_dof_access > &accessor)
static void set(typename VectorType::value_type value, const types::global_dof_index i, VectorType &V)