Reference documentation for deal.II version 9.6.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_handler.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1998 - 2024 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_handler_h
16#define dealii_dof_handler_h
17
18
19
20#include <deal.II/base/config.h>
21
27#include <deal.II/base/types.h>
28
35
37
38#include <boost/serialization/split_member.hpp>
39
40#include <map>
41#include <memory>
42#include <set>
43#include <vector>
44
46
47// Forward declarations
48#ifndef DOXYGEN
49template <int dim, int spacedim>
50class FiniteElement;
51template <int dim, int spacedim>
53class Triangulation;
54
55namespace internal
56{
57 namespace DoFHandlerImplementation
58 {
59 struct Implementation;
60
61 namespace Policy
62 {
63 template <int dim, int spacedim>
64 class PolicyBase;
65 struct Implementation;
66 } // namespace Policy
67 } // namespace DoFHandlerImplementation
68
69 namespace DoFAccessorImplementation
70 {
71 struct Implementation;
72 }
73
74 namespace DoFCellAccessorImplementation
75 {
76 struct Implementation;
77 }
78
79 namespace hp
80 {
81 namespace DoFHandlerImplementation
82 {
83 struct Implementation;
84 }
85 } // namespace hp
86} // namespace internal
87
88namespace parallel
89{
90 namespace distributed
91 {
92 template <int dim, int spacedim, typename VectorType>
93 class CellDataTransfer;
94 }
95} // namespace parallel
96#endif
97
314template <int dim, int spacedim = dim>
317{
322
323public:
336 using cell_accessor = typename ActiveSelector::CellAccessor;
337
350 using face_accessor = typename ActiveSelector::FaceAccessor;
351
359 using line_iterator = typename ActiveSelector::line_iterator;
360
374 using active_line_iterator = typename ActiveSelector::active_line_iterator;
375
383 using quad_iterator = typename ActiveSelector::quad_iterator;
384
398 using active_quad_iterator = typename ActiveSelector::active_quad_iterator;
399
407 using hex_iterator = typename ActiveSelector::hex_iterator;
408
418 using active_hex_iterator = typename ActiveSelector::active_hex_iterator;
419
440 using active_cell_iterator = typename ActiveSelector::active_cell_iterator;
441
468 using cell_iterator = typename ActiveSelector::cell_iterator;
469
486 using face_iterator = typename ActiveSelector::face_iterator;
487
499 using active_face_iterator = typename ActiveSelector::active_face_iterator;
500
501 using level_cell_accessor = typename LevelSelector::CellAccessor;
502 using level_face_accessor = typename LevelSelector::FaceAccessor;
503
504 using level_cell_iterator = typename LevelSelector::cell_iterator;
505 using level_face_iterator = typename LevelSelector::face_iterator;
506
507
511 static constexpr unsigned int dimension = dim;
512
516 static constexpr unsigned int space_dimension = spacedim;
517
521 static const types::fe_index default_fe_index = 0;
522
528 static const unsigned int invalid_fe_index DEAL_II_DEPRECATED =
530
537
541 using offset_type = unsigned int;
542
549 static const types::fe_index invalid_active_fe_index DEAL_II_DEPRECATED =
551
557
562
569 DoFHandler(const DoFHandler &) = delete;
570
574 virtual ~DoFHandler() override;
575
582 DoFHandler &
583 operator=(const DoFHandler &) = delete;
584
599 void
600 set_active_fe_indices(const std::vector<types::fe_index> &active_fe_indices);
601
608 void
609 set_active_fe_indices(const std::vector<unsigned int> &active_fe_indices);
610
624 std::vector<types::fe_index>
626
644 void
645 get_active_fe_indices(std::vector<unsigned int> &active_fe_indices) const;
646
660 void
661 set_future_fe_indices(const std::vector<types::fe_index> &future_fe_indices);
662
673 std::vector<types::fe_index>
675
683 void
685
715 void
717
721 void
723
729 void
731
735 bool
737
743 bool
745
760 bool
762
769 void
771
775 void
777
830 void
831 renumber_dofs(const std::vector<types::global_dof_index> &new_numbers);
832
837 void
838 renumber_dofs(const unsigned int level,
839 const std::vector<types::global_dof_index> &new_numbers);
840
869 unsigned int
871
884 unsigned int
886
887 /*--------------------------------------*/
888
899 begin(const unsigned int level = 0) const;
900
918 begin_active(const unsigned int level = 0) const;
919
925 end() const;
926
932 end(const unsigned int level) const;
933
940 end_active(const unsigned int level) const;
941
947 begin_mg(const unsigned int level = 0) const;
948
954 end_mg(const unsigned int level) const;
955
961 end_mg() const;
981
1024
1038
1055 cell_iterators_on_level(const unsigned int level) const;
1056
1073 active_cell_iterators_on_level(const unsigned int level) const;
1074
1091 mg_cell_iterators_on_level(const unsigned int level) const;
1092
1096 /*---------------------------------------*/
1097
1098
1123 n_dofs() const;
1124
1133 n_dofs(const unsigned int level) const;
1134
1141
1153 template <typename number>
1156 const std::map<types::boundary_id, const Function<spacedim, number> *>
1157 &boundary_ids) const;
1158
1164 n_boundary_dofs(const std::set<types::boundary_id> &boundary_ids) const;
1165
1181 const BlockInfo &
1182 block_info() const;
1183
1205
1211 const IndexSet &
1213
1218 const IndexSet &
1219 locally_owned_mg_dofs(const unsigned int level) const;
1220
1226 get_fe(const types::fe_index index = 0) const;
1227
1234
1240
1244 MPI_Comm
1246
1263 void
1265
1281 void
1283
1292 virtual std::size_t
1294
1300 template <class Archive>
1301 void
1302 save(Archive &ar, const unsigned int version) const;
1303
1309 template <class Archive>
1310 void
1311 load(Archive &ar, const unsigned int version);
1312
1313#ifdef DOXYGEN
1319 template <class Archive>
1320 void
1321 serialize(Archive &archive, const unsigned int version);
1322#else
1323 // This macro defines the serialize() method that is compatible with
1324 // the templated save() and load() method that have been implemented.
1325 BOOST_SERIALIZATION_SPLIT_MEMBER()
1326#endif
1327
1331 DeclException0(ExcNoFESelected);
1336 DeclException0(ExcInvalidBoundaryIndicator);
1341 DeclException1(ExcInvalidLevel,
1342 int,
1343 << "The given level " << arg1
1344 << " is not in the valid range!");
1349 DeclException1(ExcNewNumbersNotConsecutive,
1351 << "The given list of new dof indices is not consecutive: "
1352 << "the index " << arg1 << " does not exist.");
1356 DeclException2(ExcInvalidFEIndex,
1357 int,
1358 int,
1359 << "The mesh contains a cell with an active FE index of "
1360 << arg1 << ", but the finite element collection only has "
1361 << arg2 << " elements");
1362
1367 DeclExceptionMsg(ExcOnlyAvailableWithHP,
1368 "The current function doesn't make sense when used with a "
1369 "DoFHandler without hp-capabilities.");
1370
1375 DeclExceptionMsg(ExcNotImplementedWithHP,
1376 "The current function has not yet been implemented for a "
1377 "DoFHandler with hp-capabilities.");
1378
1379private:
1387 {
1388 public:
1393
1399 void
1400 init(const unsigned int coarsest_level,
1401 const unsigned int finest_level,
1402 const unsigned int dofs_per_vertex);
1403
1407 unsigned int
1409
1413 unsigned int
1415
1421 access_index(const unsigned int level,
1422 const unsigned int dof_number,
1423 const unsigned int dofs_per_vertex);
1424
1425 private:
1429 unsigned int coarsest_level;
1430
1434 unsigned int finest_level;
1435
1445 std::unique_ptr<types::global_dof_index[]> indices;
1446 };
1447
1455 {
1460 std::map<const cell_iterator, const types::fe_index>
1462
1467 std::map<const cell_iterator, const types::fe_index> refined_cells_fe_index;
1468
1473 std::map<const cell_iterator, const types::fe_index>
1475
1481 std::vector<types::fe_index> active_fe_indices;
1482
1488 std::unique_ptr<
1489 parallel::distributed::
1490 CellDataTransfer<dim, spacedim, std::vector<types::fe_index>>>
1492 };
1493
1498
1504
1510
1517
1522 std::unique_ptr<::internal::DoFHandlerImplementation::Policy::
1523 PolicyBase<dim, spacedim>>
1525
1534
1538 std::vector<::internal::DoFHandlerImplementation::NumberCache>
1540
1546 mutable std::vector<std::array<std::vector<types::global_dof_index>, dim + 1>>
1548
1557 mutable std::vector<std::array<std::vector<offset_type>, dim + 1>>
1559
1565 mutable std::array<std::vector<types::fe_index>, dim + 1>
1567
1571 mutable std::array<std::vector<offset_type>, dim + 1> hp_object_fe_ptr;
1572
1577 mutable std::vector<std::vector<types::fe_index>> hp_cell_active_fe_indices;
1578
1583 mutable std::vector<std::vector<types::fe_index>> hp_cell_future_fe_indices;
1584
1589 std::vector<MGVertexDoFs> mg_vertex_dofs;
1590
1594 std::vector<
1595 std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel<dim>>>
1597
1601 std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces<dim>>
1603
1608 std::unique_ptr<ActiveFEIndexTransfer> active_fe_index_transfer;
1609
1614 std::vector<boost::signals2::connection> tria_listeners;
1615
1621 std::vector<boost::signals2::connection> tria_listeners_for_transfer;
1622
1626 void
1628
1632 void
1634
1638 void
1640
1644 void
1646
1658 void
1660
1673 void
1675
1685 void
1687
1696 void
1698
1707 void
1709
1718 void
1720
1721
1722 // Make accessor objects friends.
1723 template <int, int, int, bool>
1724 friend class ::DoFAccessor;
1725 template <int, int, bool>
1726 friend class ::DoFCellAccessor;
1727 friend struct ::internal::DoFAccessorImplementation::Implementation;
1728 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
1729
1730 // Likewise for DoFLevel objects since they need to access the vertex dofs
1731 // in the functions that set and retrieve vertex dof indices.
1732 friend struct ::internal::DoFHandlerImplementation::Implementation;
1733 friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
1734 friend struct ::internal::DoFHandlerImplementation::Policy::
1735 Implementation;
1736
1737 // explicitly check for sensible template arguments, but not on windows
1738 // because MSVC creates bogus warnings during normal compilation
1739#ifndef DEAL_II_MSVC
1740 static_assert(dim <= spacedim,
1741 "The dimension <dim> of a DoFHandler must be less than or "
1742 "equal to the space dimension <spacedim> in which it lives.");
1743#endif
1744};
1745
1746namespace internal
1747{
1748 namespace hp
1749 {
1750 namespace DoFHandlerImplementation
1751 {
1763 template <int dim, int spacedim>
1764 void
1766
1791 template <int dim, int spacedim = dim>
1792 unsigned int
1794 const typename DoFHandler<dim, spacedim>::cell_iterator &parent);
1795
1801 "No FiniteElement has been found in your FECollection that is "
1802 "dominated by all children of a cell you are trying to coarsen!");
1803 } // namespace DoFHandlerImplementation
1804 } // namespace hp
1805} // namespace internal
1806
1807#ifndef DOXYGEN
1808
1809/* ----------------------- Inline functions ----------------------------------
1810 */
1811
1812
1813template <int dim, int spacedim>
1816{
1817 return hp_capability_enabled;
1818}
1819
1820
1821
1822template <int dim, int spacedim>
1825{
1826 return mg_number_cache.size() > 0;
1827}
1828
1829
1830
1831template <int dim, int spacedim>
1834{
1835 return number_cache.n_global_dofs > 0;
1836}
1837
1838
1839
1840template <int dim, int spacedim>
1843{
1844 return number_cache.n_global_dofs;
1845}
1846
1847
1848
1849template <int dim, int spacedim>
1852 DoFHandler<dim, spacedim>::n_dofs(const unsigned int level) const
1853{
1854 Assert(has_level_dofs(),
1855 ExcMessage(
1856 "n_dofs(level) can only be called after distribute_mg_dofs()"));
1857 Assert(level < mg_number_cache.size(), ExcInvalidLevel(level));
1858 return mg_number_cache[level].n_global_dofs;
1859}
1860
1861
1862
1863template <int dim, int spacedim>
1866{
1867 return number_cache.n_locally_owned_dofs;
1868}
1869
1870
1871
1872template <int dim, int spacedim>
1875{
1876 return number_cache.locally_owned_dofs;
1877}
1878
1879
1880
1881template <int dim, int spacedim>
1884 const unsigned int level) const
1885{
1886 Assert(level < this->get_triangulation().n_global_levels(),
1887 ExcMessage("The given level index exceeds the number of levels "
1888 "present in the triangulation"));
1889 Assert(
1890 mg_number_cache.size() == this->get_triangulation().n_global_levels(),
1891 ExcMessage(
1892 "The level dofs are not set up properly! Did you call distribute_mg_dofs()?"));
1893 return mg_number_cache[level].locally_owned_dofs;
1894}
1895
1896
1897
1898template <int dim, int spacedim>
1901 const types::fe_index number) const
1902{
1903 Assert(fe_collection.size() > 0,
1904 ExcMessage("No finite element collection is associated with "
1905 "this DoFHandler"));
1906 return fe_collection[number];
1907}
1908
1909
1910
1911template <int dim, int spacedim>
1915{
1916 return fe_collection;
1917}
1918
1919
1920
1921template <int dim, int spacedim>
1925{
1926 Assert(tria != nullptr,
1927 ExcMessage("This DoFHandler object has not been associated "
1928 "with a triangulation."));
1929 return *tria;
1930}
1931
1932
1933
1934template <int dim, int spacedim>
1937{
1938 Assert(tria != nullptr,
1939 ExcMessage("This DoFHandler object has not been associated "
1940 "with a triangulation."));
1941 return tria->get_communicator();
1942}
1943
1944
1945
1946template <int dim, int spacedim>
1949{
1950 Assert(this->hp_capability_enabled == false, ExcNotImplementedWithHP());
1951
1952 return block_info_object;
1953}
1954
1955
1956
1957template <int dim, int spacedim>
1959template <typename number>
1961 const std::map<types::boundary_id, const Function<spacedim, number> *>
1962 &boundary_ids) const
1963{
1964 Assert(!(dim == 2 && spacedim == 3) || this->hp_capability_enabled == false,
1965 ExcNotImplementedWithHP());
1966
1967 // extract the set of boundary ids and forget about the function object
1968 // pointers
1969 std::set<types::boundary_id> boundary_ids_only;
1970 for (typename std::map<types::boundary_id,
1971 const Function<spacedim, number> *>::const_iterator p =
1972 boundary_ids.begin();
1973 p != boundary_ids.end();
1974 ++p)
1975 boundary_ids_only.insert(p->first);
1976
1977 // then just hand everything over to the other function that does the work
1978 return n_boundary_dofs(boundary_ids_only);
1979}
1980
1981
1982
1983namespace internal
1984{
1992 template <int dim, int spacedim>
1993 std::string
1994 policy_to_string(const ::internal::DoFHandlerImplementation::Policy::
1995 PolicyBase<dim, spacedim> &policy);
1996} // namespace internal
1997
1998
1999
2000template <int dim, int spacedim>
2002template <class Archive>
2003void DoFHandler<dim, spacedim>::save(Archive &ar, const unsigned int) const
2004{
2005 if (this->hp_capability_enabled)
2006 {
2007 ar &this->object_dof_indices;
2008 ar &this->object_dof_ptr;
2009
2010 ar &this->hp_cell_active_fe_indices;
2011 ar &this->hp_cell_future_fe_indices;
2012
2013 ar &hp_object_fe_ptr;
2014 ar &hp_object_fe_indices;
2015
2016 ar &number_cache;
2017
2018 ar &mg_number_cache;
2019
2020 // write out the number of triangulation cells and later check during
2021 // loading that this number is indeed correct; same with something that
2022 // identifies the policy
2023 const unsigned int n_cells = this->tria->n_cells();
2024 std::string policy_name =
2025 ::internal::policy_to_string(*this->policy);
2026
2027 ar &n_cells &policy_name;
2028 }
2029 else
2030 {
2031 ar &this->block_info_object;
2032 ar &number_cache;
2033
2034 ar &this->object_dof_indices;
2035 ar &this->object_dof_ptr;
2036
2037 // write out the number of triangulation cells and later check during
2038 // loading that this number is indeed correct; same with something that
2039 // identifies the FE and the policy
2040 unsigned int n_cells = this->tria->n_cells();
2041 std::string fe_name = this->get_fe(0).get_name();
2042 std::string policy_name = internal::policy_to_string(*this->policy);
2043
2044 ar &n_cells &fe_name &policy_name;
2045 }
2046}
2047
2048
2049
2050template <int dim, int spacedim>
2052template <class Archive>
2053void DoFHandler<dim, spacedim>::load(Archive &ar, const unsigned int)
2054{
2055 if (this->hp_capability_enabled)
2056 {
2057 ar &this->object_dof_indices;
2058 ar &this->object_dof_ptr;
2059
2060 ar &this->hp_cell_active_fe_indices;
2061 ar &this->hp_cell_future_fe_indices;
2062
2063 ar &hp_object_fe_ptr;
2064 ar &hp_object_fe_indices;
2065
2066 ar &number_cache;
2067
2068 ar &mg_number_cache;
2069
2070 // these are the checks that correspond to the last block in the save()
2071 // function
2072 unsigned int n_cells;
2073 std::string policy_name;
2074
2075 ar &n_cells &policy_name;
2076
2078 n_cells == this->tria->n_cells(),
2079 ExcMessage(
2080 "The object being loaded into does not match the triangulation "
2081 "that has been stored previously."));
2083 policy_name == ::internal::policy_to_string(*this->policy),
2084 ExcMessage("The policy currently associated with this DoFHandler (" +
2085 ::internal::policy_to_string(*this->policy) +
2086 ") does not match the one that was associated with the "
2087 "DoFHandler previously stored (" +
2088 policy_name + ")."));
2089 }
2090 else
2091 {
2092 ar &this->block_info_object;
2093 ar &number_cache;
2094
2095 object_dof_indices.clear();
2096
2097 object_dof_ptr.clear();
2098
2099 ar &this->object_dof_indices;
2100 ar &this->object_dof_ptr;
2101
2102 // these are the checks that correspond to the last block in the save()
2103 // function
2104 unsigned int n_cells;
2105 std::string fe_name;
2106 std::string policy_name;
2107
2108 ar &n_cells &fe_name &policy_name;
2109
2111 n_cells == this->tria->n_cells(),
2112 ExcMessage(
2113 "The object being loaded into does not match the triangulation "
2114 "that has been stored previously."));
2116 fe_name == this->get_fe(0).get_name(),
2117 ExcMessage(
2118 "The finite element associated with this DoFHandler does not match "
2119 "the one that was associated with the DoFHandler previously stored."));
2120 AssertThrow(policy_name == internal::policy_to_string(*this->policy),
2121 ExcMessage(
2122 "The policy currently associated with this DoFHandler (" +
2123 internal::policy_to_string(*this->policy) +
2124 ") does not match the one that was associated with the "
2125 "DoFHandler previously stored (" +
2126 policy_name + ")."));
2127 }
2128}
2129
2130
2131
2132template <int dim, int spacedim>
2136 const unsigned int level,
2137 const unsigned int dof_number,
2138 const unsigned int dofs_per_vertex)
2139{
2142 return indices[dofs_per_vertex * (level - coarsest_level) + dof_number];
2143}
2144
2145
2146
2147extern template class DoFHandler<1, 1>;
2148extern template class DoFHandler<1, 2>;
2149extern template class DoFHandler<1, 3>;
2150extern template class DoFHandler<2, 2>;
2151extern template class DoFHandler<2, 3>;
2152extern template class DoFHandler<3, 3>;
2153
2154
2155#endif // DOXYGEN
2156
2158
2159#endif
A small class collecting the different BlockIndices involved in global, multilevel and local computat...
Definition block_info.h:95
types::global_dof_index & access_index(const unsigned int level, const unsigned int dof_number, const unsigned int dofs_per_vertex)
void init(const unsigned int coarsest_level, const unsigned int finest_level, const unsigned int dofs_per_vertex)
unsigned int get_finest_level() const
std::unique_ptr< types::global_dof_index[]> indices
unsigned int get_coarsest_level() const
cell_iterator end() const
void pre_distributed_transfer_action()
void post_transfer_action()
virtual std::size_t memory_consumption() const
std::vector< types::fe_index > get_future_fe_indices() const
unsigned int max_couplings_between_dofs() const
std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel< dim > > > mg_levels
std::unique_ptr< ActiveFEIndexTransfer > active_fe_index_transfer
hp::FECollection< dim, spacedim > fe_collection
void create_active_fe_table()
std::vector< std::vector< types::fe_index > > hp_cell_future_fe_indices
SmartPointer< const Triangulation< dim, spacedim >, DoFHandler< dim, spacedim > > tria
const hp::FECollection< dim, spacedim > & get_fe_collection() const
void save(Archive &ar, const unsigned int version) const
void set_active_fe_indices(const std::vector< unsigned int > &active_fe_indices)
BlockInfo block_info_object
bool has_active_dofs() const
std::vector< boost::signals2::connection > tria_listeners_for_transfer
level_cell_iterator end_mg() const
void renumber_dofs(const std::vector< types::global_dof_index > &new_numbers)
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
void renumber_dofs(const unsigned int level, const std::vector< types::global_dof_index > &new_numbers)
level_cell_iterator begin_mg(const unsigned int level=0) const
void distribute_dofs(const FiniteElement< dim, spacedim > &fe)
types::global_dof_index n_boundary_dofs() const
DoFHandler(const Triangulation< dim, spacedim > &tria)
void clear_mg_space()
void connect_to_triangulation_signals()
std::vector< std::vector< types::fe_index > > hp_cell_active_fe_indices
std::vector< MGVertexDoFs > mg_vertex_dofs
std::vector< boost::signals2::connection > tria_listeners
void post_distributed_transfer_action()
std::vector< std::array< std::vector< types::global_dof_index >, dim+1 > > object_dof_indices
const Triangulation< dim, spacedim > & get_triangulation() const
typename LevelSelector::CellAccessor level_cell_accessor
void pre_transfer_action()
void clear_space()
const IndexSet & locally_owned_dofs() const
void reinit(const Triangulation< dim, spacedim > &tria)
std::vector< std::array< std::vector< offset_type >, dim+1 > > object_dof_ptr
std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces< dim > > mg_faces
void serialize(Archive &archive, const unsigned int version)
types::fe_index active_fe_index_type
DoFHandler & operator=(const DoFHandler &)=delete
DoFHandler(const DoFHandler &)=delete
active_cell_iterator begin_active(const unsigned int level=0) const
std::vector< types::fe_index > get_active_fe_indices() const
void distribute_mg_dofs()
active_cell_iterator end_active(const unsigned int level) const
bool hp_capability_enabled
bool has_hp_capabilities() const
void initialize_local_block_info()
types::global_dof_index n_dofs() const
void update_active_fe_table()
level_cell_iterator end_mg(const unsigned int level) const
cell_iterator end(const unsigned int level) const
typename LevelSelector::cell_iterator level_cell_iterator
void prepare_for_serialization_of_active_fe_indices()
std::vector<::internal::DoFHandlerImplementation::NumberCache > mg_number_cache
void load(Archive &ar, const unsigned int version)
typename LevelSelector::face_iterator level_face_iterator
cell_iterator begin(const unsigned int level=0) const
std::array< std::vector< offset_type >, dim+1 > hp_object_fe_ptr
MPI_Comm get_communicator() const
void clear()
void get_active_fe_indices(std::vector< unsigned int > &active_fe_indices) const
void distribute_dofs(const hp::FECollection< dim, spacedim > &fe)
types::global_dof_index n_boundary_dofs(const std::map< types::boundary_id, const Function< spacedim, number > * > &boundary_ids) const
bool has_level_dofs() const
std::unique_ptr<::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > > policy
typename LevelSelector::FaceAccessor level_face_accessor
types::global_dof_index n_dofs(const unsigned int level) const
std::array< std::vector< types::fe_index >, dim+1 > hp_object_fe_indices
void set_active_fe_indices(const std::vector< types::fe_index > &active_fe_indices)
unsigned int max_couplings_between_boundary_dofs() const
void setup_policy()
types::global_dof_index n_locally_owned_dofs() const
void deserialize_active_fe_indices()
types::global_dof_index n_boundary_dofs(const std::set< types::boundary_id > &boundary_ids) const
const BlockInfo & block_info() const
::internal::DoFHandlerImplementation::NumberCache number_cache
void set_future_fe_indices(const std::vector< types::fe_index > &future_fe_indices)
virtual ~DoFHandler() override
#define DEAL_II_DEPRECATED
Definition config.h:207
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:177
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
unsigned int level
Definition grid_out.cc:4626
IteratorRange< active_cell_iterator > active_cell_iterators() const
IteratorRange< level_cell_iterator > mg_cell_iterators() const
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
IteratorRange< cell_iterator > cell_iterators() const
IteratorRange< active_cell_iterator > active_cell_iterators_on_level(const unsigned int level) const
IteratorRange< level_cell_iterator > mg_cell_iterators_on_level(const unsigned int level) const
#define DeclException0(Exception0)
Definition exceptions.h:471
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition exceptions.h:539
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:494
static ::ExceptionBase & ExcNoDominatedFiniteElementOnChildren()
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:516
static ::ExceptionBase & ExcInvalidLevel(int arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
typename ActiveSelector::FaceAccessor face_accessor
typename ActiveSelector::cell_iterator cell_iterator
typename ActiveSelector::line_iterator line_iterator
typename ActiveSelector::face_iterator face_iterator
typename ActiveSelector::active_line_iterator active_line_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
typename ActiveSelector::CellAccessor cell_accessor
typename ActiveSelector::quad_iterator quad_iterator
typename ActiveSelector::active_hex_iterator active_hex_iterator
typename ActiveSelector::active_quad_iterator active_quad_iterator
typename ActiveSelector::hex_iterator hex_iterator
typename ActiveSelector::active_face_iterator active_face_iterator
Definition hp.h:117
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
Definition tria.cc:14883
unsigned int dominated_future_fe_on_children(const typename DoFHandler< dim, spacedim >::cell_iterator &parent)
void communicate_future_fe_indices(DoFHandler< dim, spacedim > &dof_handler)
std::string policy_to_string(const ::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > &policy)
const types::fe_index invalid_fe_index
Definition types.h:243
unsigned short int fe_index
Definition types.h:59
std::map< const cell_iterator, const types::fe_index > refined_cells_fe_index
std::map< const cell_iterator, const types::fe_index > persisting_cells_fe_index
std::vector< types::fe_index > active_fe_indices
std::map< const cell_iterator, const types::fe_index > coarsened_cells_fe_index
std::unique_ptr< parallel::distributed::CellDataTransfer< dim, spacedim, std::vector< types::fe_index > > > cell_data_transfer