Reference documentation for deal.II version Git f40be01994 2020-04-09 07:13:12 +0200
\(\newcommand{\vcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\vcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
dof_handler.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2019 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/geometry_info.h>
17 #include <deal.II/base/memory_consumption.h>
18 #include <deal.II/base/std_cxx14/memory.h>
19 #include <deal.II/base/thread_management.h>
20 
21 #include <deal.II/distributed/cell_data_transfer.templates.h>
22 #include <deal.II/distributed/shared_tria.h>
23 #include <deal.II/distributed/tria.h>
24 
25 #include <deal.II/dofs/dof_accessor.h>
26 #include <deal.II/dofs/dof_handler_policy.h>
27 
28 #include <deal.II/fe/fe.h>
29 
30 #include <deal.II/grid/grid_tools.h>
31 #include <deal.II/grid/tria.h>
32 #include <deal.II/grid/tria_accessor.h>
33 #include <deal.II/grid/tria_iterator.h>
34 #include <deal.II/grid/tria_levels.h>
35 
36 #include <deal.II/hp/dof_faces.h>
37 #include <deal.II/hp/dof_handler.h>
38 #include <deal.II/hp/dof_level.h>
39 
40 #include <boost/serialization/array.hpp>
41 
42 #include <algorithm>
43 #include <functional>
44 #include <set>
45 #include <unordered_set>
46 
47 DEAL_II_NAMESPACE_OPEN
48 
49 // The following is necessary for compilation under Visual Studio which is
50 // unable to correctly distinguish between ::DoFHandler and
51 // ::hp::DoFHandler. Plus it makes code in dof_handler.cc easier to read.
52 #if defined(_MSC_VER) && (_MSC_VER >= 1800)
53 template <int dim, int spacedim>
54 using HpDoFHandler = ::hp::DoFHandler<dim, spacedim>;
55 #else
56 // When using older Visual Studio or a different compiler just fall back.
57 # define HpDoFHandler DoFHandler
58 #endif
59 
60 namespace parallel
61 {
62  namespace distributed
63  {
64  template <int, int>
65  class Triangulation;
66  }
67 } // namespace parallel
68 
69 
70 
71 namespace internal
72 {
73  namespace hp
74  {
75  namespace DoFHandlerImplementation
76  {
77  // access class ::hp::DoFHandler instead of namespace
78  // internal::hp::DoFHandler, etc
79  using ::hp::DoFHandler;
80 
86  {
92  template <int dim, int spacedim>
93  static void
95  DoFHandler<dim, spacedim> &dof_handler)
96  {
97  (void)dof_handler;
98  for (const auto &cell : dof_handler.active_cell_iterators())
99  if (cell->is_locally_owned())
100  Assert(
101  !cell->future_fe_index_set(),
102  ExcMessage(
103  "There shouldn't be any cells flagged for p-adaptation when partitioning."));
104  }
105 
106 
107 
112  template <int dim, int spacedim>
113  static void
115  {
116  // Release all space except the fields for active_fe_indices and
117  // refinement flags which we have to back up before
118  {
119  std::vector<std::vector<DoFLevel::active_fe_index_type>>
120  active_fe_backup(dof_handler.levels.size()),
121  future_fe_backup(dof_handler.levels.size());
122  for (unsigned int level = 0; level < dof_handler.levels.size();
123  ++level)
124  {
125  active_fe_backup[level] =
126  std::move(dof_handler.levels[level]->active_fe_indices);
127  future_fe_backup[level] =
128  std::move(dof_handler.levels[level]->future_fe_indices);
129  }
130 
131  // delete all levels and set them up newly, since vectors
132  // are troublesome if you want to change their size
133  dof_handler.clear_space();
134 
135  for (unsigned int level = 0; level < dof_handler.tria->n_levels();
136  ++level)
137  {
138  dof_handler.levels.emplace_back(new internal::hp::DoFLevel);
139  // recover backups
140  dof_handler.levels[level]->active_fe_indices =
141  std::move(active_fe_backup[level]);
142  dof_handler.levels[level]->future_fe_indices =
143  std::move(future_fe_backup[level]);
144  }
145 
146  if (dim > 1)
147  dof_handler.faces =
148  std_cxx14::make_unique<internal::hp::DoFIndicesOnFaces<dim>>();
149  }
150  }
151 
152 
153 
158  template <int dim, int spacedim>
159  static void
161  {
162  // The final step in all of the reserve_space() functions is to set
163  // up vertex dof information. since vertices are sequentially
164  // numbered, what we do first is to set up an array in which
165  // we record whether a vertex is associated with any of the
166  // given fe's, by setting a bit. in a later step, we then
167  // actually allocate memory for the required dofs
168  //
169  // in the following, we only need to consider vertices that are
170  // adjacent to either a locally owned or a ghost cell; we never
171  // store anything on vertices that are only surrounded by
172  // artificial cells. so figure out that subset of vertices
173  // first
174  std::vector<bool> locally_used_vertices(
175  dof_handler.tria->n_vertices(), false);
176  for (const auto &cell : dof_handler.active_cell_iterators())
177  if (!cell->is_artificial())
178  for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
179  locally_used_vertices[cell->vertex_index(v)] = true;
180 
181  std::vector<std::vector<bool>> vertex_fe_association(
182  dof_handler.fe_collection.size(),
183  std::vector<bool>(dof_handler.tria->n_vertices(), false));
184 
185  for (const auto &cell : dof_handler.active_cell_iterators())
186  if (!cell->is_artificial())
187  for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
188  vertex_fe_association[cell->active_fe_index()]
189  [cell->vertex_index(v)] = true;
190 
191  // in debug mode, make sure that each vertex is associated
192  // with at least one fe (note that except for unused
193  // vertices, all vertices are actually active). this is of
194  // course only true for vertices that are part of either
195  // ghost or locally owned cells
196 #ifdef DEBUG
197  for (unsigned int v = 0; v < dof_handler.tria->n_vertices(); ++v)
198  if (locally_used_vertices[v] == true)
199  if (dof_handler.tria->vertex_used(v) == true)
200  {
201  unsigned int fe = 0;
202  for (; fe < dof_handler.fe_collection.size(); ++fe)
203  if (vertex_fe_association[fe][v] == true)
204  break;
205  Assert(fe != dof_handler.fe_collection.size(),
206  ExcInternalError());
207  }
208 #endif
209 
210  // next count how much memory we actually need. for each
211  // vertex, we need one slot per fe to store the fe_index,
212  // plus dofs_per_vertex for this fe. in addition, we need
213  // one slot as the end marker for the fe_indices. at the
214  // same time already fill the vertex_dof_offsets field
215  dof_handler.vertex_dof_offsets.resize(dof_handler.tria->n_vertices(),
217 
218  unsigned int vertex_slots_needed = 0;
219  for (unsigned int v = 0; v < dof_handler.tria->n_vertices(); ++v)
220  if (dof_handler.tria->vertex_used(v) == true)
221  if (locally_used_vertices[v] == true)
222  {
223  dof_handler.vertex_dof_offsets[v] = vertex_slots_needed;
224 
225  for (unsigned int fe = 0;
226  fe < dof_handler.fe_collection.size();
227  ++fe)
228  if (vertex_fe_association[fe][v] == true)
229  vertex_slots_needed +=
230  dof_handler.get_fe(fe).dofs_per_vertex + 1;
231 
232  // don't forget the end_marker:
233  ++vertex_slots_needed;
234  }
235 
236  // now allocate the space we have determined we need, and
237  // set up the linked lists for each of the vertices
238  dof_handler.vertex_dofs.resize(vertex_slots_needed,
240  for (unsigned int v = 0; v < dof_handler.tria->n_vertices(); ++v)
241  if (dof_handler.tria->vertex_used(v) == true)
242  if (locally_used_vertices[v] == true)
243  {
244  unsigned int current_index =
245  dof_handler.vertex_dof_offsets[v];
246  for (unsigned int fe = 0;
247  fe < dof_handler.fe_collection.size();
248  ++fe)
249  if (vertex_fe_association[fe][v] == true)
250  {
251  // if this vertex uses this fe, then set the
252  // fe_index and move the pointer ahead
253  dof_handler.vertex_dofs[current_index] = fe;
254  current_index +=
255  dof_handler.get_fe(fe).dofs_per_vertex + 1;
256  }
257  // finally place the end marker
258  dof_handler.vertex_dofs[current_index] =
260  }
261  }
262 
263 
264 
269  template <int dim, int spacedim>
270  static void
272  {
273  // count how much space we need on each level for the cell
274  // dofs and set the dof_*_offsets data. initially set the
275  // latter to an invalid index, and only later set it to
276  // something reasonable for active dof_handler.cells
277  //
278  // note that for dof_handler.cells, the situation is simpler
279  // than for other (lower dimensional) objects since exactly
280  // one finite element is used for it
281  for (unsigned int level = 0; level < dof_handler.tria->n_levels();
282  ++level)
283  {
284  dof_handler.levels[level]->dof_offsets =
285  std::vector<DoFLevel::offset_type>(
286  dof_handler.tria->n_raw_cells(level),
287  static_cast<DoFLevel::offset_type>(-1));
288  dof_handler.levels[level]->cell_cache_offsets =
289  std::vector<DoFLevel::offset_type>(
290  dof_handler.tria->n_raw_cells(level),
291  static_cast<DoFLevel::offset_type>(-1));
292 
293  types::global_dof_index next_free_dof = 0;
294  types::global_dof_index cache_size = 0;
295  typename HpDoFHandler<dim, spacedim>::active_cell_iterator
296  cell = dof_handler.begin_active(level),
297  endc = dof_handler.end_active(level);
298  for (; cell != endc; ++cell)
299  if (cell->is_active() && !cell->is_artificial())
300  {
301  dof_handler.levels[level]->dof_offsets[cell->index()] =
302  next_free_dof;
303  next_free_dof +=
304  cell->get_fe().template n_dofs_per_object<dim>();
305 
306  dof_handler.levels[level]
307  ->cell_cache_offsets[cell->index()] = cache_size;
308  cache_size += cell->get_fe().dofs_per_cell;
309  }
310 
311  dof_handler.levels[level]->dof_indices =
312  std::vector<types::global_dof_index>(
313  next_free_dof, numbers::invalid_dof_index);
314  dof_handler.levels[level]->cell_dof_indices_cache =
315  std::vector<types::global_dof_index>(
316  cache_size, numbers::invalid_dof_index);
317  }
318 
319  // safety check: make sure that the number of DoFs we
320  // allocated is actually correct (above we have also set the
321  // dof_*_offsets field, so we couldn't use this simpler
322  // algorithm)
323 #ifdef DEBUG
324  for (unsigned int level = 0; level < dof_handler.tria->n_levels();
325  ++level)
326  {
327  types::global_dof_index counter = 0;
328  for (const auto &cell :
329  dof_handler.active_cell_iterators_on_level(level))
330  if (cell->is_active() && !cell->is_artificial())
331  counter += cell->get_fe().template n_dofs_per_object<dim>();
332 
333  Assert(dof_handler.levels[level]->dof_indices.size() == counter,
334  ExcInternalError());
335 
336  // also check that the number of unassigned slots in the
337  // dof_offsets equals the number of cells on that level minus the
338  // number of active, non-artificial cells (because these are
339  // exactly the cells on which we do something)
340  unsigned int n_active_non_artificial_cells = 0;
341  for (const auto &cell :
342  dof_handler.active_cell_iterators_on_level(level))
343  if (cell->is_active() && !cell->is_artificial())
344  ++n_active_non_artificial_cells;
345 
346  Assert(static_cast<unsigned int>(std::count(
347  dof_handler.levels[level]->dof_offsets.begin(),
348  dof_handler.levels[level]->dof_offsets.end(),
349  static_cast<DoFLevel::offset_type>(-1))) ==
350  dof_handler.tria->n_raw_cells(level) -
351  n_active_non_artificial_cells,
352  ExcInternalError());
353  }
354 #endif
355  }
356 
357 
358 
363  template <int dim, int spacedim>
364  static void
366  {
367  // make the code generic between lines and quads
368  std::vector<unsigned int> &face_dof_offsets =
369  (dim == 2 ?
370  reinterpret_cast<::internal::hp::DoFIndicesOnFaces<2> &>(
371  *dof_handler.faces)
372  .lines.dof_offsets :
373  reinterpret_cast<::internal::hp::DoFIndicesOnFaces<3> &>(
374  *dof_handler.faces)
375  .quads.dof_offsets);
376 
377  std::vector<types::global_dof_index> &face_dof_indices =
378  (dim == 2 ?
379  reinterpret_cast<::internal::hp::DoFIndicesOnFaces<2> &>(
380  *dof_handler.faces)
381  .lines.dofs :
382  reinterpret_cast<::internal::hp::DoFIndicesOnFaces<3> &>(
383  *dof_handler.faces)
384  .quads.dofs);
385 
386  // FACE DOFS
387  //
388  // Count face dofs, then allocate as much space
389  // as we need and prime the linked list for faces (see the
390  // description in hp::DoFLevel) with the indices we will
391  // need. Note that our task is more complicated than for the
392  // cell case above since two adjacent cells may have different
393  // active_fe_indices, in which case we need to allocate
394  // *two* sets of face dofs for the same face. But they don't
395  // *have* to be different, and so we need to prepare for this
396  // as well.
397  //
398  // The way we do things is that we loop over all active
399  // cells (these are the only ones that have DoFs
400  // anyway) and all their faces. We note in the
401  // user flags whether we have previously visited a face and
402  // if so skip it (consequently, we have to save and later
403  // restore the face flags)
404  {
405  std::vector<bool> saved_face_user_flags;
406  switch (dim)
407  {
408  case 2:
409  {
410  const_cast<::Triangulation<dim, spacedim> &>(
411  *dof_handler.tria)
412  .save_user_flags_line(saved_face_user_flags);
413  const_cast<::Triangulation<dim, spacedim> &>(
414  *dof_handler.tria)
415  .clear_user_flags_line();
416 
417  break;
418  }
419 
420  case 3:
421  {
422  const_cast<::Triangulation<dim, spacedim> &>(
423  *dof_handler.tria)
424  .save_user_flags_quad(saved_face_user_flags);
425  const_cast<::Triangulation<dim, spacedim> &>(
426  *dof_handler.tria)
427  .clear_user_flags_quad();
428 
429  break;
430  }
431 
432  default:
433  Assert(false, ExcNotImplemented());
434  }
435 
436  // An array to hold how many slots (see the hp::DoFLevel
437  // class) we will have to store on each level
438  unsigned int n_face_slots = 0;
439 
440  for (const auto &cell : dof_handler.active_cell_iterators())
441  if (!cell->is_artificial())
442  for (const unsigned int face :
444  if (cell->face(face)->user_flag_set() == false)
445  {
446  // Ok, face has not been visited. So we need to
447  // allocate space for it. Let's see how much we
448  // need: we need one set if a) there is no
449  // neighbor behind this face, or b) the neighbor
450  // is either coarser or finer than we are, or c)
451  // the neighbor is artificial, or d) the neighbor
452  // is neither coarser nor finer, nor is artificial,
453  // and just so happens to have the same active_fe_index :
454  if (cell->at_boundary(face) ||
455  cell->face(face)->has_children() ||
456  cell->neighbor_is_coarser(face) ||
457  (!cell->at_boundary(face) &&
458  cell->neighbor(face)->is_artificial()) ||
459  (!cell->at_boundary(face) &&
460  !cell->neighbor(face)->is_artificial() &&
461  (cell->active_fe_index() ==
462  cell->neighbor(face)->active_fe_index())))
463  // Ok, one set of dofs. that makes one active_fe_index,
464  // 1 times dofs_per_face dofs, and one stop index
465  n_face_slots +=
466  1 +
467  dof_handler.get_fe(cell->active_fe_index())
468  .template n_dofs_per_object<dim - 1>() +
469  1;
470 
471  // Otherwise we do indeed need two sets, i.e. two
472  // active_fe_indices, two sets of dofs, and one stop
473  // index:
474  else
475  n_face_slots +=
476  (1 + // the active_fe_index
477  dof_handler.get_fe(cell->active_fe_index())
478  .template n_dofs_per_object<
479  dim - 1>() // actual DoF indices
480  + 1 // the second active_fe_index
481  + dof_handler
482  .get_fe(cell->neighbor(face)->active_fe_index())
483  .template n_dofs_per_object<
484  dim - 1>() // actual DoF indices
485  + 1); // stop marker
486 
487  // mark this face as visited
488  cell->face(face)->set_user_flag();
489  }
490 
491  // Now that we know how many face dofs we will have to
492  // have, allocate the memory. Note that we
493  // allocate offsets for all faces, though only the active
494  // ones will have a non-invalid value later on
495  face_dof_offsets =
496  std::vector<unsigned int>(dof_handler.tria->n_raw_faces(),
498  face_dof_indices =
499  std::vector<types::global_dof_index>(n_face_slots,
501 
502  // With the memory now allocated, loop over the
503  // dof_handler cells again and prime the _offset values as
504  // well as the fe_index fields
505  switch (dim)
506  {
507  case 2:
508  {
509  const_cast<::Triangulation<dim, spacedim> &>(
510  *dof_handler.tria)
511  .clear_user_flags_line();
512 
513  break;
514  }
515 
516  case 3:
517  {
518  const_cast<::Triangulation<dim, spacedim> &>(
519  *dof_handler.tria)
520  .clear_user_flags_quad();
521 
522  break;
523  }
524 
525  default:
526  Assert(false, ExcNotImplemented());
527  }
528 
529  unsigned int next_free_face_slot = 0;
530 
531  for (const auto &cell : dof_handler.active_cell_iterators())
532  if (!cell->is_artificial())
533  for (const unsigned int face :
535  if (!cell->face(face)->user_flag_set())
536  {
537  // Same decision tree as before
538  if (cell->at_boundary(face) ||
539  cell->face(face)->has_children() ||
540  cell->neighbor_is_coarser(face) ||
541  (!cell->at_boundary(face) &&
542  cell->neighbor(face)->is_artificial()) ||
543  (!cell->at_boundary(face) &&
544  !cell->neighbor(face)->is_artificial() &&
545  (cell->active_fe_index() ==
546  cell->neighbor(face)->active_fe_index())))
547  {
548  // Only one active_fe_index lives on this
549  // face
550  face_dof_offsets[cell->face(face)->index()] =
551  next_free_face_slot;
552 
553  // Set the first and only slot for this face to
554  // active_fe_index of this face
555  face_dof_indices[next_free_face_slot] =
556  cell->active_fe_index();
557 
558  // The next dofs_per_face indices remain unset
559  // for the moment (i.e. at invalid_dof_index).
560  // Following this comes the stop index, which
561  // also is invalid_dof_index and therefore
562  // does not have to be explicitly set
563 
564  // Finally, move the current marker forward:
565  next_free_face_slot +=
566  1 // the active_fe_index
567  + dof_handler.get_fe(cell->active_fe_index())
568  .template n_dofs_per_object<
569  dim - 1>() // actual DoF indices
570  + 1; // the end marker
571  }
572  else
573  {
574  // There are two active_fe_indices that live on this
575  // face.
576  face_dof_offsets[cell->face(face)->index()] =
577  next_free_face_slot;
578 
579  // Store the two indices we will have to deal with.
580  // We sort these so that it does not matter which of
581  // the cells adjacent to this face we visit first. In
582  // sequential computations, this does not matter
583  // because the order in which we visit these cells is
584  // deterministic and always the same. But in parallel
585  // computations, we can get into trouble because two
586  // processors visit the cells in different order
587  // (because the mesh creation history on the two
588  // processors is different), and in that case it can
589  // happen that the order of active_fe_index values for
590  // a given face is different on the two processes,
591  // even though they agree on which two values need to
592  // be stored. Since the DoF unification on faces
593  // takes into account the order of the
594  // active_fe_indices, this leads to quite subtle bugs.
595  // We could fix this in the place where we do the DoF
596  // unification on cells, but it is better to just make
597  // sure that every process stores the exact same
598  // information (and in the same order) on each face.
599  unsigned int active_fe_indices[2] = {
600  cell->active_fe_index(),
601  cell->neighbor(face)->active_fe_index()};
602  if (active_fe_indices[1] < active_fe_indices[0])
603  std::swap(active_fe_indices[0],
604  active_fe_indices[1]);
605 
606  // Set first slot for this face to
607  // active_fe_index of this face
608  face_dof_indices[next_free_face_slot] =
609  active_fe_indices[0];
610 
611  // The next dofs_per_line/quad indices remain unset
612  // for the moment (i.e. at invalid_dof_index).
613  //
614  // Then comes the fe_index for the second slot:
615  face_dof_indices
616  [next_free_face_slot + 1 +
617  dof_handler.get_fe(active_fe_indices[0])
618  .template n_dofs_per_object<dim - 1>()] =
619  active_fe_indices[1];
620  // Then again a set of dofs that we need not
621  // set right now
622  //
623  // Following this comes the stop index, which
624  // also is invalid_dof_index and therefore
625  // does not have to be explicitly set
626 
627  // Finally, move the current marker forward:
628  next_free_face_slot +=
629  (1 // first active_fe_index
630  + dof_handler.get_fe(active_fe_indices[0])
631  .template n_dofs_per_object<
632  dim - 1>() // actual DoF indices
633  + 1 // second active_fe_index
634  + dof_handler.get_fe(active_fe_indices[1])
635  .template n_dofs_per_object<
636  dim - 1>() // actual DoF indices
637  + 1); // end marker
638  }
639 
640  // mark this face as visited
641  cell->face(face)->set_user_flag();
642  }
643 
644  // we should have moved the cursor for each level to the
645  // total number of dofs on that level. check that
646  Assert(next_free_face_slot == n_face_slots, ExcInternalError());
647 
648  // at the end, restore the user flags for the faces
649  switch (dim)
650  {
651  case 2:
652  {
653  const_cast<::Triangulation<dim, spacedim> &>(
654  *dof_handler.tria)
655  .load_user_flags_line(saved_face_user_flags);
656 
657  break;
658  }
659 
660  case 3:
661  {
662  const_cast<::Triangulation<dim, spacedim> &>(
663  *dof_handler.tria)
664  .load_user_flags_quad(saved_face_user_flags);
665 
666  break;
667  }
668 
669  default:
670  Assert(false, ExcNotImplemented());
671  }
672  }
673  }
674 
675 
676 
683  template <int spacedim>
684  static void reserve_space(DoFHandler<1, spacedim> &dof_handler)
685  {
686  Assert(dof_handler.fe_collection.size() > 0,
688  Assert(dof_handler.tria->n_levels() > 0,
689  ExcMessage("The current Triangulation must not be empty."));
690  Assert(dof_handler.tria->n_levels() == dof_handler.levels.size(),
691  ExcInternalError());
692 
693  reserve_space_release_space(dof_handler);
694 
695  Threads::TaskGroup<> tasks;
696  tasks +=
697  Threads::new_task(&reserve_space_cells<1, spacedim>, dof_handler);
698  tasks += Threads::new_task(&reserve_space_vertices<1, spacedim>,
699  dof_handler);
700  tasks.join_all();
701  }
702 
703 
704 
705  template <int spacedim>
706  static void reserve_space(DoFHandler<2, spacedim> &dof_handler)
707  {
708  Assert(dof_handler.fe_collection.size() > 0,
710  Assert(dof_handler.tria->n_levels() > 0,
711  ExcMessage("The current Triangulation must not be empty."));
712  Assert(dof_handler.tria->n_levels() == dof_handler.levels.size(),
713  ExcInternalError());
714 
715  reserve_space_release_space(dof_handler);
716 
717  Threads::TaskGroup<> tasks;
718  tasks +=
719  Threads::new_task(&reserve_space_cells<2, spacedim>, dof_handler);
720  tasks +=
721  Threads::new_task(&reserve_space_faces<2, spacedim>, dof_handler);
722  tasks += Threads::new_task(&reserve_space_vertices<2, spacedim>,
723  dof_handler);
724  tasks.join_all();
725  }
726 
727 
728 
729  template <int spacedim>
730  static void reserve_space(DoFHandler<3, spacedim> &dof_handler)
731  {
732  const unsigned int dim = 3;
733 
734  Assert(dof_handler.fe_collection.size() > 0,
736  Assert(dof_handler.tria->n_levels() > 0,
737  ExcMessage("The current Triangulation must not be empty."));
738  Assert(dof_handler.tria->n_levels() == dof_handler.levels.size(),
739  ExcInternalError());
740 
741  reserve_space_release_space(dof_handler);
742 
743  Threads::TaskGroup<> tasks;
744  tasks +=
745  Threads::new_task(&reserve_space_cells<3, spacedim>, dof_handler);
746  tasks +=
747  Threads::new_task(&reserve_space_faces<3, spacedim>, dof_handler);
748  tasks += Threads::new_task(&reserve_space_vertices<3, spacedim>,
749  dof_handler);
750 
751  // While the tasks above are running, we can turn to line dofs
752 
753  // the situation here is pretty much like with vertices:
754  // there can be an arbitrary number of finite elements
755  // associated with each line.
756  //
757  // the algorithm we use is somewhat similar to what we do in
758  // reserve_space_vertices()
759  {
760  // what we do first is to set up an array in which we
761  // record whether a line is associated with any of the
762  // given fe's, by setting a bit. in a later step, we
763  // then actually allocate memory for the required dofs
764  std::vector<std::vector<bool>> line_fe_association(
765  dof_handler.fe_collection.size(),
766  std::vector<bool>(dof_handler.tria->n_raw_lines(), false));
767 
768  for (const auto &cell : dof_handler.active_cell_iterators())
769  if (!cell->is_artificial())
770  for (unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell;
771  ++l)
772  line_fe_association[cell->active_fe_index()]
773  [cell->line_index(l)] = true;
774 
775  // first check which of the lines is used at all,
776  // i.e. is associated with a finite element. we do this
777  // since not all lines may actually be used, in which
778  // case we do not have to allocate any memory at all
779  std::vector<bool> line_is_used(dof_handler.tria->n_raw_lines(),
780  false);
781  for (unsigned int line = 0; line < dof_handler.tria->n_raw_lines();
782  ++line)
783  for (unsigned int fe = 0; fe < dof_handler.fe_collection.size();
784  ++fe)
785  if (line_fe_association[fe][line] == true)
786  {
787  line_is_used[line] = true;
788  break;
789  }
790 
791  // next count how much memory we actually need. for each
792  // line, we need one slot per fe to store the fe_index,
793  // plus dofs_per_line for this fe. in addition, we need
794  // one slot as the end marker for the fe_indices. at the
795  // same time already fill the line_dofs_offsets field
796  dof_handler.faces->lines.dof_offsets.resize(
797  dof_handler.tria->n_raw_lines(), numbers::invalid_unsigned_int);
798 
799  unsigned int line_slots_needed = 0;
800  for (unsigned int line = 0; line < dof_handler.tria->n_raw_lines();
801  ++line)
802  if (line_is_used[line] == true)
803  {
804  dof_handler.faces->lines.dof_offsets[line] =
805  line_slots_needed;
806 
807  for (unsigned int fe = 0;
808  fe < dof_handler.fe_collection.size();
809  ++fe)
810  if (line_fe_association[fe][line] == true)
811  line_slots_needed +=
812  dof_handler.get_fe(fe).dofs_per_line + 1;
813  ++line_slots_needed;
814  }
815 
816  // now allocate the space we have determined we need,
817  // and set up the linked lists for each of the lines
818  dof_handler.faces->lines.dofs.resize(line_slots_needed,
820  for (unsigned int line = 0; line < dof_handler.tria->n_raw_lines();
821  ++line)
822  if (line_is_used[line] == true)
823  {
824  unsigned int pointer =
825  dof_handler.faces->lines.dof_offsets[line];
826  for (unsigned int fe = 0;
827  fe < dof_handler.fe_collection.size();
828  ++fe)
829  if (line_fe_association[fe][line] == true)
830  {
831  // if this line uses this fe, then set the
832  // fe_index and move the pointer ahead
833  dof_handler.faces->lines.dofs[pointer] = fe;
834  pointer += dof_handler.get_fe(fe).dofs_per_line + 1;
835  }
836  // finally place the end marker
837  dof_handler.faces->lines.dofs[pointer] =
839  }
840  }
841 
842  // Ensure that everything is done at this point.
843  tasks.join_all();
844  }
845 
846 
847 
851  template <int spacedim>
852  static unsigned int
854  {
855  return std::min(static_cast<types::global_dof_index>(
856  3 *
857  dof_handler.fe_collection.max_dofs_per_vertex() +
858  2 * dof_handler.fe_collection.max_dofs_per_line()),
859  dof_handler.n_dofs());
860  }
861 
862 
863 
864  template <int spacedim>
865  static unsigned int
866  max_couplings_between_dofs(const DoFHandler<2, spacedim> &dof_handler)
867  {
868  // get these numbers by drawing pictures and counting...
869  // example:
870  // | | |
871  // --x-----x--x--X--
872  // | | | |
873  // | x--x--x
874  // | | | |
875  // --x--x--*--x--x--
876  // | | | |
877  // x--x--x |
878  // | | | |
879  // --X--x--x-----x--
880  // | | |
881  // x = vertices connected with center vertex *;
882  // = total of 19
883  //
884  // (the X vertices are connected with * if the vertices
885  // adjacent to X are hanging nodes) count lines -> 28 (don't
886  // forget to count mother and children separately!)
887  types::global_dof_index max_couplings;
888  switch (dof_handler.tria->max_adjacent_cells())
889  {
890  case 4:
891  max_couplings =
892  19 * dof_handler.fe_collection.max_dofs_per_vertex() +
893  28 * dof_handler.fe_collection.max_dofs_per_line() +
894  8 * dof_handler.fe_collection.max_dofs_per_quad();
895  break;
896  case 5:
897  max_couplings =
898  21 * dof_handler.fe_collection.max_dofs_per_vertex() +
899  31 * dof_handler.fe_collection.max_dofs_per_line() +
900  9 * dof_handler.fe_collection.max_dofs_per_quad();
901  break;
902  case 6:
903  max_couplings =
904  28 * dof_handler.fe_collection.max_dofs_per_vertex() +
905  42 * dof_handler.fe_collection.max_dofs_per_line() +
906  12 * dof_handler.fe_collection.max_dofs_per_quad();
907  break;
908  case 7:
909  max_couplings =
910  30 * dof_handler.fe_collection.max_dofs_per_vertex() +
911  45 * dof_handler.fe_collection.max_dofs_per_line() +
912  13 * dof_handler.fe_collection.max_dofs_per_quad();
913  break;
914  case 8:
915  max_couplings =
916  37 * dof_handler.fe_collection.max_dofs_per_vertex() +
917  56 * dof_handler.fe_collection.max_dofs_per_line() +
918  16 * dof_handler.fe_collection.max_dofs_per_quad();
919  break;
920  default:
921  Assert(false, ExcNotImplemented());
922  max_couplings = 0;
923  }
924  return std::min(max_couplings, dof_handler.n_dofs());
925  }
926 
927 
928 
929  template <int spacedim>
930  static unsigned int
931  max_couplings_between_dofs(const DoFHandler<3, spacedim> &dof_handler)
932  {
933  // TODO:[?] Invent significantly better estimates than the ones in
934  // this function
935  // doing the same thing here is a rather complicated thing,
936  // compared to the 2d case, since it is hard to draw
937  // pictures with several refined hexahedra :-) so I
938  // presently only give a coarse estimate for the case that
939  // at most 8 hexes meet at each vertex
940  //
941  // can anyone give better estimate here?
942  const unsigned int max_adjacent_cells =
943  dof_handler.tria->max_adjacent_cells();
944 
945  types::global_dof_index max_couplings;
946  if (max_adjacent_cells <= 8)
947  max_couplings =
948  7 * 7 * 7 * dof_handler.fe_collection.max_dofs_per_vertex() +
949  7 * 6 * 7 * 3 * dof_handler.fe_collection.max_dofs_per_line() +
950  9 * 4 * 7 * 3 * dof_handler.fe_collection.max_dofs_per_quad() +
951  27 * dof_handler.fe_collection.max_dofs_per_hex();
952  else
953  {
954  Assert(false, ExcNotImplemented());
955  max_couplings = 0;
956  }
957 
958  return std::min(max_couplings, dof_handler.n_dofs());
959  }
960 
961 
962 
974  template <int dim, int spacedim>
975  static void
977  ::hp::DoFHandler<dim, spacedim> &dof_handler)
978  {
979  if (const ::parallel::shared::Triangulation<dim, spacedim> *tr =
980  dynamic_cast<
981  const ::parallel::shared::Triangulation<dim, spacedim>
982  *>(&dof_handler.get_triangulation()))
983  {
984  // we have a shared triangulation. in this case, every processor
985  // knows about all cells, but every processor only has knowledge
986  // about the active_fe_index on the cells it owns.
987  //
988  // we can create a complete set of active_fe_indices by letting
989  // every processor create a vector of indices for all cells,
990  // filling only those on the cells it owns and setting the indices
991  // on the other cells to zero. then we add all of these vectors
992  // up, and because every vector entry has exactly one processor
993  // that owns it, the sum is correct
994  std::vector<unsigned int> active_fe_indices(tr->n_active_cells(),
995  0u);
996  for (const auto &cell : dof_handler.active_cell_iterators())
997  if (cell->is_locally_owned())
998  active_fe_indices[cell->active_cell_index()] =
999  cell->active_fe_index();
1000 
1001  Utilities::MPI::sum(active_fe_indices,
1002  tr->get_communicator(),
1003  active_fe_indices);
1004 
1005  // now go back and fill the active_fe_index on all other
1006  // cells. we would like to call cell->set_active_fe_index(),
1007  // but that function does not allow setting these indices on
1008  // non-locally_owned cells. so we have to work around the
1009  // issue a little bit by accessing the underlying data
1010  // structures directly
1011  for (const auto &cell : dof_handler.active_cell_iterators())
1012  if (!cell->is_locally_owned())
1013  dof_handler.levels[cell->level()]->set_active_fe_index(
1014  cell->index(),
1015  active_fe_indices[cell->active_cell_index()]);
1016  }
1017  else if (const ::parallel::distributed::Triangulation<dim,
1018  spacedim>
1019  *tr = dynamic_cast<const ::parallel::distributed::
1020  Triangulation<dim, spacedim> *>(
1021  &dof_handler.get_triangulation()))
1022  {
1023  // For completely distributed meshes, use the function that is
1024  // able to move data from locally owned cells on one processor to
1025  // the corresponding ghost cells on others. To this end, we need
1026  // to have functions that can pack and unpack the data we want to
1027  // transport -- namely, the single unsigned int active_fe_index
1028  // objects
1029  auto pack =
1030  [](const typename ::hp::DoFHandler<dim, spacedim>::
1031  active_cell_iterator &cell) -> unsigned int {
1032  return cell->active_fe_index();
1033  };
1034 
1035  auto unpack =
1036  [&dof_handler](
1037  const typename ::hp::DoFHandler<dim, spacedim>::
1038  active_cell_iterator &cell,
1039  const unsigned int active_fe_index) -> void {
1040  // we would like to say
1041  // cell->set_active_fe_index(active_fe_index);
1042  // but this is not allowed on cells that are not
1043  // locally owned, and we are on a ghost cell
1044  dof_handler.levels[cell->level()]->set_active_fe_index(
1045  cell->index(), active_fe_index);
1046  };
1047 
1049  unsigned int,
1050  ::hp::DoFHandler<dim, spacedim>>(dof_handler,
1051  pack,
1052  unpack);
1053  }
1054  else
1055  {
1056  // a sequential triangulation. there is nothing we need to do here
1057  Assert(
1058  (dynamic_cast<
1059  const ::parallel::TriangulationBase<dim, spacedim> *>(
1060  &dof_handler.get_triangulation()) == nullptr),
1061  ExcInternalError());
1062  }
1063  }
1064 
1065 
1066 
1087  template <int dim, int spacedim>
1088  static void
1090  ::hp::DoFHandler<dim, spacedim> &dof_handler)
1091  {
1092  const auto &fe_transfer = dof_handler.active_fe_index_transfer;
1093 
1094  for (const auto &cell : dof_handler.active_cell_iterators())
1095  if (cell->is_locally_owned())
1096  {
1097  if (cell->refine_flag_set())
1098  {
1099  // Store the active_fe_index of each cell that will be
1100  // refined to and distribute it later on its children.
1101  // Pick their future index if flagged for p-refinement.
1102  fe_transfer->refined_cells_fe_index.insert(
1103  {cell, cell->future_fe_index()});
1104  }
1105  else if (cell->coarsen_flag_set())
1106  {
1107  // From all cells that will be coarsened, determine their
1108  // parent and calculate its proper active_fe_index, so that
1109  // it can be set after refinement. But first, check if that
1110  // particular cell has a parent at all.
1111  Assert(cell->level() > 0, ExcInternalError());
1112  const auto &parent = cell->parent();
1113 
1114  // Check if the active_fe_index for the current cell has
1115  // been determined already.
1116  if (fe_transfer->coarsened_cells_fe_index.find(parent) ==
1117  fe_transfer->coarsened_cells_fe_index.end())
1118  {
1119  // Find a suitable active_fe_index for the parent cell
1120  // based on the 'least dominant finite element' of its
1121  // children. Consider the childrens' hypothetical future
1122  // index when they have been flagged for p-refinement.
1123  std::set<unsigned int> fe_indices_children;
1124  for (unsigned int child_index = 0;
1125  child_index < parent->n_children();
1126  ++child_index)
1127  {
1128  const auto sibling = parent->child(child_index);
1129  Assert(sibling->is_active() &&
1130  sibling->coarsen_flag_set(),
1131  typename ::Triangulation<
1132  dim>::ExcInconsistentCoarseningFlags());
1133 
1134  fe_indices_children.insert(
1135  sibling->future_fe_index());
1136  }
1137  Assert(!fe_indices_children.empty(),
1138  ExcInternalError());
1139 
1140  const unsigned int fe_index =
1141  dof_handler.fe_collection.find_dominated_fe_extended(
1142  fe_indices_children, /*codim=*/0);
1143 
1145  typename ::hp::FECollection<dim>::
1146  ExcNoDominatedFiniteElementAmongstChildren());
1147 
1148  fe_transfer->coarsened_cells_fe_index.insert(
1149  {parent, fe_index});
1150  }
1151  }
1152  else
1153  {
1154  // No h-refinement is scheduled for this cell.
1155  // However, it may have p-refinement indicators, so we
1156  // choose a new active_fe_index based on its flags.
1157  if (cell->future_fe_index_set() == true)
1158  fe_transfer->persisting_cells_fe_index.insert(
1159  {cell, cell->future_fe_index()});
1160  }
1161  }
1162  }
1163 
1164 
1165 
1170  template <int dim, int spacedim>
1171  static void
1173  ::hp::DoFHandler<dim, spacedim> &dof_handler)
1174  {
1175  const auto &fe_transfer = dof_handler.active_fe_index_transfer;
1176 
1177  // Set active_fe_indices on persisting cells.
1178  for (const auto &persist : fe_transfer->persisting_cells_fe_index)
1179  {
1180  const auto &cell = persist.first;
1181 
1182  if (cell->is_locally_owned())
1183  {
1184  Assert(cell->is_active(), ExcInternalError());
1185  cell->set_active_fe_index(persist.second);
1186  }
1187  }
1188 
1189  // Distribute active_fe_indices from all refined cells on their
1190  // respective children.
1191  for (const auto &refine : fe_transfer->refined_cells_fe_index)
1192  {
1193  const auto &parent = refine.first;
1194 
1195  for (unsigned int child_index = 0;
1196  child_index < parent->n_children();
1197  ++child_index)
1198  {
1199  const auto &child = parent->child(child_index);
1200  Assert(child->is_locally_owned() && child->is_active(),
1201  ExcInternalError());
1202  child->set_active_fe_index(refine.second);
1203  }
1204  }
1205 
1206  // Set active_fe_indices on coarsened cells that have been determined
1207  // before the actual coarsening happened.
1208  for (const auto &coarsen : fe_transfer->coarsened_cells_fe_index)
1209  {
1210  const auto &cell = coarsen.first;
1211  Assert(cell->is_locally_owned() && cell->is_active(),
1212  ExcInternalError());
1213  cell->set_active_fe_index(coarsen.second);
1214  }
1215  }
1216 
1217 
1228  template <int dim, int spacedim>
1229  static unsigned int
1231  const std::vector<unsigned int> & children_fe_indices,
1232  ::hp::FECollection<dim, spacedim> &fe_collection)
1233  {
1234  Assert(!children_fe_indices.empty(), ExcInternalError());
1235 
1236  // convert vector to set
1237  const std::set<unsigned int> children_fe_indices_set(
1238  children_fe_indices.begin(), children_fe_indices.end());
1239 
1240  const unsigned int dominated_fe_index =
1241  fe_collection.find_dominated_fe_extended(children_fe_indices_set,
1242  /*codim=*/0);
1243 
1244  Assert(dominated_fe_index != numbers::invalid_unsigned_int,
1245  typename ::hp::FECollection<
1246  dim>::ExcNoDominatedFiniteElementAmongstChildren());
1247 
1248  return dominated_fe_index;
1249  }
1250  };
1251  } // namespace DoFHandlerImplementation
1252  } // namespace hp
1253 } // namespace internal
1254 
1255 
1256 
1257 namespace hp
1258 {
1259  template <int dim, int spacedim>
1260  const unsigned int DoFHandler<dim, spacedim>::dimension;
1261 
1262  template <int dim, int spacedim>
1264 
1265 
1266 
1267  template <int dim, int spacedim>
1269  : tria(nullptr, typeid(*this).name())
1270  , faces(nullptr)
1271  {}
1272 
1273 
1274 
1275  template <int dim, int spacedim>
1278  : tria(&tria, typeid(*this).name())
1279  , faces(nullptr)
1280  {
1281  setup_policy_and_listeners();
1282 
1283  create_active_fe_table();
1284  }
1285 
1286 
1287 
1288  template <int dim, int spacedim>
1290  {
1291  // unsubscribe as a listener to refinement of the underlying
1292  // triangulation
1293  for (auto &connection : tria_listeners)
1294  connection.disconnect();
1295  tria_listeners.clear();
1296 
1297  // ...and release allocated memory
1298  // virtual functions called in constructors and destructors never use the
1299  // override in a derived class
1300  // for clarity be explicit on which function is called
1302  }
1303 
1304 
1305 
1306  /*------------------------ Cell iterator functions ------------------------*/
1307 
1308 
1309  template <int dim, int spacedim>
1311  DoFHandler<dim, spacedim>::begin(const unsigned int level) const
1312  {
1313  return cell_iterator(*this->get_triangulation().begin(level), this);
1314  }
1315 
1316 
1317 
1318  template <int dim, int spacedim>
1320  DoFHandler<dim, spacedim>::begin_active(const unsigned int level) const
1321  {
1322  // level is checked in begin
1323  cell_iterator i = begin(level);
1324  if (i.state() != IteratorState::valid)
1325  return i;
1326  while (i->has_children())
1327  if ((++i).state() != IteratorState::valid)
1328  return i;
1329  return i;
1330  }
1331 
1332 
1333 
1334  template <int dim, int spacedim>
1337  {
1338  return cell_iterator(&this->get_triangulation(), -1, -1, this);
1339  }
1340 
1341 
1342 
1343  template <int dim, int spacedim>
1345  DoFHandler<dim, spacedim>::end(const unsigned int level) const
1346  {
1347  return (level == this->get_triangulation().n_levels() - 1 ?
1348  end() :
1349  begin(level + 1));
1350  }
1351 
1352 
1353 
1354  template <int dim, int spacedim>
1356  DoFHandler<dim, spacedim>::end_active(const unsigned int level) const
1357  {
1358  return (level == this->get_triangulation().n_levels() - 1 ?
1360  begin_active(level + 1));
1361  }
1362 
1363 
1364 
1365  template <int dim, int spacedim>
1368  {
1370  begin(), end());
1371  }
1372 
1373 
1374 
1375  template <int dim, int spacedim>
1378  {
1379  return IteratorRange<
1381  end());
1382  }
1383 
1384 
1385 
1386  template <int dim, int spacedim>
1389  const unsigned int level) const
1390  {
1392  begin(level), end(level));
1393  }
1394 
1395 
1396 
1397  template <int dim, int spacedim>
1400  const unsigned int level) const
1401  {
1402  return IteratorRange<
1404  begin_active(level), end_active(level));
1405  }
1406 
1407 
1408 
1409  //------------------------------------------------------------------
1410 
1411 
1412  template <int dim, int spacedim>
1415  {
1416  Assert(fe_collection.size() > 0, ExcNoFESelected());
1417 
1418  std::unordered_set<types::global_dof_index> boundary_dofs;
1419  std::vector<types::global_dof_index> dofs_on_face;
1420  dofs_on_face.reserve(this->get_fe_collection().max_dofs_per_face());
1421 
1422  const IndexSet &owned_dofs = locally_owned_dofs();
1423 
1424  // loop over all faces to check whether they are at a
1425  // boundary. note that we need not take special care of single
1426  // lines in 3d (using @p{cell->has_boundary_lines}), since we do
1427  // not support boundaries of dimension dim-2, and so every
1428  // boundary line is also part of a boundary face.
1429  typename HpDoFHandler<dim, spacedim>::active_cell_iterator
1430  cell = this->begin_active(),
1431  endc = this->end();
1432  for (; cell != endc; ++cell)
1433  if (cell->is_locally_owned() && cell->at_boundary())
1434  {
1435  for (auto f : GeometryInfo<dim>::face_indices())
1436  if (cell->at_boundary(f))
1437  {
1438  const unsigned int dofs_per_face = cell->get_fe().dofs_per_face;
1439  dofs_on_face.resize(dofs_per_face);
1440 
1441  cell->face(f)->get_dof_indices(dofs_on_face,
1442  cell->active_fe_index());
1443  for (unsigned int i = 0; i < dofs_per_face; ++i)
1444  {
1445  const unsigned int global_idof_index = dofs_on_face[i];
1446  if (owned_dofs.is_element(global_idof_index))
1447  {
1448  boundary_dofs.insert(global_idof_index);
1449  }
1450  }
1451  }
1452  }
1453  return boundary_dofs.size();
1454  }
1455 
1456 
1457 
1458  template <int dim, int spacedim>
1461  const std::set<types::boundary_id> &boundary_ids) const
1462  {
1463  Assert(fe_collection.size() > 0, ExcNoFESelected());
1464  Assert(boundary_ids.find(numbers::internal_face_boundary_id) ==
1465  boundary_ids.end(),
1467 
1468  // same as above, but with additional checks for set of boundary
1469  // indicators
1470  std::unordered_set<types::global_dof_index> boundary_dofs;
1471  std::vector<types::global_dof_index> dofs_on_face;
1472  dofs_on_face.reserve(this->get_fe_collection().max_dofs_per_face());
1473 
1474  const IndexSet &owned_dofs = locally_owned_dofs();
1475 
1476  typename HpDoFHandler<dim, spacedim>::active_cell_iterator
1477  cell = this->begin_active(),
1478  endc = this->end();
1479  for (; cell != endc; ++cell)
1480  if (cell->is_locally_owned() && cell->at_boundary())
1481  {
1482  for (auto f : GeometryInfo<dim>::face_indices())
1483  if (cell->at_boundary(f) &&
1484  (boundary_ids.find(cell->face(f)->boundary_id()) !=
1485  boundary_ids.end()))
1486  {
1487  const unsigned int dofs_per_face = cell->get_fe().dofs_per_face;
1488  dofs_on_face.resize(dofs_per_face);
1489 
1490  cell->face(f)->get_dof_indices(dofs_on_face,
1491  cell->active_fe_index());
1492  for (unsigned int i = 0; i < dofs_per_face; ++i)
1493  {
1494  const unsigned int global_idof_index = dofs_on_face[i];
1495  if (owned_dofs.is_element(global_idof_index))
1496  {
1497  boundary_dofs.insert(global_idof_index);
1498  }
1499  }
1500  }
1501  }
1502  return boundary_dofs.size();
1503  }
1504 
1505 
1506 
1507  template <int dim, int spacedim>
1508  std::size_t
1510  {
1511  std::size_t mem =
1519  MemoryConsumption::memory_consumption(vertex_dof_offsets));
1520  for (const auto &level : levels)
1523 
1524  return mem;
1525  }
1526 
1527 
1528 
1529  template <int dim, int spacedim>
1530  void
1532  const std::vector<unsigned int> &active_fe_indices)
1533  {
1534  Assert(active_fe_indices.size() == get_triangulation().n_active_cells(),
1535  ExcDimensionMismatch(active_fe_indices.size(),
1536  get_triangulation().n_active_cells()));
1537 
1538  create_active_fe_table();
1539  // we could set the values directly, since they are stored as
1540  // protected data of this object, but for simplicity we use the
1541  // cell-wise access. this way we also have to pass some debug-mode
1542  // tests which we would have to duplicate ourselves otherwise
1543  for (const auto &cell : active_cell_iterators())
1544  if (cell->is_locally_owned())
1545  cell->set_active_fe_index(active_fe_indices[cell->active_cell_index()]);
1546  }
1547 
1548 
1549 
1550  template <int dim, int spacedim>
1551  void
1553  std::vector<unsigned int> &active_fe_indices) const
1554  {
1555  active_fe_indices.resize(get_triangulation().n_active_cells());
1556 
1557  // we could try to extract the values directly, since they are
1558  // stored as protected data of this object, but for simplicity we
1559  // use the cell-wise access.
1560  for (const auto &cell : active_cell_iterators())
1561  if (!cell->is_artificial())
1562  active_fe_indices[cell->active_cell_index()] = cell->active_fe_index();
1563  }
1564 
1565 
1566 
1567  template <int dim, int spacedim>
1568  void
1572  {
1573  clear();
1574 
1575  if (this->tria != &tria)
1576  {
1577  for (auto &connection : tria_listeners)
1578  connection.disconnect();
1579  tria_listeners.clear();
1580 
1581  this->tria = &tria;
1582 
1583  setup_policy_and_listeners();
1584  }
1585 
1586  create_active_fe_table();
1587 
1588  distribute_dofs(fe);
1589  }
1590 
1591 
1592 
1593  template <int dim, int spacedim>
1594  void
1596  {
1597  Assert(
1598  tria != nullptr,
1599  ExcMessage(
1600  "You need to set the Triangulation in the DoFHandler using initialize() or "
1601  "in the constructor before you can distribute DoFs."));
1602  Assert(tria->n_levels() > 0,
1603  ExcMessage("The Triangulation you are using is empty!"));
1604  Assert(ff.size() > 0, ExcMessage("The hp::FECollection given is empty!"));
1605 
1606  // don't create a new object if the one we have is already appropriate
1607  if (fe_collection != ff)
1609 
1610  // ensure that the active_fe_indices vectors are initialized correctly
1611  create_active_fe_table();
1612 
1613  // make sure every processor knows the active_fe_indices
1614  // on both its own cells and all ghost cells
1617 
1618  // make sure that the fe collection is large enough to
1619  // cover all fe indices presently in use on the mesh
1620  for (const auto &cell : active_cell_iterators())
1621  if (!cell->is_artificial())
1622  Assert(cell->active_fe_index() < fe_collection.size(),
1623  ExcInvalidFEIndex(cell->active_fe_index(),
1624  fe_collection.size()));
1625  }
1626 
1627 
1628 
1629  template <int dim, int spacedim>
1630  void
1633  {
1634  // assign the fe_collection and initialize all active_fe_indices
1635  set_fe(ff);
1636 
1637  // If an underlying shared::Tria allows artificial cells,
1638  // then save the current set of subdomain ids, and set
1639  // subdomain ids to the "true" owner of each cell. we later
1640  // restore these flags
1641  std::vector<types::subdomain_id> saved_subdomain_ids;
1643  (dynamic_cast<const parallel::shared::Triangulation<dim, spacedim> *>(
1644  &get_triangulation()));
1645  if (shared_tria != nullptr && shared_tria->with_artificial_cells())
1646  {
1647  saved_subdomain_ids.resize(shared_tria->n_active_cells());
1648 
1649  const std::vector<types::subdomain_id> &true_subdomain_ids =
1650  shared_tria->get_true_subdomain_ids_of_cells();
1651 
1652  for (const auto &cell : shared_tria->active_cell_iterators())
1653  {
1654  const unsigned int index = cell->active_cell_index();
1655  saved_subdomain_ids[index] = cell->subdomain_id();
1656  cell->set_subdomain_id(true_subdomain_ids[index]);
1657  }
1658  }
1659 
1660  // then allocate space for all the other tables
1662  reserve_space(*this);
1663 
1664  // now undo the subdomain modification
1665  if (shared_tria != nullptr && shared_tria->with_artificial_cells())
1666  for (const auto &cell : shared_tria->active_cell_iterators())
1667  cell->set_subdomain_id(saved_subdomain_ids[cell->active_cell_index()]);
1668 
1669 
1670  // Clear user flags because we will need them. But first we save
1671  // them and make sure that we restore them later such that at the
1672  // end of this function the Triangulation will be in the same
1673  // state as it was at the beginning of this function.
1674  std::vector<bool> user_flags;
1675  tria->save_user_flags(user_flags);
1676  const_cast<Triangulation<dim, spacedim> &>(*tria).clear_user_flags();
1677 
1678 
1680 
1681  // Now for the real work:
1682  number_cache = policy->distribute_dofs();
1683 
1685 
1686  // do some housekeeping: compress indices
1687  {
1689  for (int level = levels.size() - 1; level >= 0; --level)
1690  tg += Threads::new_task(
1691  &::internal::hp::DoFLevel::compress_data<dim, spacedim>,
1692  *levels[level],
1693  fe_collection);
1694  tg.join_all();
1695  }
1696 
1697  // finally restore the user flags
1698  const_cast<Triangulation<dim, spacedim> &>(*tria).load_user_flags(
1699  user_flags);
1700  }
1701 
1702 
1703 
1704  template <int dim, int spacedim>
1705  void
1707  {
1708  // connect functions to signals of the underlying triangulation
1709  tria_listeners.push_back(this->tria->signals.pre_refinement.connect(
1710  [this]() { this->pre_refinement_action(); }));
1711  tria_listeners.push_back(this->tria->signals.post_refinement.connect(
1712  [this]() { this->post_refinement_action(); }));
1713  tria_listeners.push_back(this->tria->signals.create.connect(
1714  [this]() { this->post_refinement_action(); }));
1715 
1716  // decide whether we need a sequential or a parallel shared/distributed
1717  // policy and attach corresponding callback functions dealing with the
1718  // transfer of active_fe_indices
1720  *>(&this->get_triangulation()))
1721  {
1722  policy = std_cxx14::make_unique<
1724  DoFHandler<dim, spacedim>>>(*this);
1725 
1726  // repartitioning signals
1727  tria_listeners.push_back(
1728  this->tria->signals.pre_distributed_repartition.connect([this]() {
1729  internal::hp::DoFHandlerImplementation::Implementation::
1730  ensure_absence_of_future_fe_indices<dim, spacedim>(*this);
1731  }));
1732  tria_listeners.push_back(
1733  this->tria->signals.pre_distributed_repartition.connect(
1734  [this]() { this->pre_distributed_active_fe_index_transfer(); }));
1735  tria_listeners.push_back(
1736  this->tria->signals.post_distributed_repartition.connect(
1737  [this] { this->post_distributed_active_fe_index_transfer(); }));
1738 
1739  // refinement signals
1740  tria_listeners.push_back(
1741  this->tria->signals.pre_distributed_refinement.connect(
1742  [this]() { this->pre_distributed_active_fe_index_transfer(); }));
1743  tria_listeners.push_back(
1744  this->tria->signals.post_distributed_refinement.connect(
1745  [this]() { this->post_distributed_active_fe_index_transfer(); }));
1746 
1747  // serialization signals
1748  tria_listeners.push_back(
1749  this->tria->signals.post_distributed_save.connect([this]() {
1750  this->post_distributed_serialization_of_active_fe_indices();
1751  }));
1752  }
1753  else if (dynamic_cast<const parallel::shared::Triangulation<dim, spacedim>
1754  *>(&this->get_triangulation()) != nullptr)
1755  {
1756  policy =
1757  std_cxx14::make_unique<internal::DoFHandlerImplementation::Policy::
1758  ParallelShared<DoFHandler<dim, spacedim>>>(
1759  *this);
1760 
1761  // partitioning signals
1762  tria_listeners.push_back(
1763  this->tria->signals.pre_partition.connect([this]() {
1764  internal::hp::DoFHandlerImplementation::Implementation::
1765  ensure_absence_of_future_fe_indices(*this);
1766  }));
1767 
1768  // refinement signals
1769  tria_listeners.push_back(this->tria->signals.pre_refinement.connect(
1770  [this] { this->pre_active_fe_index_transfer(); }));
1771  tria_listeners.push_back(this->tria->signals.post_refinement.connect(
1772  [this] { this->post_active_fe_index_transfer(); }));
1773  }
1774  else
1775  {
1776  policy =
1777  std_cxx14::make_unique<internal::DoFHandlerImplementation::Policy::
1778  Sequential<DoFHandler<dim, spacedim>>>(
1779  *this);
1780 
1781  // refinement signals
1782  tria_listeners.push_back(this->tria->signals.pre_refinement.connect(
1783  [this] { this->pre_active_fe_index_transfer(); }));
1784  tria_listeners.push_back(this->tria->signals.post_refinement.connect(
1785  [this] { this->post_active_fe_index_transfer(); }));
1786  }
1787  }
1788 
1789 
1790 
1791  template <int dim, int spacedim>
1792  void
1794  {
1795  // release memory
1796  clear_space();
1797  }
1798 
1799 
1800 
1801  template <int dim, int spacedim>
1802  void
1804  const std::vector<types::global_dof_index> &new_numbers)
1805  {
1806  Assert(levels.size() > 0,
1807  ExcMessage(
1808  "You need to distribute DoFs before you can renumber them."));
1809 
1810  AssertDimension(new_numbers.size(), n_locally_owned_dofs());
1811 
1812 #ifdef DEBUG
1813  // assert that the new indices are
1814  // consecutively numbered if we are
1815  // working on a single
1816  // processor. this doesn't need to
1817  // hold in the case of a parallel
1818  // mesh since we map the interval
1819  // [0...n_dofs()) into itself but
1820  // only globally, not on each
1821  // processor
1822  if (n_locally_owned_dofs() == n_dofs())
1823  {
1824  std::vector<types::global_dof_index> tmp(new_numbers);
1825  std::sort(tmp.begin(), tmp.end());
1826  std::vector<types::global_dof_index>::const_iterator p = tmp.begin();
1828  for (; p != tmp.end(); ++p, ++i)
1829  Assert(*p == i, ExcNewNumbersNotConsecutive(i));
1830  }
1831  else
1832  for (const auto new_number : new_numbers)
1833  Assert(new_number < n_dofs(),
1834  ExcMessage(
1835  "New DoF index is not less than the total number of dofs."));
1836 #endif
1837 
1838  // uncompress the internal storage scheme of dofs on cells so that
1839  // we can access dofs in turns. uncompress in parallel, starting
1840  // with the most expensive levels (the highest ones)
1841  {
1843  for (int level = levels.size() - 1; level >= 0; --level)
1844  tg += Threads::new_task(
1845  &::internal::hp::DoFLevel::uncompress_data<dim, spacedim>,
1846  *levels[level],
1847  fe_collection);
1848  tg.join_all();
1849  }
1850 
1851  // do the renumbering
1852  number_cache = policy->renumber_dofs(new_numbers);
1853 
1854  // now re-compress the dof indices
1855  {
1857  for (int level = levels.size() - 1; level >= 0; --level)
1858  tg += Threads::new_task(
1859  &::internal::hp::DoFLevel::compress_data<dim, spacedim>,
1860  *levels[level],
1861  fe_collection);
1862  tg.join_all();
1863  }
1864  }
1865 
1866 
1867 
1868  template <int dim, int spacedim>
1869  unsigned int
1871  {
1872  Assert(fe_collection.size() > 0, ExcNoFESelected());
1873  return ::internal::hp::DoFHandlerImplementation::Implementation::
1874  max_couplings_between_dofs(*this);
1875  }
1876 
1877 
1878 
1879  template <int dim, int spacedim>
1880  unsigned int
1882  {
1883  Assert(fe_collection.size() > 0, ExcNoFESelected());
1884 
1885  switch (dim)
1886  {
1887  case 1:
1888  return fe_collection.max_dofs_per_vertex();
1889  case 2:
1890  return (3 * fe_collection.max_dofs_per_vertex() +
1891  2 * fe_collection.max_dofs_per_line());
1892  case 3:
1893  // we need to take refinement of one boundary face into
1894  // consideration here; in fact, this function returns what
1895  // #max_coupling_between_dofs<2> returns
1896  //
1897  // we assume here, that only four faces meet at the boundary;
1898  // this assumption is not justified and needs to be fixed some
1899  // time. fortunately, omitting it for now does no harm since
1900  // the matrix will cry foul if its requirements are not
1901  // satisfied
1902  return (19 * fe_collection.max_dofs_per_vertex() +
1903  28 * fe_collection.max_dofs_per_line() +
1904  8 * fe_collection.max_dofs_per_quad());
1905  default:
1906  Assert(false, ExcNotImplemented());
1907  return 0;
1908  }
1909  }
1910 
1911 
1912 
1913  template <int dim, int spacedim>
1914  void
1916  {
1917  // Create sufficiently many hp::DoFLevels.
1918  while (levels.size() < tria->n_levels())
1919  levels.emplace_back(new ::internal::hp::DoFLevel);
1920 
1921  // then make sure that on each level we have the appropriate size
1922  // of active_fe_indices; preset them to zero, i.e. the default FE
1923  for (unsigned int level = 0; level < levels.size(); ++level)
1924  {
1925  if (levels[level]->active_fe_indices.size() == 0 &&
1926  levels[level]->future_fe_indices.size() == 0)
1927  {
1928  levels[level]->active_fe_indices.resize(tria->n_raw_cells(level),
1929  0);
1930  levels[level]->future_fe_indices.resize(
1931  tria->n_raw_cells(level),
1933  }
1934  else
1935  {
1936  // Either the active_fe_indices have size zero because
1937  // they were just created, or the correct size. Other
1938  // sizes indicate that something went wrong.
1939  Assert(levels[level]->active_fe_indices.size() ==
1940  tria->n_raw_cells(level) &&
1941  levels[level]->future_fe_indices.size() ==
1942  tria->n_raw_cells(level),
1943  ExcInternalError());
1944  }
1945 
1946  // it may be that the previous table was compressed; in that
1947  // case, restore the correct active_fe_index. the fact that
1948  // this no longer matches the indices in the table is of no
1949  // importance because the current function is called at a
1950  // point where we have to recreate the dof_indices tables in
1951  // the levels anyway
1952  levels[level]->normalize_active_fe_indices();
1953  }
1954  }
1955 
1956 
1957 
1958  template <int dim, int spacedim>
1959  void
1961  {
1962  create_active_fe_table();
1963  }
1964 
1965 
1966 
1967  template <int dim, int spacedim>
1968  void
1970  {
1971  // Normally only one level is added, but if this Triangulation
1972  // is created by copy_triangulation, it can be more than one level.
1973  while (levels.size() < tria->n_levels())
1974  levels.emplace_back(new ::internal::hp::DoFLevel);
1975 
1976  // Coarsening can lead to the loss of levels. Hence remove them.
1977  while (levels.size() > tria->n_levels())
1978  {
1979  // drop the last element. that also releases the memory pointed to
1980  levels.pop_back();
1981  }
1982 
1983  Assert(levels.size() == tria->n_levels(), ExcInternalError());
1984  for (unsigned int i = 0; i < levels.size(); ++i)
1985  {
1986  // Resize active_fe_indices vectors. Use zero indicator to extend.
1987  levels[i]->active_fe_indices.resize(tria->n_raw_cells(i), 0);
1988 
1989  // Resize future_fe_indices vectors. Make sure that all
1990  // future_fe_indices have been cleared after refinement happened.
1991  //
1992  // We have used future_fe_indices to update all active_fe_indices
1993  // before refinement happened, thus we are safe to clear them now.
1994  levels[i]->future_fe_indices.assign(
1995  tria->n_raw_cells(i),
1997  }
1998  }
1999 
2000 
2001 
2002  template <int dim, int spacedim>
2003  void
2005  {
2006  // Finite elements need to be assigned to each cell by calling
2007  // distribute_dofs() first to make this functionality available.
2008  if (fe_collection.size() > 0)
2009  {
2010  Assert(active_fe_index_transfer == nullptr, ExcInternalError());
2011 
2012  active_fe_index_transfer =
2013  std_cxx14::make_unique<ActiveFEIndexTransfer>();
2014 
2017  }
2018  }
2019 
2020 
2021 
2022  template <int dim, int spacedim>
2023  void
2025  {
2026 #ifndef DEAL_II_WITH_P4EST
2027  Assert(false, ExcInternalError());
2028 #else
2029  // Finite elements need to be assigned to each cell by calling
2030  // distribute_dofs() first to make this functionality available.
2031  if (fe_collection.size() > 0)
2032  {
2033  Assert(active_fe_index_transfer == nullptr, ExcInternalError());
2034 
2035  active_fe_index_transfer =
2036  std_cxx14::make_unique<ActiveFEIndexTransfer>();
2037 
2038  // If we work on a p::d::Triangulation, we have to transfer all
2039  // active_fe_indices since ownership of cells may change. We will
2040  // use our p::d::CellDataTransfer member to achieve this. Further,
2041  // we prepare the values in such a way that they will correspond to
2042  // the active_fe_indices on the new mesh.
2043 
2044  // Gather all current future_fe_indices.
2045  active_fe_index_transfer->active_fe_indices.resize(
2046  get_triangulation().n_active_cells(), numbers::invalid_unsigned_int);
2047 
2048  for (const auto &cell : active_cell_iterators())
2049  if (cell->is_locally_owned())
2050  active_fe_index_transfer
2051  ->active_fe_indices[cell->active_cell_index()] =
2052  cell->future_fe_index();
2053 
2054  // Create transfer object and attach to it.
2055  const auto *distributed_tria = dynamic_cast<
2057  &this->get_triangulation());
2058 
2059  active_fe_index_transfer->cell_data_transfer = std_cxx14::make_unique<
2060  parallel::distributed::
2061  CellDataTransfer<dim, spacedim, std::vector<unsigned int>>>(
2062  *distributed_tria,
2063  /*transfer_variable_size_data=*/false,
2064  [this](const std::vector<unsigned int> &children_fe_indices) {
2065  return ::internal::hp::DoFHandlerImplementation::
2066  Implementation::determine_fe_from_children<dim, spacedim>(
2067  children_fe_indices, fe_collection);
2068  });
2069 
2070  active_fe_index_transfer->cell_data_transfer
2071  ->prepare_for_coarsening_and_refinement(
2072  active_fe_index_transfer->active_fe_indices);
2073  }
2074 #endif
2075  }
2076 
2077 
2078 
2079  template <int dim, int spacedim>
2080  void
2082  {
2083  // Finite elements need to be assigned to each cell by calling
2084  // distribute_dofs() first to make this functionality available.
2085  if (fe_collection.size() > 0)
2086  {
2087  Assert(active_fe_index_transfer != nullptr, ExcInternalError());
2088 
2091 
2092  // We have to distribute the information about active_fe_indices
2093  // of all cells (including the artificial ones) on all processors,
2094  // if a parallel::shared::Triangulation has been used.
2097 
2098  // Free memory.
2099  active_fe_index_transfer.reset();
2100  }
2101  }
2102 
2103 
2104 
2105  template <int dim, int spacedim>
2106  void
2108  {
2109 #ifndef DEAL_II_WITH_P4EST
2110  Assert(false, ExcInternalError());
2111 #else
2112  // Finite elements need to be assigned to each cell by calling
2113  // distribute_dofs() first to make this functionality available.
2114  if (fe_collection.size() > 0)
2115  {
2116  Assert(active_fe_index_transfer != nullptr, ExcInternalError());
2117 
2118  // Unpack active_fe_indices.
2119  active_fe_index_transfer->active_fe_indices.resize(
2120  get_triangulation().n_active_cells(), numbers::invalid_unsigned_int);
2121  active_fe_index_transfer->cell_data_transfer->unpack(
2122  active_fe_index_transfer->active_fe_indices);
2123 
2124  // Update all locally owned active_fe_indices.
2125  set_active_fe_indices(active_fe_index_transfer->active_fe_indices);
2126 
2127  // Update active_fe_indices on ghost cells.
2130 
2131  // Free memory.
2132  active_fe_index_transfer.reset();
2133  }
2134 #endif
2135  }
2136 
2137 
2138 
2139  template <int dim, int spacedim>
2140  void
2142  {
2143 #ifndef DEAL_II_WITH_P4EST
2144  Assert(false,
2145  ExcMessage(
2146  "You are attempting to use a functionality that is only available "
2147  "if deal.II was configured to use p4est, but cmake did not find a "
2148  "valid p4est library."));
2149 #else
2150  Assert(
2152  *>(&this->get_triangulation()) != nullptr),
2153  ExcMessage(
2154  "This functionality requires a parallel::distributed::Triangulation object."));
2155 
2156  // Finite elements need to be assigned to each cell by calling
2157  // distribute_dofs() first to make this functionality available.
2158  if (fe_collection.size() > 0)
2159  {
2160  Assert(active_fe_index_transfer == nullptr, ExcInternalError());
2161 
2162  active_fe_index_transfer =
2163  std_cxx14::make_unique<ActiveFEIndexTransfer>();
2164 
2165  // Create transfer object and attach to it.
2166  const auto *distributed_tria = dynamic_cast<
2168  &this->get_triangulation());
2169 
2170  active_fe_index_transfer->cell_data_transfer = std_cxx14::make_unique<
2171  parallel::distributed::
2172  CellDataTransfer<dim, spacedim, std::vector<unsigned int>>>(
2173  *distributed_tria,
2174  /*transfer_variable_size_data=*/false,
2175  [this](const std::vector<unsigned int> &children_fe_indices) {
2176  return ::internal::hp::DoFHandlerImplementation::
2177  Implementation::determine_fe_from_children<dim, spacedim>(
2178  children_fe_indices, fe_collection);
2179  });
2180 
2181  // If we work on a p::d::Triangulation, we have to transfer all
2182  // active fe indices since ownership of cells may change.
2183 
2184  // Gather all current active_fe_indices
2185  get_active_fe_indices(active_fe_index_transfer->active_fe_indices);
2186 
2187  // Attach to transfer object
2188  active_fe_index_transfer->cell_data_transfer->prepare_for_serialization(
2189  active_fe_index_transfer->active_fe_indices);
2190  }
2191 #endif
2192  }
2193 
2194 
2195  template <int dim, int spacedim>
2196  void
2197  DoFHandler<dim,
2198  spacedim>::post_distributed_serialization_of_active_fe_indices()
2199  {
2200 #ifndef DEAL_II_WITH_P4EST
2201  Assert(false,
2202  ExcMessage(
2203  "You are attempting to use a functionality that is only available "
2204  "if deal.II was configured to use p4est, but cmake did not find a "
2205  "valid p4est library."));
2206 #else
2207  if (fe_collection.size() > 0)
2208  {
2209  Assert(active_fe_index_transfer != nullptr, ExcInternalError());
2210 
2211  // Free memory.
2212  active_fe_index_transfer.reset();
2213  }
2214 #endif
2215  }
2216 
2217 
2218 
2219  template <int dim, int spacedim>
2220  void
2222  {
2223 #ifndef DEAL_II_WITH_P4EST
2224  Assert(false,
2225  ExcMessage(
2226  "You are attempting to use a functionality that is only available "
2227  "if deal.II was configured to use p4est, but cmake did not find a "
2228  "valid p4est library."));
2229 #else
2230  Assert(
2232  *>(&this->get_triangulation()) != nullptr),
2233  ExcMessage(
2234  "This functionality requires a parallel::distributed::Triangulation object."));
2235 
2236  // Finite elements need to be assigned to each cell by calling
2237  // distribute_dofs() first to make this functionality available.
2238  if (fe_collection.size() > 0)
2239  {
2240  Assert(active_fe_index_transfer == nullptr, ExcInternalError());
2241 
2242  active_fe_index_transfer =
2243  std_cxx14::make_unique<ActiveFEIndexTransfer>();
2244 
2245  // Create transfer object and attach to it.
2246  const auto *distributed_tria = dynamic_cast<
2248  &this->get_triangulation());
2249 
2250  active_fe_index_transfer->cell_data_transfer = std_cxx14::make_unique<
2251  parallel::distributed::
2252  CellDataTransfer<dim, spacedim, std::vector<unsigned int>>>(
2253  *distributed_tria,
2254  /*transfer_variable_size_data=*/false,
2255  [this](const std::vector<unsigned int> &children_fe_indices) {
2256  return ::internal::hp::DoFHandlerImplementation::
2257  Implementation::determine_fe_from_children<dim, spacedim>(
2258  children_fe_indices, fe_collection);
2259  });
2260 
2261  // Unpack active_fe_indices.
2262  active_fe_index_transfer->active_fe_indices.resize(
2263  get_triangulation().n_active_cells(), numbers::invalid_unsigned_int);
2264  active_fe_index_transfer->cell_data_transfer->deserialize(
2265  active_fe_index_transfer->active_fe_indices);
2266 
2267  // Update all locally owned active_fe_indices.
2268  set_active_fe_indices(active_fe_index_transfer->active_fe_indices);
2269 
2270  // Update active_fe_indices on ghost cells.
2273 
2274  // Free memory.
2275  active_fe_index_transfer.reset();
2276  }
2277 #endif
2278  }
2279 
2280 
2281 
2282  template <int dim, int spacedim>
2283  template <int structdim>
2285  DoFHandler<dim, spacedim>::get_dof_index(const unsigned int,
2286  const unsigned int,
2287  const unsigned int,
2288  const unsigned int) const
2289  {
2290  Assert(false, ExcNotImplemented());
2292  }
2293 
2294 
2295 
2296  template <int dim, int spacedim>
2297  template <int structdim>
2298  void
2299  DoFHandler<dim, spacedim>::set_dof_index(const unsigned int,
2300  const unsigned int,
2301  const unsigned int,
2302  const unsigned int,
2303  const types::global_dof_index) const
2304  {
2305  Assert(false, ExcNotImplemented());
2306  }
2307 
2308 
2309 
2310  template <int dim, int spacedim>
2311  void
2313  {
2314  levels.clear();
2315  faces.reset();
2316 
2317  vertex_dofs = std::vector<types::global_dof_index>();
2318  vertex_dof_offsets = std::vector<unsigned int>();
2319  }
2320 } // namespace hp
2321 
2322 
2323 
2324 /*-------------- Explicit Instantiations -------------------------------*/
2325 #include "dof_handler.inst"
2326 
2327 
2328 DEAL_II_NAMESPACE_CLOSE
unsigned int find_dominated_fe_extended(const std::set< unsigned int > &fes, const unsigned int codim=0) const
IteratorRange< active_cell_iterator > active_cell_iterators_on_level(const unsigned int level) const
unsigned int max_couplings_between_dofs() const
unsigned int n_active_cells() const
Definition: tria.cc:12673
const Triangulation< dim, spacedim > & get_triangulation() const
static const unsigned int invalid_unsigned_int
Definition: types.h:190
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
void clear_user_flags()
Definition: tria.cc:11172
cell_iterator begin(const unsigned int level=0) const
Definition: dof_handler.cc:922
unsigned int offset_type
Definition: dof_level.h:115
Task< RT > new_task(const std::function< RT()> &function)
static void distribute_fe_indices_on_refined_cells(::hp::DoFHandler< dim, spacedim > &dof_handler)
const std::vector< types::subdomain_id > & get_true_subdomain_ids_of_cells() const
Definition: shared_tria.cc:324
IteratorRange< active_cell_iterator > active_cell_iterators() const
Definition: tria.cc:12183
std::vector< std::unique_ptr<::internal::hp::DoFLevel > > levels
Definition: dof_handler.h:1226
void load_user_flags(std::istream &in)
Definition: tria.cc:11232
cell_iterator end() const
Definition: dof_handler.cc:951
::internal::DoFHandlerImplementation::NumberCache number_cache
Definition: dof_handler.h:1254
virtual void clear()
unsigned int max_dofs_per_face(const DoFHandler< dim, spacedim > &dh)
static unsigned int determine_fe_from_children(const std::vector< unsigned int > &children_fe_indices, ::hp::FECollection< dim, spacedim > &fe_collection)
unsigned int size() const
virtual std::size_t memory_consumption() const
std::unique_ptr< ActiveFEIndexTransfer > active_fe_index_transfer
Definition: dof_handler.h:1329
static void collect_fe_indices_on_cells_to_be_refined(::hp::DoFHandler< dim, spacedim > &dof_handler)
static void reserve_space_cells(DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_handler.cc:271
void initialize(const Triangulation< dim, spacedim > &tria, const FiniteElement< dim, spacedim > &fe)
Definition: dof_handler.cc:887
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
typename ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:312
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: dof_handler.cc:935
types::global_dof_index n_locally_owned_dofs() const
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:325
T sum(const T &t, const MPI_Comm &mpi_communicator)
#define Assert(cond, exc)
Definition: exceptions.h:1419
hp::FECollection< dim, spacedim > fe_collection
Definition: dof_handler.h:1237
IteratorRange< active_cell_iterator > active_cell_iterators() const
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::vector< types::global_dof_index > dofs
Definition: dof_faces.h:111
active_cell_iterator end_active(const unsigned int level) const
Definition: dof_handler.cc:971
types::global_dof_index n_boundary_dofs() const
static void reserve_space_release_space(DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_handler.cc:114
SmartPointer< const Triangulation< dim, spacedim >, DoFHandler< dim, spacedim > > tria
Definition: dof_handler.h:1230
types::global_dof_index n_dofs() const
std::unique_ptr<::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > > policy
Definition: dof_handler.h:1245
const hp::FECollection< dim, spacedim > & get_fe_collection() const
virtual void set_fe(const FiniteElement< dim, spacedim > &fe)
static unsigned int max_couplings_between_dofs(const DoFHandler< 1, spacedim > &dof_handler)
Definition: dof_handler.cc:853
std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces< dim > > faces
Definition: dof_handler.h:1396
std::vector< unsigned int > dof_offsets
Definition: dof_faces.h:105
Definition: hp.h:117
std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel< dim > > > levels
Definition: dof_handler.h:1384
internal::hp::DoFIndicesOnFacesOrEdges< 2 > quads
Definition: dof_faces.h:320
void renumber_dofs(const std::vector< types::global_dof_index > &new_numbers)
void swap(Vector< Number > &u, Vector< Number > &v)
Definition: vector.h:1390
unsigned int global_dof_index
Definition: types.h:77
hp::FECollection< dim, spacedim > fe_collection
Definition: dof_handler.h:1097
void clear_space()
static void ensure_absence_of_future_fe_indices(DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_handler.cc:94
const Triangulation< dim, spacedim > & get_triangulation() const
virtual void distribute_dofs(const FiniteElement< dim, spacedim > &fe)
static void reserve_space_vertices(DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_handler.cc:160
static const active_fe_index_type invalid_active_fe_index
Definition: dof_level.h:132
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
static ::ExceptionBase & ExcNotImplemented()
Iterator points to a valid object.
bool is_element(const size_type index) const
Definition: index_set.h:1766
static ::ExceptionBase & ExcNewNumbersNotConsecutive(types::global_dof_index arg1)
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:357
const types::boundary_id internal_face_boundary_id
Definition: types.h:247
unsigned int max_couplings_between_boundary_dofs() const
static void communicate_active_fe_indices(::hp::DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_handler.cc:976
static void reserve_space(DoFHandler< 1, spacedim > &dof_handler)
Definition: dof_handler.cc:684
typename ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:329
const types::global_dof_index invalid_dof_index
Definition: types.h:205
virtual ~DoFHandler() override
Definition: dof_handler.cc:870
void exchange_cell_data_to_ghosts(const MeshType &mesh, const std::function< std_cxx17::optional< DataType >(const typename MeshType::active_cell_iterator &)> &pack, const std::function< void(const typename MeshType::active_cell_iterator &, const DataType &)> &unpack)
const IndexSet & locally_owned_dofs() const
static ::ExceptionBase & ExcNoFESelected()
IteratorRange< cell_iterator > cell_iterators() const
static void reserve_space_faces(DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_handler.cc:365
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
std::vector< types::global_dof_index > vertex_dofs
Definition: dof_handler.h:1370
static ::ExceptionBase & ExcInternalError()