Reference documentation for deal.II version Git 39b7d3efb0 2021-05-07 15:49:09 -0400
\(\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\}}\)
dof_handler.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2021 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_dof_handler_h
17 #define dealii_dof_handler_h
18 
19 
20 
21 #include <deal.II/base/config.h>
22 
24 #include <deal.II/base/function.h>
25 #include <deal.II/base/index_set.h>
28 
30 
33 #include <deal.II/dofs/dof_faces.h>
37 
39 
40 #include <boost/serialization/split_member.hpp>
41 
42 #include <map>
43 #include <memory>
44 #include <set>
45 #include <vector>
46 
48 
49 // Forward declarations
50 #ifndef DOXYGEN
51 template <int dim, int spacedim>
52 class FiniteElement;
53 template <int dim, int spacedim>
54 class Triangulation;
55 
56 namespace internal
57 {
58  namespace DoFHandlerImplementation
59  {
60  struct Implementation;
61 
62  namespace Policy
63  {
64  template <int dim, int spacedim>
65  class PolicyBase;
66  struct Implementation;
67  } // namespace Policy
68  } // namespace DoFHandlerImplementation
69 
70  namespace DoFAccessorImplementation
71  {
72  struct Implementation;
73  }
74 
75  namespace DoFCellAccessorImplementation
76  {
77  struct Implementation;
78  }
79 
80  namespace hp
81  {
82  namespace DoFHandlerImplementation
83  {
84  struct Implementation;
85  }
86  } // namespace hp
87 } // namespace internal
88 
89 namespace parallel
90 {
91  namespace distributed
92  {
93  template <int dim, int spacedim, typename VectorType>
94  class CellDataTransfer;
95  }
96 } // namespace parallel
97 #endif
98 
313 template <int dim, int spacedim = dim>
314 class DoFHandler : public Subscriptor
315 {
316  using ActiveSelector =
318  using LevelSelector =
320 
321 public:
334  using cell_accessor = typename ActiveSelector::CellAccessor;
335 
348  using face_accessor = typename ActiveSelector::FaceAccessor;
349 
357  using line_iterator = typename ActiveSelector::line_iterator;
358 
372  using active_line_iterator = typename ActiveSelector::active_line_iterator;
373 
381  using quad_iterator = typename ActiveSelector::quad_iterator;
382 
396  using active_quad_iterator = typename ActiveSelector::active_quad_iterator;
397 
405  using hex_iterator = typename ActiveSelector::hex_iterator;
406 
416  using active_hex_iterator = typename ActiveSelector::active_hex_iterator;
417 
438  using active_cell_iterator = typename ActiveSelector::active_cell_iterator;
439 
466  using cell_iterator = typename ActiveSelector::cell_iterator;
467 
484  using face_iterator = typename ActiveSelector::face_iterator;
485 
497  using active_face_iterator = typename ActiveSelector::active_face_iterator;
498 
499  using level_cell_accessor = typename LevelSelector::CellAccessor;
500  using level_face_accessor = typename LevelSelector::FaceAccessor;
501 
502  using level_cell_iterator = typename LevelSelector::cell_iterator;
503  using level_face_iterator = typename LevelSelector::face_iterator;
504 
505 
509  static const unsigned int dimension = dim;
510 
514  static const unsigned int space_dimension = spacedim;
515 
519  static const unsigned int default_fe_index = 0;
520 
524  static const unsigned int invalid_fe_index = numbers::invalid_unsigned_int;
525 
529  using active_fe_index_type = unsigned short int;
530 
534  using offset_type = unsigned int;
535 
540  static const active_fe_index_type invalid_active_fe_index =
541  static_cast<active_fe_index_type>(-1);
542 
547  DoFHandler();
548 
552  explicit DoFHandler(const Triangulation<dim, spacedim> &tria);
553 
560  DoFHandler(const DoFHandler &) = delete;
561 
565  virtual ~DoFHandler() override;
566 
573  DoFHandler &
574  operator=(const DoFHandler &) = delete;
575 
583  void
584  initialize(const Triangulation<dim, spacedim> &tria,
585  const FiniteElement<dim, spacedim> &fe);
586 
593  void
594  initialize(const Triangulation<dim, spacedim> & tria,
596 
621  void
622  set_fe(const FiniteElement<dim, spacedim> &fe);
623 
630  void
631  set_fe(const hp::FECollection<dim, spacedim> &fe);
632 
637  void
638  set_active_fe_indices(const std::vector<unsigned int> &active_fe_indices);
639 
645  void
646  get_active_fe_indices(std::vector<unsigned int> &active_fe_indices) const;
647 
655  void
657 
687  void
688  distribute_dofs(const FiniteElement<dim, spacedim> &fe);
689 
693  void
694  distribute_dofs(const hp::FECollection<dim, spacedim> &fe);
695 
701  void
702  distribute_mg_dofs();
703 
707  bool
708  has_hp_capabilities() const;
709 
715  bool
716  has_level_dofs() const;
717 
732  bool
733  has_active_dofs() const;
734 
741  void
742  initialize_local_block_info();
743 
747  void
748  clear();
749 
802  void
803  renumber_dofs(const std::vector<types::global_dof_index> &new_numbers);
804 
809  void
810  renumber_dofs(const unsigned int level,
811  const std::vector<types::global_dof_index> &new_numbers);
812 
841  unsigned int
842  max_couplings_between_dofs() const;
843 
856  unsigned int
857  max_couplings_between_boundary_dofs() const;
858 
859  /*--------------------------------------*/
860 
865  /*
866  * @{
867  */
868 
873  begin(const unsigned int level = 0) const;
874 
892  begin_active(const unsigned int level = 0) const;
893 
899  end() const;
900 
906  end(const unsigned int level) const;
907 
914  end_active(const unsigned int level) const;
915 
921  begin_mg(const unsigned int level = 0) const;
922 
928  end_mg(const unsigned int level) const;
929 
935  end_mg() const;
936 
952  cell_iterators() const;
953 
995  active_cell_iterators() const;
996 
1009  mg_cell_iterators() const;
1010 
1027  cell_iterators_on_level(const unsigned int level) const;
1028 
1045  active_cell_iterators_on_level(const unsigned int level) const;
1046 
1063  mg_cell_iterators_on_level(const unsigned int level) const;
1064 
1065  /*
1066  * @}
1067  */
1068 
1069 
1070  /*---------------------------------------*/
1071 
1072 
1097  n_dofs() const;
1098 
1107  n_dofs(const unsigned int level) const;
1108 
1114  n_boundary_dofs() const;
1115 
1127  template <typename number>
1129  n_boundary_dofs(
1130  const std::map<types::boundary_id, const Function<spacedim, number> *>
1131  &boundary_ids) const;
1132 
1138  n_boundary_dofs(const std::set<types::boundary_id> &boundary_ids) const;
1139 
1155  const BlockInfo &
1156  block_info() const;
1157 
1178  n_locally_owned_dofs() const;
1179 
1185  const IndexSet &
1186  locally_owned_dofs() const;
1187 
1192  const IndexSet &
1193  locally_owned_mg_dofs(const unsigned int level) const;
1194 
1208  DEAL_II_DEPRECATED const std::vector<IndexSet> &
1209  locally_owned_dofs_per_processor() const;
1210 
1227  DEAL_II_DEPRECATED const std::vector<types::global_dof_index> &
1228  n_locally_owned_dofs_per_processor() const;
1229 
1244  DEAL_II_DEPRECATED const std::vector<IndexSet> &
1245  locally_owned_mg_dofs_per_processor(const unsigned int level) const;
1246 
1252  get_fe(const unsigned int index = 0) const;
1253 
1259  get_fe_collection() const;
1260 
1265  get_triangulation() const;
1266 
1270  MPI_Comm
1271  get_communicator() const;
1272 
1289  void
1290  prepare_for_serialization_of_active_fe_indices();
1291 
1307  void
1308  deserialize_active_fe_indices();
1309 
1318  virtual std::size_t
1319  memory_consumption() const;
1320 
1326  template <class Archive>
1327  void
1328  save(Archive &ar, const unsigned int version) const;
1329 
1335  template <class Archive>
1336  void
1337  load(Archive &ar, const unsigned int version);
1338 
1339 #ifdef DOXYGEN
1340 
1345  template <class Archive>
1346  void
1347  serialize(Archive &archive, const unsigned int version);
1348 #else
1349  // This macro defines the serialize() method that is compatible with
1350  // the templated save() and load() method that have been implemented.
1351  BOOST_SERIALIZATION_SPLIT_MEMBER()
1352 #endif
1353 
1367  DeclException1(ExcInvalidLevel,
1368  int,
1369  << "The given level " << arg1
1370  << " is not in the valid range!");
1375  DeclException1(ExcNewNumbersNotConsecutive,
1377  << "The given list of new dof indices is not consecutive: "
1378  << "the index " << arg1 << " does not exist.");
1382  DeclException2(ExcInvalidFEIndex,
1383  int,
1384  int,
1385  << "The mesh contains a cell with an active FE index of "
1386  << arg1 << ", but the finite element collection only has "
1387  << arg2 << " elements");
1388 
1393  DeclExceptionMsg(ExcOnlyAvailableWithHP,
1394  "The current function doesn't make sense when used with a "
1395  "DoFHandler without hp-capabilities.");
1396 
1401  DeclExceptionMsg(ExcNotImplementedWithHP,
1402  "The current function has not yet been implemented for a "
1403  "DoFHandler with hp-capabilities.");
1404 
1405 private:
1413  {
1414  public:
1418  MGVertexDoFs();
1419 
1425  void
1426  init(const unsigned int coarsest_level,
1427  const unsigned int finest_level,
1428  const unsigned int dofs_per_vertex);
1429 
1433  unsigned int
1434  get_coarsest_level() const;
1435 
1439  unsigned int
1440  get_finest_level() const;
1441 
1447  get_index(const unsigned int level,
1448  const unsigned int dof_number,
1449  const unsigned int dofs_per_vertex) const;
1450 
1455  void
1456  set_index(const unsigned int level,
1457  const unsigned int dof_number,
1458  const unsigned int dofs_per_vertex,
1459  const types::global_dof_index index);
1460 
1461  private:
1465  unsigned int coarsest_level;
1466 
1470  unsigned int finest_level;
1471 
1481  std::unique_ptr<types::global_dof_index[]> indices;
1482  };
1483 
1491  {
1496  std::map<const cell_iterator, const unsigned int> persisting_cells_fe_index;
1497 
1502  std::map<const cell_iterator, const unsigned int> refined_cells_fe_index;
1503 
1508  std::map<const cell_iterator, const unsigned int> coarsened_cells_fe_index;
1509 
1515  std::vector<unsigned int> active_fe_indices;
1516 
1522  std::unique_ptr<
1523  parallel::distributed::
1524  CellDataTransfer<dim, spacedim, std::vector<unsigned int>>>
1526  };
1527 
1532 
1538 
1544 
1551 
1556  std::unique_ptr<::internal::DoFHandlerImplementation::Policy::
1557  PolicyBase<dim, spacedim>>
1559 
1568 
1572  std::vector<::internal::DoFHandlerImplementation::NumberCache>
1574 
1580  mutable std::vector<std::vector<types::global_dof_index>>
1582 
1587  mutable std::vector<std::vector<offset_type>> cell_dof_cache_ptr;
1588 
1594  mutable std::vector<std::array<std::vector<types::global_dof_index>, dim + 1>>
1596 
1605  mutable std::vector<std::array<std::vector<offset_type>, dim + 1>>
1607 
1613  mutable std::array<std::vector<active_fe_index_type>, dim + 1>
1615 
1619  mutable std::array<std::vector<offset_type>, dim + 1> hp_object_fe_ptr;
1620 
1625  mutable std::vector<std::vector<active_fe_index_type>>
1627 
1632  mutable std::vector<std::vector<active_fe_index_type>>
1634 
1639  std::vector<MGVertexDoFs> mg_vertex_dofs;
1640 
1644  std::vector<
1645  std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel<dim>>>
1647 
1651  std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces<dim>>
1653 
1658  std::unique_ptr<ActiveFEIndexTransfer> active_fe_index_transfer;
1659 
1664  std::vector<boost::signals2::connection> tria_listeners;
1665 
1671  std::vector<boost::signals2::connection> tria_listeners_for_transfer;
1672 
1676  void
1677  clear_space();
1678 
1682  void
1683  clear_mg_space();
1684 
1688  template <int structdim>
1690  get_dof_index(const unsigned int obj_level,
1691  const unsigned int obj_index,
1692  const unsigned int fe_index,
1693  const unsigned int local_index) const;
1694 
1698  template <int structdim>
1699  void
1700  set_dof_index(const unsigned int obj_level,
1701  const unsigned int obj_index,
1702  const unsigned int fe_index,
1703  const unsigned int local_index,
1705 
1709  void
1710  setup_policy();
1711 
1715  void
1716  connect_to_triangulation_signals();
1717 
1729  void
1730  create_active_fe_table();
1731 
1744  void
1745  update_active_fe_table();
1746 
1756  void
1757  pre_transfer_action();
1758 
1767  void
1768  post_transfer_action();
1769 
1778  void
1779  pre_distributed_transfer_action();
1780 
1789  void
1790  post_distributed_transfer_action();
1791 
1792 
1793  // Make accessor objects friends.
1794  template <int, int, int, bool>
1795  friend class ::DoFAccessor;
1796  template <int, int, bool>
1797  friend class ::DoFCellAccessor;
1798  friend struct ::internal::DoFAccessorImplementation::Implementation;
1799  friend struct ::internal::DoFCellAccessorImplementation::Implementation;
1800 
1801  // Likewise for DoFLevel objects since they need to access the vertex dofs
1802  // in the functions that set and retrieve vertex dof indices.
1803  friend struct ::internal::DoFHandlerImplementation::Implementation;
1804  friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
1805  friend struct ::internal::DoFHandlerImplementation::Policy::
1806  Implementation;
1807 
1808  // explicitly check for sensible template arguments, but not on windows
1809  // because MSVC creates bogus warnings during normal compilation
1810 #ifndef DEAL_II_MSVC
1811  static_assert(dim <= spacedim,
1812  "The dimension <dim> of a DoFHandler must be less than or "
1813  "equal to the space dimension <spacedim> in which it lives.");
1814 #endif
1815 };
1816 
1817 namespace internal
1818 {
1819  namespace hp
1820  {
1821  namespace DoFHandlerImplementation
1822  {
1834  template <int dim, int spacedim>
1835  void
1837 
1862  template <int dim, int spacedim = dim>
1863  unsigned int
1865  const typename DoFHandler<dim, spacedim>::cell_iterator &parent);
1866 
1872  "No FiniteElement has been found in your FECollection that is "
1873  "dominated by all children of a cell you are trying to coarsen!");
1874  } // namespace DoFHandlerImplementation
1875  } // namespace hp
1876 } // namespace internal
1877 
1878 #ifndef DOXYGEN
1879 
1880 /* ----------------------- Inline functions ----------------------------------
1881  */
1882 
1883 
1884 template <int dim, int spacedim>
1885 inline bool
1887 {
1888  return hp_capability_enabled;
1889 }
1890 
1891 
1892 
1893 template <int dim, int spacedim>
1894 inline bool
1896 {
1897  return mg_number_cache.size() > 0;
1898 }
1899 
1900 
1901 
1902 template <int dim, int spacedim>
1903 inline bool
1905 {
1906  return number_cache.n_global_dofs > 0;
1907 }
1908 
1909 
1910 
1911 template <int dim, int spacedim>
1914 {
1915  return number_cache.n_global_dofs;
1916 }
1917 
1918 
1919 
1920 template <int dim, int spacedim>
1922 DoFHandler<dim, spacedim>::n_dofs(const unsigned int level) const
1923 {
1924  Assert(has_level_dofs(),
1925  ExcMessage(
1926  "n_dofs(level) can only be called after distribute_mg_dofs()"));
1927  Assert(level < mg_number_cache.size(), ExcInvalidLevel(level));
1928  return mg_number_cache[level].n_global_dofs;
1929 }
1930 
1931 
1932 
1933 template <int dim, int spacedim>
1936 {
1937  return number_cache.n_locally_owned_dofs;
1938 }
1939 
1940 
1941 
1942 template <int dim, int spacedim>
1943 const IndexSet &
1945 {
1946  return number_cache.locally_owned_dofs;
1947 }
1948 
1949 
1950 
1951 template <int dim, int spacedim>
1952 const IndexSet &
1953 DoFHandler<dim, spacedim>::locally_owned_mg_dofs(const unsigned int level) const
1954 {
1955  Assert(level < this->get_triangulation().n_global_levels(),
1956  ExcMessage("The given level index exceeds the number of levels "
1957  "present in the triangulation"));
1958  Assert(
1959  mg_number_cache.size() == this->get_triangulation().n_global_levels(),
1960  ExcMessage(
1961  "The level dofs are not set up properly! Did you call distribute_mg_dofs()?"));
1962  return mg_number_cache[level].locally_owned_dofs;
1963 }
1964 
1965 
1966 
1967 template <int dim, int spacedim>
1968 const std::vector<types::global_dof_index> &
1970 {
1971  if (number_cache.n_locally_owned_dofs_per_processor.empty() &&
1972  number_cache.n_global_dofs > 0)
1973  {
1975  number_cache)
1976  .n_locally_owned_dofs_per_processor =
1977  number_cache.get_n_locally_owned_dofs_per_processor(get_communicator());
1978  }
1979  return number_cache.n_locally_owned_dofs_per_processor;
1980 }
1981 
1982 
1983 
1984 template <int dim, int spacedim>
1985 const std::vector<IndexSet> &
1987 {
1988  if (number_cache.locally_owned_dofs_per_processor.empty() &&
1989  number_cache.n_global_dofs > 0)
1990  {
1992  number_cache)
1993  .locally_owned_dofs_per_processor =
1994  number_cache.get_locally_owned_dofs_per_processor(get_communicator());
1995  }
1996  return number_cache.locally_owned_dofs_per_processor;
1997 }
1998 
1999 
2000 
2001 template <int dim, int spacedim>
2002 const std::vector<IndexSet> &
2004  const unsigned int level) const
2005 {
2006  Assert(level < this->get_triangulation().n_global_levels(),
2007  ExcMessage("The given level index exceeds the number of levels "
2008  "present in the triangulation"));
2009  Assert(
2010  mg_number_cache.size() == this->get_triangulation().n_global_levels(),
2011  ExcMessage(
2012  "The level dofs are not set up properly! Did you call distribute_mg_dofs()?"));
2013  if (mg_number_cache[level].locally_owned_dofs_per_processor.empty() &&
2014  mg_number_cache[level].n_global_dofs > 0)
2015  {
2017  mg_number_cache[level])
2018  .locally_owned_dofs_per_processor =
2019  mg_number_cache[level].get_locally_owned_dofs_per_processor(
2020  get_communicator());
2021  }
2022  return mg_number_cache[level].locally_owned_dofs_per_processor;
2023 }
2024 
2025 
2026 
2027 template <int dim, int spacedim>
2028 inline const FiniteElement<dim, spacedim> &
2029 DoFHandler<dim, spacedim>::get_fe(const unsigned int number) const
2030 {
2031  Assert(fe_collection.size() > 0,
2032  ExcMessage("No finite element collection is associated with "
2033  "this DoFHandler"));
2034  return fe_collection[number];
2035 }
2036 
2037 
2038 
2039 template <int dim, int spacedim>
2040 inline const hp::FECollection<dim, spacedim> &
2042 {
2043  return fe_collection;
2044 }
2045 
2046 
2047 
2048 template <int dim, int spacedim>
2049 inline const Triangulation<dim, spacedim> &
2051 {
2052  Assert(tria != nullptr,
2053  ExcMessage("This DoFHandler object has not been associated "
2054  "with a triangulation."));
2055  return *tria;
2056 }
2057 
2058 
2059 
2060 template <int dim, int spacedim>
2061 inline MPI_Comm
2063 {
2064  Assert(tria != nullptr,
2065  ExcMessage("This DoFHandler object has not been associated "
2066  "with a triangulation."));
2067  return tria->get_communicator();
2068 }
2069 
2070 
2071 
2072 template <int dim, int spacedim>
2073 inline const BlockInfo &
2075 {
2076  Assert(this->hp_capability_enabled == false, ExcNotImplementedWithHP());
2077 
2078  return block_info_object;
2079 }
2080 
2081 
2082 
2083 template <int dim, int spacedim>
2084 template <typename number>
2087  const std::map<types::boundary_id, const Function<spacedim, number> *>
2088  &boundary_ids) const
2089 {
2090  Assert(!(dim == 2 && spacedim == 3) || this->hp_capability_enabled == false,
2091  ExcNotImplementedWithHP());
2092 
2093  // extract the set of boundary ids and forget about the function object
2094  // pointers
2095  std::set<types::boundary_id> boundary_ids_only;
2096  for (typename std::map<types::boundary_id,
2097  const Function<spacedim, number> *>::const_iterator p =
2098  boundary_ids.begin();
2099  p != boundary_ids.end();
2100  ++p)
2101  boundary_ids_only.insert(p->first);
2102 
2103  // then just hand everything over to the other function that does the work
2104  return n_boundary_dofs(boundary_ids_only);
2105 }
2106 
2107 
2108 
2109 namespace internal
2110 {
2118  template <int dim, int spacedim>
2119  std::string
2120  policy_to_string(const ::internal::DoFHandlerImplementation::Policy::
2121  PolicyBase<dim, spacedim> &policy);
2122 } // namespace internal
2123 
2124 
2125 
2126 template <int dim, int spacedim>
2127 template <class Archive>
2128 void
2129 DoFHandler<dim, spacedim>::save(Archive &ar, const unsigned int) const
2130 {
2131  if (this->hp_capability_enabled)
2132  {
2133  ar & this->object_dof_indices;
2134  ar & this->object_dof_ptr;
2135 
2136  ar & this->cell_dof_cache_indices;
2137  ar & this->cell_dof_cache_ptr;
2138 
2139  ar & this->hp_cell_active_fe_indices;
2140  ar & this->hp_cell_future_fe_indices;
2141 
2142  ar &hp_object_fe_ptr;
2143  ar &hp_object_fe_indices;
2144 
2145  ar &number_cache;
2146 
2147  ar &mg_number_cache;
2148 
2149  // write out the number of triangulation cells and later check during
2150  // loading that this number is indeed correct; same with something that
2151  // identifies the policy
2152  const unsigned int n_cells = this->tria->n_cells();
2153  std::string policy_name =
2154  ::internal::policy_to_string(*this->policy);
2155 
2156  ar &n_cells &policy_name;
2157  }
2158  else
2159  {
2160  ar & this->block_info_object;
2161  ar &number_cache;
2162 
2163  ar & this->object_dof_indices;
2164  ar & this->object_dof_ptr;
2165 
2166  ar & this->cell_dof_cache_indices;
2167  ar & this->cell_dof_cache_ptr;
2168 
2169  // write out the number of triangulation cells and later check during
2170  // loading that this number is indeed correct; same with something that
2171  // identifies the FE and the policy
2172  unsigned int n_cells = this->tria->n_cells();
2173  std::string fe_name = this->get_fe(0).get_name();
2174  std::string policy_name = internal::policy_to_string(*this->policy);
2175 
2176  ar &n_cells &fe_name &policy_name;
2177  }
2178 }
2179 
2180 
2181 
2182 template <int dim, int spacedim>
2183 template <class Archive>
2184 void
2185 DoFHandler<dim, spacedim>::load(Archive &ar, const unsigned int)
2186 {
2187  if (this->hp_capability_enabled)
2188  {
2189  ar & this->object_dof_indices;
2190  ar & this->object_dof_ptr;
2191 
2192  ar & this->cell_dof_cache_indices;
2193  ar & this->cell_dof_cache_ptr;
2194 
2195  ar & this->hp_cell_active_fe_indices;
2196  ar & this->hp_cell_future_fe_indices;
2197 
2198  ar &hp_object_fe_ptr;
2199  ar &hp_object_fe_indices;
2200 
2201  ar &number_cache;
2202 
2203  ar &mg_number_cache;
2204 
2205  // these are the checks that correspond to the last block in the save()
2206  // function
2207  unsigned int n_cells;
2208  std::string policy_name;
2209 
2210  ar &n_cells &policy_name;
2211 
2212  AssertThrow(
2213  n_cells == this->tria->n_cells(),
2214  ExcMessage(
2215  "The object being loaded into does not match the triangulation "
2216  "that has been stored previously."));
2217  AssertThrow(
2218  policy_name == ::internal::policy_to_string(*this->policy),
2219  ExcMessage("The policy currently associated with this DoFHandler (" +
2220  ::internal::policy_to_string(*this->policy) +
2221  ") does not match the one that was associated with the "
2222  "DoFHandler previously stored (" +
2223  policy_name + ")."));
2224  }
2225  else
2226  {
2227  ar & this->block_info_object;
2228  ar &number_cache;
2229 
2230  object_dof_indices.clear();
2231 
2232  object_dof_ptr.clear();
2233 
2234  ar & this->object_dof_indices;
2235  ar & this->object_dof_ptr;
2236 
2237  ar & this->cell_dof_cache_indices;
2238  ar & this->cell_dof_cache_ptr;
2239 
2240  // these are the checks that correspond to the last block in the save()
2241  // function
2242  unsigned int n_cells;
2243  std::string fe_name;
2244  std::string policy_name;
2245 
2246  ar &n_cells &fe_name &policy_name;
2247 
2248  AssertThrow(
2249  n_cells == this->tria->n_cells(),
2250  ExcMessage(
2251  "The object being loaded into does not match the triangulation "
2252  "that has been stored previously."));
2253  AssertThrow(
2254  fe_name == this->get_fe(0).get_name(),
2255  ExcMessage(
2256  "The finite element associated with this DoFHandler does not match "
2257  "the one that was associated with the DoFHandler previously stored."));
2258  AssertThrow(policy_name == internal::policy_to_string(*this->policy),
2259  ExcMessage(
2260  "The policy currently associated with this DoFHandler (" +
2261  internal::policy_to_string(*this->policy) +
2262  ") does not match the one that was associated with the "
2263  "DoFHandler previously stored (" +
2264  policy_name + ")."));
2265  }
2266 }
2267 
2268 
2269 
2270 template <int dim, int spacedim>
2273  const unsigned int level,
2274  const unsigned int dof_number,
2275  const unsigned int dofs_per_vertex) const
2276 {
2277  Assert((level >= coarsest_level) && (level <= finest_level),
2278  ExcInvalidLevel(level));
2279  return indices[dofs_per_vertex * (level - coarsest_level) + dof_number];
2280 }
2281 
2282 
2283 
2284 template <int dim, int spacedim>
2285 inline void
2287  const unsigned int level,
2288  const unsigned int dof_number,
2289  const unsigned int dofs_per_vertex,
2290  const types::global_dof_index index)
2291 {
2292  Assert((level >= coarsest_level) && (level <= finest_level),
2293  ExcInvalidLevel(level));
2294  indices[dofs_per_vertex * (level - coarsest_level) + dof_number] = index;
2295 }
2296 
2297 
2298 
2299 extern template class DoFHandler<1, 1>;
2300 extern template class DoFHandler<1, 2>;
2301 extern template class DoFHandler<1, 3>;
2302 extern template class DoFHandler<2, 2>;
2303 extern template class DoFHandler<2, 3>;
2304 extern template class DoFHandler<3, 3>;
2305 
2306 
2307 #endif // DOXYGEN
2308 
2310 
2311 #endif
std::vector< types::global_dof_index > get_n_locally_owned_dofs_per_processor(const MPI_Comm &mpi_communicator) const
Definition: number_cache.cc:87
std::vector< MGVertexDoFs > mg_vertex_dofs
Definition: dof_handler.h:1639
static const unsigned int invalid_unsigned_int
Definition: types.h:196
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:618
typename ActiveSelector::FaceAccessor face_accessor
Definition: dof_handler.h:348
std::array< std::vector< offset_type >, dim+1 > hp_object_fe_ptr
Definition: dof_handler.h:1619
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:538
typename LevelSelector::face_iterator level_face_iterator
Definition: dof_handler.h:503
std::map< const cell_iterator, const unsigned int > coarsened_cells_fe_index
Definition: dof_handler.h:1508
static ::ExceptionBase & ExcNoDominatedFiniteElementOnChildren()
typename ActiveSelector::quad_iterator quad_iterator
Definition: dof_handler.h:381
void communicate_future_fe_indices(DoFHandler< dim, spacedim > &dof_handler)
unsigned short int active_fe_index_type
Definition: dof_handler.h:529
::internal::DoFHandlerImplementation::NumberCache number_cache
Definition: dof_handler.h:1567
std::vector< std::array< std::vector< offset_type >, dim+1 > > object_dof_ptr
Definition: dof_handler.h:1606
typename ActiveSelector::CellAccessor cell_accessor
Definition: dof_handler.h:334
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel< dim > > > mg_levels
Definition: dof_handler.h:1646
void save(Archive &ar, const unsigned int version) const
std::vector< std::vector< active_fe_index_type > > hp_cell_future_fe_indices
Definition: dof_handler.h:1633
typename ActiveSelector::face_iterator face_iterator
Definition: dof_handler.h:484
typename ActiveSelector::active_face_iterator active_face_iterator
Definition: dof_handler.h:497
MPI_Comm get_communicator() const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
std::vector< std::vector< types::global_dof_index > > cell_dof_cache_indices
Definition: dof_handler.h:1581
std::map< const cell_iterator, const unsigned int > refined_cells_fe_index
Definition: dof_handler.h:1502
BlockInfo block_info_object
Definition: dof_handler.h:1531
static ::ExceptionBase & ExcMessage(std::string arg1)
std::map< const cell_iterator, const unsigned int > persisting_cells_fe_index
Definition: dof_handler.h:1496
unsigned int dominated_future_fe_on_children(const typename DoFHandler< dim, spacedim >::cell_iterator &parent)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:515
typename LevelSelector::CellAccessor level_cell_accessor
Definition: dof_handler.h:499
typename ActiveSelector::hex_iterator hex_iterator
Definition: dof_handler.h:405
std::vector<::internal::DoFHandlerImplementation::NumberCache > mg_number_cache
Definition: dof_handler.h:1573
bool has_active_dofs() const
#define Assert(cond, exc)
Definition: exceptions.h:1465
hp::FECollection< dim, spacedim > fe_collection
Definition: dof_handler.h:1550
std::vector< IndexSet > get_locally_owned_dofs_per_processor(const MPI_Comm &mpi_communicator) const
types::global_dof_index n_dofs() const
std::unique_ptr< ActiveFEIndexTransfer > active_fe_index_transfer
Definition: dof_handler.h:1658
std::vector< std::array< std::vector< types::global_dof_index >, dim+1 > > object_dof_indices
Definition: dof_handler.h:1595
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces< dim > > mg_faces
Definition: dof_handler.h:1652
#define DeclException0(Exception0)
Definition: exceptions.h:470
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:395
unsigned int level
Definition: grid_out.cc:4590
SmartPointer< const Triangulation< dim, spacedim >, DoFHandler< dim, spacedim > > tria
Definition: dof_handler.h:1543
VectorType::value_type * end(VectorType &V)
typename LevelSelector::cell_iterator level_cell_iterator
Definition: dof_handler.h:502
void load(Archive &ar, const unsigned int version)
std::unique_ptr<::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > > policy
Definition: dof_handler.h:1558
void set_index(const unsigned int level, const unsigned int dof_number, const unsigned int dofs_per_vertex, const types::global_dof_index index)
std::unique_ptr< parallel::distributed::CellDataTransfer< dim, spacedim, std::vector< unsigned int > > > cell_data_transfer
Definition: dof_handler.h:1525
Definition: hp.h:117
std::vector< std::vector< active_fe_index_type > > hp_cell_active_fe_indices
Definition: dof_handler.h:1626
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
Definition: tria.cc:12575
typename LevelSelector::FaceAccessor level_face_accessor
Definition: dof_handler.h:500
std::vector< unsigned int > active_fe_indices
Definition: dof_handler.h:1515
std::array< std::vector< active_fe_index_type >, dim+1 > hp_object_fe_indices
Definition: dof_handler.h:1614
const IndexSet & locally_owned_dofs() const
typename ActiveSelector::active_line_iterator active_line_iterator
Definition: dof_handler.h:372
std::unique_ptr< types::global_dof_index[]> indices
Definition: dof_handler.h:1481
std::string policy_to_string(const ::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > &policy)
Definition: dof_handler.cc:51
bool has_hp_capabilities() const
typename ActiveSelector::line_iterator line_iterator
Definition: dof_handler.h:357
const Triangulation< dim, spacedim > & get_triangulation() const
unsigned int coarsest_level
Definition: dof_handler.h:1465
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:394
VectorType::value_type * begin(VectorType &V)
const std::vector< IndexSet > & locally_owned_dofs_per_processor() const
types::global_dof_index get_index(const unsigned int level, const unsigned int dof_number, const unsigned int dofs_per_vertex) const
types::global_dof_index n_boundary_dofs() const
void get_active_fe_indices(const DoFHandler< dim, spacedim > &dof_handler, std::vector< unsigned int > &active_fe_indices)
Definition: dof_tools.cc:1399
std::vector< IndexSet > locally_owned_dofs_per_processor
Definition: number_cache.h:180
A small class collecting the different BlockIndices involved in global, multilevel and local computat...
Definition: block_info.h:93
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:466
const std::vector< types::global_dof_index > & n_locally_owned_dofs_per_processor() const
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
typename ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:438
std::vector< boost::signals2::connection > tria_listeners
Definition: dof_handler.h:1664
std::vector< boost::signals2::connection > tria_listeners_for_transfer
Definition: dof_handler.h:1671
bool has_level_dofs() const
#define DEAL_II_DEPRECATED
Definition: config.h:157
const BlockInfo & block_info() const
static ::ExceptionBase & ExcNoFESelected()
bool hp_capability_enabled
Definition: dof_handler.h:1537
const std::vector< IndexSet > & locally_owned_mg_dofs_per_processor(const unsigned int level) const
typename ActiveSelector::active_quad_iterator active_quad_iterator
Definition: dof_handler.h:396
types::global_dof_index n_locally_owned_dofs() const
std::vector< std::vector< offset_type > > cell_dof_cache_ptr
Definition: dof_handler.h:1587
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
typename ActiveSelector::active_hex_iterator active_hex_iterator
Definition: dof_handler.h:416