Reference documentation for deal.II version 9.1.0-pre
dof_handler.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/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/shared_tria.h>
22 #include <deal.II/distributed/tria.h>
23 
24 #include <deal.II/dofs/dof_accessor.h>
25 #include <deal.II/dofs/dof_handler_policy.h>
26 
27 #include <deal.II/fe/fe.h>
28 
29 #include <deal.II/grid/grid_tools.h>
30 #include <deal.II/grid/tria.h>
31 #include <deal.II/grid/tria_accessor.h>
32 #include <deal.II/grid/tria_iterator.h>
33 #include <deal.II/grid/tria_levels.h>
34 
35 #include <deal.II/hp/dof_faces.h>
36 #include <deal.II/hp/dof_handler.h>
37 #include <deal.II/hp/dof_level.h>
38 
39 #include <boost/serialization/array.hpp>
40 
41 #include <algorithm>
42 #include <functional>
43 #include <set>
44 
45 DEAL_II_NAMESPACE_OPEN
46 
47 // The following is necessary for compilation under Visual Studio which is
48 // unable to correctly distinguish between ::DoFHandler and
49 // ::hp::DoFHandler. Plus it makes code in dof_handler.cc easier to read.
50 #if defined(_MSC_VER) && (_MSC_VER >= 1800)
51 template <int dim, int spacedim>
52 using HpDoFHandler = ::hp::DoFHandler<dim, spacedim>;
53 #else
54 // When using older Visual Studio or a different compiler just fall back.
55 # define HpDoFHandler DoFHandler
56 #endif
57 
58 namespace parallel
59 {
60  namespace distributed
61  {
62  template <int, int>
63  class Triangulation;
64  }
65 } // namespace parallel
66 
67 
68 
69 namespace internal
70 {
71  namespace hp
72  {
73  namespace DoFHandlerImplementation
74  {
75  // access class ::hp::DoFHandler instead of namespace
76  // internal::hp::DoFHandler, etc
77  using ::hp::DoFHandler;
78 
84  {
89  template <int dim, int spacedim>
90  static void
92  {
93  // Release all space except the active_fe_indices field
94  // which we have to back up before
95  {
96  std::vector<std::vector<DoFLevel::active_fe_index_type>>
97  active_fe_backup(dof_handler.levels.size());
98  for (unsigned int level = 0; level < dof_handler.levels.size();
99  ++level)
100  active_fe_backup[level] =
101  std::move(dof_handler.levels[level]->active_fe_indices);
102 
103  // delete all levels and set them up newly, since vectors
104  // are troublesome if you want to change their size
105  dof_handler.clear_space();
106 
107  for (unsigned int level = 0; level < dof_handler.tria->n_levels();
108  ++level)
109  {
110  dof_handler.levels.emplace_back(new internal::hp::DoFLevel);
111  dof_handler.levels[level]->active_fe_indices =
112  std::move(active_fe_backup[level]);
113  }
114 
115  if (dim > 1)
116  dof_handler.faces =
117  std_cxx14::make_unique<internal::hp::DoFIndicesOnFaces<dim>>();
118  }
119  }
120 
121 
122 
127  template <int dim, int spacedim>
128  static void
130  {
131  // The final step in all of the reserve_space() functions is to set
132  // up vertex dof information. since vertices are sequentially
133  // numbered, what we do first is to set up an array in which
134  // we record whether a vertex is associated with any of the
135  // given fe's, by setting a bit. in a later step, we then
136  // actually allocate memory for the required dofs
137  //
138  // in the following, we only need to consider vertices that are
139  // adjacent to either a locally owned or a ghost cell; we never
140  // store anything on vertices that are only surrounded by
141  // artificial cells. so figure out that subset of vertices
142  // first
143  std::vector<bool> locally_used_vertices(
144  dof_handler.tria->n_vertices(), false);
145  for (typename HpDoFHandler<dim, spacedim>::active_cell_iterator cell =
146  dof_handler.begin_active();
147  cell != dof_handler.end();
148  ++cell)
149  if (!cell->is_artificial())
150  for (unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell;
151  ++v)
152  locally_used_vertices[cell->vertex_index(v)] = true;
153 
154  std::vector<std::vector<bool>> vertex_fe_association(
155  dof_handler.fe_collection.size(),
156  std::vector<bool>(dof_handler.tria->n_vertices(), false));
157 
158  for (typename HpDoFHandler<dim, spacedim>::active_cell_iterator cell =
159  dof_handler.begin_active();
160  cell != dof_handler.end();
161  ++cell)
162  if (!cell->is_artificial())
163  for (unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell;
164  ++v)
165  vertex_fe_association[cell->active_fe_index()]
166  [cell->vertex_index(v)] = true;
167 
168  // in debug mode, make sure that each vertex is associated
169  // with at least one fe (note that except for unused
170  // vertices, all vertices are actually active). this is of
171  // course only true for vertices that are part of either
172  // ghost or locally owned cells
173 #ifdef DEBUG
174  for (unsigned int v = 0; v < dof_handler.tria->n_vertices(); ++v)
175  if (locally_used_vertices[v] == true)
176  if (dof_handler.tria->vertex_used(v) == true)
177  {
178  unsigned int fe = 0;
179  for (; fe < dof_handler.fe_collection.size(); ++fe)
180  if (vertex_fe_association[fe][v] == true)
181  break;
182  Assert(fe != dof_handler.fe_collection.size(),
183  ExcInternalError());
184  }
185 #endif
186 
187  // next count how much memory we actually need. for each
188  // vertex, we need one slot per fe to store the fe_index,
189  // plus dofs_per_vertex for this fe. in addition, we need
190  // one slot as the end marker for the fe_indices. at the
191  // same time already fill the vertex_dof_offsets field
192  dof_handler.vertex_dof_offsets.resize(dof_handler.tria->n_vertices(),
194 
195  unsigned int vertex_slots_needed = 0;
196  for (unsigned int v = 0; v < dof_handler.tria->n_vertices(); ++v)
197  if (dof_handler.tria->vertex_used(v) == true)
198  if (locally_used_vertices[v] == true)
199  {
200  dof_handler.vertex_dof_offsets[v] = vertex_slots_needed;
201 
202  for (unsigned int fe = 0;
203  fe < dof_handler.fe_collection.size();
204  ++fe)
205  if (vertex_fe_association[fe][v] == true)
206  vertex_slots_needed +=
207  dof_handler.get_fe(fe).dofs_per_vertex + 1;
208 
209  // don't forget the end_marker:
210  ++vertex_slots_needed;
211  }
212 
213  // now allocate the space we have determined we need, and
214  // set up the linked lists for each of the vertices
215  dof_handler.vertex_dofs.resize(vertex_slots_needed,
217  for (unsigned int v = 0; v < dof_handler.tria->n_vertices(); ++v)
218  if (dof_handler.tria->vertex_used(v) == true)
219  if (locally_used_vertices[v] == true)
220  {
221  unsigned int current_index =
222  dof_handler.vertex_dof_offsets[v];
223  for (unsigned int fe = 0;
224  fe < dof_handler.fe_collection.size();
225  ++fe)
226  if (vertex_fe_association[fe][v] == true)
227  {
228  // if this vertex uses this fe, then set the
229  // fe_index and move the pointer ahead
230  dof_handler.vertex_dofs[current_index] = fe;
231  current_index +=
232  dof_handler.get_fe(fe).dofs_per_vertex + 1;
233  }
234  // finally place the end marker
235  dof_handler.vertex_dofs[current_index] =
237  }
238  }
239 
240 
245  template <int dim, int spacedim>
246  static void
248  {
249  // count how much space we need on each level for the cell
250  // dofs and set the dof_*_offsets data. initially set the
251  // latter to an invalid index, and only later set it to
252  // something reasonable for active dof_handler.cells
253  //
254  // note that for dof_handler.cells, the situation is simpler
255  // than for other (lower dimensional) objects since exactly
256  // one finite element is used for it
257  for (unsigned int level = 0; level < dof_handler.tria->n_levels();
258  ++level)
259  {
260  dof_handler.levels[level]->dof_offsets =
261  std::vector<DoFLevel::offset_type>(
262  dof_handler.tria->n_raw_cells(level),
263  (DoFLevel::offset_type)(-1));
264  dof_handler.levels[level]->cell_cache_offsets =
265  std::vector<DoFLevel::offset_type>(
266  dof_handler.tria->n_raw_cells(level),
267  (DoFLevel::offset_type)(-1));
268 
269  types::global_dof_index next_free_dof = 0;
270  types::global_dof_index cache_size = 0;
271  typename HpDoFHandler<dim, spacedim>::active_cell_iterator
272  cell = dof_handler.begin_active(level),
273  endc = dof_handler.end_active(level);
274  for (; cell != endc; ++cell)
275  if (!cell->has_children() && !cell->is_artificial())
276  {
277  dof_handler.levels[level]->dof_offsets[cell->index()] =
278  next_free_dof;
279  next_free_dof +=
280  cell->get_fe().template n_dofs_per_object<dim>();
281 
282  dof_handler.levels[level]
283  ->cell_cache_offsets[cell->index()] = cache_size;
284  cache_size += cell->get_fe().dofs_per_cell;
285  }
286 
287  dof_handler.levels[level]->dof_indices =
288  std::vector<types::global_dof_index>(
289  next_free_dof, numbers::invalid_dof_index);
290  dof_handler.levels[level]->cell_dof_indices_cache =
291  std::vector<types::global_dof_index>(
292  cache_size, numbers::invalid_dof_index);
293  }
294 
295  // safety check: make sure that the number of DoFs we
296  // allocated is actually correct (above we have also set the
297  // dof_*_offsets field, so we couldn't use this simpler
298  // algorithm)
299 #ifdef DEBUG
300  for (unsigned int level = 0; level < dof_handler.tria->n_levels();
301  ++level)
302  {
303  types::global_dof_index counter = 0;
304  typename HpDoFHandler<dim, spacedim>::active_cell_iterator
305  cell = dof_handler.begin_active(level),
306  endc = dof_handler.end_active(level);
307  for (; cell != endc; ++cell)
308  if (!cell->has_children() && !cell->is_artificial())
309  counter += cell->get_fe().template n_dofs_per_object<dim>();
310 
311  Assert(dof_handler.levels[level]->dof_indices.size() == counter,
312  ExcInternalError());
313 
314  // also check that the number of unassigned slots in the
315  // dof_offsets equals the number of cells on that level minus the
316  // number of active, non-artificial cells (because these are
317  // exactly the cells on which we do something)
318  unsigned int n_active_non_artificial_cells = 0;
319  for (cell = dof_handler.begin_active(level); cell != endc; ++cell)
320  if (!cell->has_children() && !cell->is_artificial())
321  ++n_active_non_artificial_cells;
322 
323  Assert(static_cast<unsigned int>(std::count(
324  dof_handler.levels[level]->dof_offsets.begin(),
325  dof_handler.levels[level]->dof_offsets.end(),
326  (DoFLevel::offset_type)(-1))) ==
327  dof_handler.tria->n_raw_cells(level) -
328  n_active_non_artificial_cells,
329  ExcInternalError());
330  }
331 #endif
332  }
333 
334 
335 
340  template <int dim, int spacedim>
341  static void
343  {
344  // make the code generic between lines and quads
345  std::vector<unsigned int> &face_dof_offsets =
346  (dim == 2 ?
347  reinterpret_cast<::internal::hp::DoFIndicesOnFaces<2> &>(
348  *dof_handler.faces)
349  .lines.dof_offsets :
350  reinterpret_cast<::internal::hp::DoFIndicesOnFaces<3> &>(
351  *dof_handler.faces)
352  .quads.dof_offsets);
353 
354  std::vector<types::global_dof_index> &face_dof_indices =
355  (dim == 2 ?
356  reinterpret_cast<::internal::hp::DoFIndicesOnFaces<2> &>(
357  *dof_handler.faces)
358  .lines.dofs :
359  reinterpret_cast<::internal::hp::DoFIndicesOnFaces<3> &>(
360  *dof_handler.faces)
361  .quads.dofs);
362 
363  // FACE DOFS
364  //
365  // count face dofs, then allocate as much space
366  // as we need and prime the linked list for faces (see the
367  // description in hp::DoFLevel) with the indices we will
368  // need. note that our task is more complicated than for the
369  // cell case above since two adjacent cells may have different
370  // active_fe_indices, in which case we need to allocate
371  // *two* sets of face dofs for the same face
372  //
373  // the way we do things is that we loop over all active
374  // cells (these are the ones that have DoFs only
375  // anyway) and all their faces. We note in the
376  // user flags whether we have previously visited a face and
377  // if so skip it (consequently, we have to save and later
378  // restore the face flags)
379  {
380  std::vector<bool> saved_face_user_flags;
381  switch (dim)
382  {
383  case 2:
384  {
385  const_cast<::Triangulation<dim, spacedim> &>(
386  *dof_handler.tria)
387  .save_user_flags_line(saved_face_user_flags);
388  const_cast<::Triangulation<dim, spacedim> &>(
389  *dof_handler.tria)
390  .clear_user_flags_line();
391 
392  break;
393  }
394 
395  case 3:
396  {
397  const_cast<::Triangulation<dim, spacedim> &>(
398  *dof_handler.tria)
399  .save_user_flags_quad(saved_face_user_flags);
400  const_cast<::Triangulation<dim, spacedim> &>(
401  *dof_handler.tria)
402  .clear_user_flags_quad();
403 
404  break;
405  }
406 
407  default:
408  Assert(false, ExcNotImplemented());
409  }
410 
411  // an array to hold how many slots (see the hp::DoFLevel
412  // class) we will have to store on each level
413  unsigned int n_face_slots = 0;
414 
415  for (typename HpDoFHandler<dim, spacedim>::active_cell_iterator
416  cell = dof_handler.begin_active();
417  cell != dof_handler.end();
418  ++cell)
419  if (!cell->is_artificial())
420  for (unsigned int face = 0;
421  face < GeometryInfo<dim>::faces_per_cell;
422  ++face)
423  if (!cell->face(face)->user_flag_set())
424  {
425  // ok, face has not been visited. so we need to
426  // allocate space for it. let's see how much we
427  // need: we need one set if a) there is no
428  // neighbor behind this face, or b) the neighbor
429  // is either coarser or finer than we are, or c)
430  // the neighbor is artificial, or d) the neighbor
431  // is neither coarser nor finer, nor is artificial,
432  // and just so happens to have the same active_fe_index :
433  if (cell->at_boundary(face) ||
434  cell->face(face)->has_children() ||
435  cell->neighbor_is_coarser(face) ||
436  (!cell->at_boundary(face) &&
437  cell->neighbor(face)->is_artificial()) ||
438  (!cell->at_boundary(face) &&
439  !cell->neighbor(face)->is_artificial() &&
440  (cell->active_fe_index() ==
441  cell->neighbor(face)->active_fe_index())))
442  // ok, one set of dofs. that makes one active_fe_index,
443  // 1 times dofs_per_face dofs, and one stop index
444  n_face_slots +=
445  1 +
446  dof_handler.get_fe(cell->active_fe_index())
447  .template n_dofs_per_object<dim - 1>() +
448  1;
449 
450  // otherwise we do indeed need two sets, i.e. two
451  // active_fe_indices, two sets of dofs, and one stop
452  // index:
453  else
454  n_face_slots +=
455  (2 +
456  dof_handler.get_fe(cell->active_fe_index())
457  .template n_dofs_per_object<dim - 1>() +
458  dof_handler
459  .get_fe(cell->neighbor(face)->active_fe_index())
460  .template n_dofs_per_object<dim - 1>() +
461  1);
462 
463  // mark this face as visited
464  cell->face(face)->set_user_flag();
465  }
466 
467  // now that we know how many face dofs we will have to
468  // have on each level, allocate the memory. note that we
469  // allocate offsets for all faces, though only the active
470  // ones will have a non-invalid value later on
471  face_dof_offsets =
472  std::vector<unsigned int>(dof_handler.tria->n_raw_faces(),
473  (unsigned int)(-1));
474  face_dof_indices =
475  std::vector<types::global_dof_index>(n_face_slots,
477 
478  // with the memory now allocated, loop over the
479  // dof_handler.cells again and prime the _offset values as
480  // well as the fe_index fields
481  switch (dim)
482  {
483  case 2:
484  {
485  const_cast<::Triangulation<dim, spacedim> &>(
486  *dof_handler.tria)
487  .clear_user_flags_line();
488 
489  break;
490  }
491 
492  case 3:
493  {
494  const_cast<::Triangulation<dim, spacedim> &>(
495  *dof_handler.tria)
496  .clear_user_flags_quad();
497 
498  break;
499  }
500 
501  default:
502  Assert(false, ExcNotImplemented());
503  }
504 
505  unsigned int next_free_face_slot = 0;
506 
507  for (typename HpDoFHandler<dim, spacedim>::active_cell_iterator
508  cell = dof_handler.begin_active();
509  cell != dof_handler.end();
510  ++cell)
511  if (!cell->is_artificial())
512  for (unsigned int face = 0;
513  face < GeometryInfo<dim>::faces_per_cell;
514  ++face)
515  if (!cell->face(face)->user_flag_set())
516  {
517  // same decision tree as before
518  if (cell->at_boundary(face) ||
519  cell->face(face)->has_children() ||
520  cell->neighbor_is_coarser(face) ||
521  (!cell->at_boundary(face) &&
522  cell->neighbor(face)->is_artificial()) ||
523  (!cell->at_boundary(face) &&
524  !cell->neighbor(face)->is_artificial() &&
525  (cell->active_fe_index() ==
526  cell->neighbor(face)->active_fe_index())))
527  {
528  face_dof_offsets[cell->face(face)->index()] =
529  next_free_face_slot;
530 
531  // set first slot for this face to
532  // active_fe_index of this face
533  face_dof_indices[next_free_face_slot] =
534  cell->active_fe_index();
535 
536  // the next dofs_per_face indices remain unset
537  // for the moment (i.e. at invalid_dof_index).
538  // following this comes the stop index, which
539  // also is invalid_dof_index and therefore
540  // does not have to be explicitly set
541 
542  // finally, mark those slots as used
543  next_free_face_slot +=
544  dof_handler.get_fe(cell->active_fe_index())
545  .template n_dofs_per_object<dim - 1>() +
546  2;
547  }
548  else
549  {
550  face_dof_offsets[cell->face(face)->index()] =
551  next_free_face_slot;
552 
553  // set first 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_line/quad indices remain unset
559  // for the moment (i.e. at invalid_dof_index).
560  //
561  // then comes the fe_index for the neighboring
562  // cell:
563  face_dof_indices
564  [next_free_face_slot +
565  dof_handler.get_fe(cell->active_fe_index())
566  .template n_dofs_per_object<dim - 1>() +
567  1] = cell->neighbor(face)->active_fe_index();
568  // then again a set of dofs that we need not
569  // set right now
570  //
571  // following this comes the stop index, which
572  // also is invalid_dof_index and therefore
573  // does not have to be explicitly set
574 
575  // finally, mark those slots as used
576  next_free_face_slot +=
577  (dof_handler.get_fe(cell->active_fe_index())
578  .template n_dofs_per_object<dim - 1>() +
579  dof_handler
580  .get_fe(cell->neighbor(face)->active_fe_index())
581  .template n_dofs_per_object<dim - 1>() +
582  3);
583  }
584 
585  // mark this face as visited
586  cell->face(face)->set_user_flag();
587  }
588 
589  // we should have moved the cursor for each level to the
590  // total number of dofs on that level. check that
591  Assert(next_free_face_slot == n_face_slots, ExcInternalError());
592 
593  // at the end, restore the user flags for the faces
594  switch (dim)
595  {
596  case 2:
597  {
598  const_cast<::Triangulation<dim, spacedim> &>(
599  *dof_handler.tria)
600  .load_user_flags_line(saved_face_user_flags);
601 
602  break;
603  }
604 
605  case 3:
606  {
607  const_cast<::Triangulation<dim, spacedim> &>(
608  *dof_handler.tria)
609  .load_user_flags_quad(saved_face_user_flags);
610 
611  break;
612  }
613 
614  default:
615  Assert(false, ExcNotImplemented());
616  }
617  }
618  }
619 
620 
627  template <int spacedim>
628  static void reserve_space(DoFHandler<1, spacedim> &dof_handler)
629  {
630  Assert(dof_handler.fe_collection.size() > 0,
632  Assert(dof_handler.tria->n_levels() > 0,
633  ExcMessage("The current Triangulation must not be empty."));
634  Assert(dof_handler.tria->n_levels() == dof_handler.levels.size(),
635  ExcInternalError());
636 
637  reserve_space_release_space(dof_handler);
638 
639  reserve_space_cells(dof_handler);
640  reserve_space_vertices(dof_handler);
641  }
642 
643 
644 
645  template <int spacedim>
646  static void reserve_space(DoFHandler<2, spacedim> &dof_handler)
647  {
648  Assert(dof_handler.fe_collection.size() > 0,
650  Assert(dof_handler.tria->n_levels() > 0,
651  ExcMessage("The current Triangulation must not be empty."));
652  Assert(dof_handler.tria->n_levels() == dof_handler.levels.size(),
653  ExcInternalError());
654 
655  reserve_space_release_space(dof_handler);
656 
657  // FIRST CELLS AND FACES
658  reserve_space_cells(dof_handler);
659  reserve_space_faces(dof_handler);
660 
661  // VERTEX DOFS
662  reserve_space_vertices(dof_handler);
663  }
664 
665 
666 
667  template <int spacedim>
668  static void reserve_space(DoFHandler<3, spacedim> &dof_handler)
669  {
670  const unsigned int dim = 3;
671 
672  Assert(dof_handler.fe_collection.size() > 0,
674  Assert(dof_handler.tria->n_levels() > 0,
675  ExcMessage("The current Triangulation must not be empty."));
676  Assert(dof_handler.tria->n_levels() == dof_handler.levels.size(),
677  ExcInternalError());
678 
679  reserve_space_release_space(dof_handler);
680 
681  // FIRST CELLS AND FACES
682  reserve_space_cells(dof_handler);
683  reserve_space_faces(dof_handler);
684 
685 
686  // LINE DOFS
687 
688  // the situation here is pretty much like with vertices:
689  // there can be an arbitrary number of finite elements
690  // associated with each line.
691  //
692  // the algorithm we use is somewhat similar to what we do in
693  // reserve_space_vertices()
694  if (true)
695  {
696  // what we do first is to set up an array in which we
697  // record whether a line is associated with any of the
698  // given fe's, by setting a bit. in a later step, we
699  // then actually allocate memory for the required dofs
700  std::vector<std::vector<bool>> line_fe_association(
701  dof_handler.fe_collection.size(),
702  std::vector<bool>(dof_handler.tria->n_raw_lines(), false));
703 
704  for (typename HpDoFHandler<dim, spacedim>::active_cell_iterator
705  cell = dof_handler.begin_active();
706  cell != dof_handler.end();
707  ++cell)
708  if (!cell->is_artificial())
709  for (unsigned int l = 0;
710  l < GeometryInfo<dim>::lines_per_cell;
711  ++l)
712  line_fe_association[cell->active_fe_index()]
713  [cell->line_index(l)] = true;
714 
715  // first check which of the lines is used at all,
716  // i.e. is associated with a finite element. we do this
717  // since not all lines may actually be used, in which
718  // case we do not have to allocate any memory at all
719  std::vector<bool> line_is_used(dof_handler.tria->n_raw_lines(),
720  false);
721  for (unsigned int line = 0;
722  line < dof_handler.tria->n_raw_lines();
723  ++line)
724  for (unsigned int fe = 0; fe < dof_handler.fe_collection.size();
725  ++fe)
726  if (line_fe_association[fe][line] == true)
727  {
728  line_is_used[line] = true;
729  break;
730  }
731 
732  // next count how much memory we actually need. for each
733  // line, we need one slot per fe to store the fe_index,
734  // plus dofs_per_line for this fe. in addition, we need
735  // one slot as the end marker for the fe_indices. at the
736  // same time already fill the line_dofs_offsets field
737  dof_handler.faces->lines.dof_offsets.resize(
738  dof_handler.tria->n_raw_lines(), numbers::invalid_unsigned_int);
739 
740  unsigned int line_slots_needed = 0;
741  for (unsigned int line = 0;
742  line < dof_handler.tria->n_raw_lines();
743  ++line)
744  if (line_is_used[line] == true)
745  {
746  dof_handler.faces->lines.dof_offsets[line] =
747  line_slots_needed;
748 
749  for (unsigned int fe = 0;
750  fe < dof_handler.fe_collection.size();
751  ++fe)
752  if (line_fe_association[fe][line] == true)
753  line_slots_needed +=
754  dof_handler.get_fe(fe).dofs_per_line + 1;
755  ++line_slots_needed;
756  }
757 
758  // now allocate the space we have determined we need,
759  // and set up the linked lists for each of the lines
760  dof_handler.faces->lines.dofs.resize(line_slots_needed,
762  for (unsigned int line = 0;
763  line < dof_handler.tria->n_raw_lines();
764  ++line)
765  if (line_is_used[line] == true)
766  {
767  unsigned int pointer =
768  dof_handler.faces->lines.dof_offsets[line];
769  for (unsigned int fe = 0;
770  fe < dof_handler.fe_collection.size();
771  ++fe)
772  if (line_fe_association[fe][line] == true)
773  {
774  // if this line uses this fe, then set the
775  // fe_index and move the pointer ahead
776  dof_handler.faces->lines.dofs[pointer] = fe;
777  pointer += dof_handler.get_fe(fe).dofs_per_line + 1;
778  }
779  // finally place the end marker
780  dof_handler.faces->lines.dofs[pointer] =
782  }
783  }
784 
785 
786  // VERTEX DOFS
787  reserve_space_vertices(dof_handler);
788  }
789 
790 
794  template <int spacedim>
795  static unsigned int
797  {
798  return std::min(static_cast<types::global_dof_index>(
799  3 *
800  dof_handler.fe_collection.max_dofs_per_vertex() +
801  2 * dof_handler.fe_collection.max_dofs_per_line()),
802  dof_handler.n_dofs());
803  }
804 
805 
806 
807  template <int spacedim>
808  static unsigned int
809  max_couplings_between_dofs(const DoFHandler<2, spacedim> &dof_handler)
810  {
811  // get these numbers by drawing pictures and counting...
812  // example:
813  // | | |
814  // --x-----x--x--X--
815  // | | | |
816  // | x--x--x
817  // | | | |
818  // --x--x--*--x--x--
819  // | | | |
820  // x--x--x |
821  // | | | |
822  // --X--x--x-----x--
823  // | | |
824  // x = vertices connected with center vertex *;
825  // = total of 19
826  //
827  // (the X vertices are connected with * if the vertices
828  // adjacent to X are hanging nodes) count lines -> 28 (don't
829  // forget to count mother and children separately!)
830  types::global_dof_index max_couplings;
831  switch (dof_handler.tria->max_adjacent_cells())
832  {
833  case 4:
834  max_couplings =
835  19 * dof_handler.fe_collection.max_dofs_per_vertex() +
836  28 * dof_handler.fe_collection.max_dofs_per_line() +
837  8 * dof_handler.fe_collection.max_dofs_per_quad();
838  break;
839  case 5:
840  max_couplings =
841  21 * dof_handler.fe_collection.max_dofs_per_vertex() +
842  31 * dof_handler.fe_collection.max_dofs_per_line() +
843  9 * dof_handler.fe_collection.max_dofs_per_quad();
844  break;
845  case 6:
846  max_couplings =
847  28 * dof_handler.fe_collection.max_dofs_per_vertex() +
848  42 * dof_handler.fe_collection.max_dofs_per_line() +
849  12 * dof_handler.fe_collection.max_dofs_per_quad();
850  break;
851  case 7:
852  max_couplings =
853  30 * dof_handler.fe_collection.max_dofs_per_vertex() +
854  45 * dof_handler.fe_collection.max_dofs_per_line() +
855  13 * dof_handler.fe_collection.max_dofs_per_quad();
856  break;
857  case 8:
858  max_couplings =
859  37 * dof_handler.fe_collection.max_dofs_per_vertex() +
860  56 * dof_handler.fe_collection.max_dofs_per_line() +
861  16 * dof_handler.fe_collection.max_dofs_per_quad();
862  break;
863  default:
864  Assert(false, ExcNotImplemented());
865  max_couplings = 0;
866  };
867  return std::min(max_couplings, dof_handler.n_dofs());
868  }
869 
870 
871  template <int spacedim>
872  static unsigned int
873  max_couplings_between_dofs(const DoFHandler<3, spacedim> &dof_handler)
874  {
875  // TODO:[?] Invent significantly better estimates than the ones in
876  // this function
877  // doing the same thing here is a rather complicated thing,
878  // compared to the 2d case, since it is hard to draw
879  // pictures with several refined hexahedra :-) so I
880  // presently only give a coarse estimate for the case that
881  // at most 8 hexes meet at each vertex
882  //
883  // can anyone give better estimate here?
884  const unsigned int max_adjacent_cells =
885  dof_handler.tria->max_adjacent_cells();
886 
887  types::global_dof_index max_couplings;
888  if (max_adjacent_cells <= 8)
889  max_couplings =
890  7 * 7 * 7 * dof_handler.fe_collection.max_dofs_per_vertex() +
891  7 * 6 * 7 * 3 * dof_handler.fe_collection.max_dofs_per_line() +
892  9 * 4 * 7 * 3 * dof_handler.fe_collection.max_dofs_per_quad() +
893  27 * dof_handler.fe_collection.max_dofs_per_hex();
894  else
895  {
896  Assert(false, ExcNotImplemented());
897  max_couplings = 0;
898  }
899 
900  return std::min(max_couplings, dof_handler.n_dofs());
901  }
902 
903 
909  template <int dim, int spacedim>
910  static void
912  ::hp::DoFHandler<dim, spacedim> &dof_handler)
913  {
915  dynamic_cast<
917  &dof_handler.get_triangulation()))
918  {
919  // we have a shared triangulation. in this case, every processor
920  // knows about all cells, but every processor only has knowledge
921  // about the active_fe_index on the cells it owns.
922  //
923  // we can create a complete set of active_fe_indices by letting
924  // every processor create a vector of indices for all cells,
925  // filling only those on the cells it owns and setting the indices
926  // on the other cells to zero. then we add all of these vectors
927  // up, and because every vector entry has exactly one processor
928  // that owns it, the sum is correct
929  std::vector<unsigned int> active_fe_indices(tr->n_active_cells(),
930  (unsigned int)0);
931  for (const auto &cell : dof_handler.active_cell_iterators())
932  if (cell->is_locally_owned())
933  active_fe_indices[cell->active_cell_index()] =
934  cell->active_fe_index();
935 
936  Utilities::MPI::sum(active_fe_indices,
937  tr->get_communicator(),
938  active_fe_indices);
939 
940  // now go back and fill the active_fe_index on ghost
941  // cells. we would like to call cell->set_active_fe_index(),
942  // but that function does not allow setting these indices on
943  // non-locally_owned cells. so we have to work around the
944  // issue a little bit by accessing the underlying data
945  // structures directly
946  for (auto cell : dof_handler.active_cell_iterators())
947  if (cell->is_ghost())
948  dof_handler.levels[cell->level()]->set_active_fe_index(
949  cell->index(),
950  active_fe_indices[cell->active_cell_index()]);
951  }
953  *tr = dynamic_cast<
955  *>(&dof_handler.get_triangulation()))
956  {
957  // For completely distributed meshes, use the function that is
958  // able to move data from locally owned cells on one processor to
959  // the corresponding ghost cells on others. To this end, we need
960  // to have functions that can pack and unpack the data we want to
961  // transport -- namely, the single unsigned int active_fe_index
962  // objects
963  auto pack =
964  [](const typename ::hp::DoFHandler<dim, spacedim>::
965  active_cell_iterator &cell) -> unsigned int {
966  return cell->active_fe_index();
967  };
968 
969  auto unpack =
970  [&dof_handler](
971  const typename ::hp::DoFHandler<dim, spacedim>::
972  active_cell_iterator &cell,
973  const unsigned int & active_fe_index) -> void {
974  // we would like to say
975  // cell->set_active_fe_index(active_fe_index);
976  // but this is not allowed on cells that are not
977  // locally owned, and we are on a ghost cell
978  dof_handler.levels[cell->level()]->set_active_fe_index(
979  cell->index(), active_fe_index);
980  };
981 
983  unsigned int,
984  ::hp::DoFHandler<dim, spacedim>>(dof_handler,
985  pack,
986  unpack);
987  }
988  else
989  {
990  // a sequential triangulation. there is nothing we need to do here
991  Assert(
992  (dynamic_cast<const parallel::Triangulation<dim, spacedim> *>(
993  &dof_handler.get_triangulation()) == nullptr),
994  ExcInternalError());
995  }
996  }
997  };
998  } // namespace DoFHandlerImplementation
999  } // namespace hp
1000 } // namespace internal
1001 
1002 
1003 namespace hp
1004 {
1005  template <int dim, int spacedim>
1006  const unsigned int DoFHandler<dim, spacedim>::dimension;
1007 
1008  template <int dim, int spacedim>
1010 
1011  template <int dim, int spacedim>
1013 
1014 
1015 
1016  template <int dim, int spacedim>
1018  : tria(nullptr, typeid(*this).name())
1019  , faces(nullptr)
1020  {}
1021 
1022 
1023  template <int dim, int spacedim>
1026  : tria(&tria, typeid(*this).name())
1027  , faces(nullptr)
1028  {
1030 
1032  }
1033 
1034 
1035  template <int dim, int spacedim>
1037  {
1038  // unsubscribe as a listener to refinement of the underlying
1039  // triangulation
1040  for (unsigned int i = 0; i < tria_listeners.size(); ++i)
1041  tria_listeners[i].disconnect();
1042  tria_listeners.clear();
1043 
1044  // ...and release allocated memory
1045  // virtual functions called in constructors and destructors never use the
1046  // override in a derived class
1047  // for clarity be explicit on which function is called
1049  }
1050 
1051 
1052  /*------------------------ Cell iterator functions ------------------------*/
1053 
1054  template <int dim, int spacedim>
1056  DoFHandler<dim, spacedim>::begin(const unsigned int level) const
1057  {
1058  return cell_iterator(*this->get_triangulation().begin(level), this);
1059  }
1060 
1061 
1062 
1063  template <int dim, int spacedim>
1065  DoFHandler<dim, spacedim>::begin_active(const unsigned int level) const
1066  {
1067  // level is checked in begin
1068  cell_iterator i = begin(level);
1069  if (i.state() != IteratorState::valid)
1070  return i;
1071  while (i->has_children())
1072  if ((++i).state() != IteratorState::valid)
1073  return i;
1074  return i;
1075  }
1076 
1077 
1078 
1079  template <int dim, int spacedim>
1082  {
1083  return cell_iterator(&this->get_triangulation(), -1, -1, this);
1084  }
1085 
1086 
1087  template <int dim, int spacedim>
1089  DoFHandler<dim, spacedim>::end(const unsigned int level) const
1090  {
1091  return (level == this->get_triangulation().n_levels() - 1 ?
1092  end() :
1093  begin(level + 1));
1094  }
1095 
1096 
1097  template <int dim, int spacedim>
1099  DoFHandler<dim, spacedim>::end_active(const unsigned int level) const
1100  {
1101  return (level == this->get_triangulation().n_levels() - 1 ?
1103  begin_active(level + 1));
1104  }
1105 
1106 
1107 
1108  template <int dim, int spacedim>
1111  {
1113  begin(), end());
1114  }
1115 
1116 
1117  template <int dim, int spacedim>
1120  {
1121  return IteratorRange<
1123  end());
1124  }
1125 
1126 
1127 
1128  template <int dim, int spacedim>
1131  const unsigned int level) const
1132  {
1134  begin(level), end(level));
1135  }
1136 
1137 
1138 
1139  template <int dim, int spacedim>
1142  const unsigned int level) const
1143  {
1144  return IteratorRange<
1146  begin_active(level), end_active(level));
1147  }
1148 
1149 
1150 
1151  //------------------------------------------------------------------
1152 
1153 
1154  template <int dim, int spacedim>
1157  {
1158  Assert(fe_collection.size() > 0, ExcNoFESelected());
1159 
1160  std::set<types::global_dof_index> boundary_dofs;
1161  std::vector<types::global_dof_index> dofs_on_face;
1162  dofs_on_face.reserve(this->get_fe_collection().max_dofs_per_face());
1163 
1164  // loop over all faces to check whether they are at a
1165  // boundary. note that we need not take special care of single
1166  // lines in 3d (using @p{cell->has_boundary_lines}), since we do
1167  // not support boundaries of dimension dim-2, and so every
1168  // boundary line is also part of a boundary face.
1169  typename HpDoFHandler<dim, spacedim>::active_cell_iterator
1170  cell = this->begin_active(),
1171  endc = this->end();
1172  for (; cell != endc; ++cell)
1173  for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
1174  if (cell->at_boundary(f))
1175  {
1176  const unsigned int dofs_per_face = cell->get_fe().dofs_per_face;
1177  dofs_on_face.resize(dofs_per_face);
1178 
1179  cell->face(f)->get_dof_indices(dofs_on_face,
1180  cell->active_fe_index());
1181  for (unsigned int i = 0; i < dofs_per_face; ++i)
1182  boundary_dofs.insert(dofs_on_face[i]);
1183  }
1184  return boundary_dofs.size();
1185  }
1186 
1187 
1188 
1189  template <int dim, int spacedim>
1192  const std::set<types::boundary_id> &boundary_ids) const
1193  {
1194  Assert(fe_collection.size() > 0, ExcNoFESelected());
1195  Assert(boundary_ids.find(numbers::internal_face_boundary_id) ==
1196  boundary_ids.end(),
1198 
1199  // same as above, but with additional checks for set of boundary
1200  // indicators
1201  std::set<types::global_dof_index> boundary_dofs;
1202  std::vector<types::global_dof_index> dofs_on_face;
1203  dofs_on_face.reserve(this->get_fe_collection().max_dofs_per_face());
1204 
1205  typename HpDoFHandler<dim, spacedim>::active_cell_iterator
1206  cell = this->begin_active(),
1207  endc = this->end();
1208  for (; cell != endc; ++cell)
1209  for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
1210  if (cell->at_boundary(f) &&
1211  (boundary_ids.find(cell->face(f)->boundary_id()) !=
1212  boundary_ids.end()))
1213  {
1214  const unsigned int dofs_per_face = cell->get_fe().dofs_per_face;
1215  dofs_on_face.resize(dofs_per_face);
1216 
1217  cell->face(f)->get_dof_indices(dofs_on_face,
1218  cell->active_fe_index());
1219  for (unsigned int i = 0; i < dofs_per_face; ++i)
1220  boundary_dofs.insert(dofs_on_face[i]);
1221  }
1222  return boundary_dofs.size();
1223  }
1224 
1225 
1226 
1227  template <>
1230  {
1231  Assert(false, ExcNotImplemented());
1232  return 0;
1233  }
1234 
1235 
1236 
1237  template <>
1238  template <typename number>
1241  const std::map<types::boundary_id, const Function<3, number> *> &) const
1242  {
1243  Assert(false, ExcNotImplemented());
1244  return 0;
1245  }
1246 
1247 
1248 
1249  template <>
1251  DoFHandler<2, 3>::n_boundary_dofs(const std::set<types::boundary_id> &) const
1252  {
1253  Assert(false, ExcNotImplemented());
1254  return 0;
1255  }
1256 
1257 
1258 
1259  template <int dim, int spacedim>
1260  std::size_t
1262  {
1263  std::size_t mem =
1273  for (unsigned int i = 0; i < levels.size(); ++i)
1276 
1277  return mem;
1278  }
1279 
1280 
1281 
1282  template <int dim, int spacedim>
1283  void
1285  const std::vector<unsigned int> &active_fe_indices)
1286  {
1287  Assert(active_fe_indices.size() == get_triangulation().n_active_cells(),
1288  ExcDimensionMismatch(active_fe_indices.size(),
1289  get_triangulation().n_active_cells()));
1290 
1292  // we could set the values directly, since they are stored as
1293  // protected data of this object, but for simplicity we use the
1294  // cell-wise access. this way we also have to pass some debug-mode
1295  // tests which we would have to duplicate ourselves otherwise
1296  active_cell_iterator cell = begin_active(), endc = end();
1297  for (unsigned int i = 0; cell != endc; ++cell, ++i)
1298  if (cell->is_locally_owned())
1299  cell->set_active_fe_index(active_fe_indices[i]);
1300  }
1301 
1302 
1303 
1304  template <int dim, int spacedim>
1305  void
1307  std::vector<unsigned int> &active_fe_indices) const
1308  {
1309  active_fe_indices.resize(get_triangulation().n_active_cells());
1310 
1311  // we could try to extract the values directly, since they are
1312  // stored as protected data of this object, but for simplicity we
1313  // use the cell-wise access.
1314  active_cell_iterator cell = begin_active(), endc = end();
1315  for (unsigned int i = 0; cell != endc; ++cell, ++i)
1316  active_fe_indices[i] = cell->active_fe_index();
1317  }
1318 
1319  template <int dim, int spacedim>
1320  void
1324  {
1325  if (this->tria != &tria)
1326  {
1327  for (unsigned int i = 0; i < tria_listeners.size(); ++i)
1328  tria_listeners[i].disconnect();
1329  tria_listeners.clear();
1330 
1331  this->tria = &tria;
1332 
1334  }
1335 
1337 
1338  distribute_dofs(fe);
1339  }
1340 
1341 
1342  template <int dim, int spacedim>
1343  void
1346  {
1347  Assert(
1348  tria != nullptr,
1349  ExcMessage(
1350  "You need to set the Triangulation in the DoFHandler using initialize() or "
1351  "in the constructor before you can distribute DoFs."));
1352  Assert(tria->n_levels() > 0,
1353  ExcMessage("The Triangulation you are using is empty!"));
1354  Assert(ff.size() > 0, ExcMessage("The hp::FECollection given is empty!"));
1355 
1356  // don't create a new object if the one we have is already appropriate
1357  if (fe_collection != ff)
1359 
1360  // at the beginning, make sure every processor knows the
1361  // active_fe_indices on both its own cells and all ghost cells
1364 
1365  // If an underlying shared::Tria allows artificial cells,
1366  // then save the current set of subdomain ids, and set
1367  // subdomain ids to the "true" owner of each cell. we later
1368  // restore these flags
1369  std::vector<types::subdomain_id> saved_subdomain_ids;
1370  if (const parallel::shared::Triangulation<dim, spacedim> *shared_tria =
1371  (dynamic_cast<const parallel::shared::Triangulation<dim, spacedim> *>(
1372  &get_triangulation())))
1373  if (shared_tria->with_artificial_cells())
1374  {
1375  saved_subdomain_ids.resize(shared_tria->n_active_cells());
1376 
1378  active_cell_iterator cell = get_triangulation().begin_active(),
1379  endc = get_triangulation().end();
1380 
1381  const std::vector<types::subdomain_id> &true_subdomain_ids =
1382  shared_tria->get_true_subdomain_ids_of_cells();
1383 
1384  for (unsigned int index = 0; cell != endc; ++cell, ++index)
1385  {
1386  saved_subdomain_ids[index] = cell->subdomain_id();
1387  cell->set_subdomain_id(true_subdomain_ids[index]);
1388  }
1389  }
1390 
1391  // ensure that the active_fe_indices vectors are
1392  // initialized correctly.
1394 
1395  // up front make sure that the fe collection is large enough to
1396  // cover all fe indices presently in use on the mesh
1397  for (active_cell_iterator cell = begin_active(); cell != end(); ++cell)
1398  if (cell->is_locally_owned())
1399  Assert(cell->active_fe_index() < fe_collection.size(),
1400  ExcInvalidFEIndex(cell->active_fe_index(),
1401  fe_collection.size()));
1402 
1403 
1404  // then allocate space for all the other tables
1406  reserve_space(*this);
1407 
1408 
1409  // now undo the subdomain modification
1410  if (const parallel::shared::Triangulation<dim, spacedim> *shared_tria =
1411  (dynamic_cast<const parallel::shared::Triangulation<dim, spacedim> *>(
1412  &get_triangulation())))
1413  if (shared_tria->with_artificial_cells())
1414  {
1416  active_cell_iterator cell = get_triangulation().begin_active(),
1417  endc = get_triangulation().end();
1418 
1419  for (unsigned int index = 0; cell != endc; ++cell, ++index)
1420  cell->set_subdomain_id(saved_subdomain_ids[index]);
1421  }
1422 
1423 
1424  // Clear user flags because we will need them. But first we save
1425  // them and make sure that we restore them later such that at the
1426  // end of this function the Triangulation will be in the same
1427  // state as it was at the beginning of this function.
1428  std::vector<bool> user_flags;
1429  tria->save_user_flags(user_flags);
1430  const_cast<Triangulation<dim, spacedim> &>(*tria).clear_user_flags();
1431 
1432 
1434 
1435  // Now for the real work:
1436  number_cache = policy->distribute_dofs();
1437 
1439 
1440  // do some housekeeping: compress indices
1441  {
1443  for (int level = levels.size() - 1; level >= 0; --level)
1444  tg += Threads::new_task(
1445  &::internal::hp::DoFLevel::compress_data<dim, spacedim>,
1446  *levels[level],
1447  fe_collection);
1448  tg.join_all();
1449  }
1450 
1451  // finally restore the user flags
1452  const_cast<Triangulation<dim, spacedim> &>(*tria).load_user_flags(
1453  user_flags);
1454  }
1455 
1456 
1457  template <int dim, int spacedim>
1458  void
1460  {
1461  // decide whether we need a sequential or a parallel shared/distributed
1462  // policy
1463  if (dynamic_cast<const parallel::shared::Triangulation<dim, spacedim> *>(
1464  &*this->tria) != nullptr)
1465  policy =
1466  std_cxx14::make_unique<internal::DoFHandlerImplementation::Policy::
1468  *this);
1469  else if (dynamic_cast<
1471  &*this->tria) != nullptr)
1472  policy = std_cxx14::make_unique<
1474  DoFHandler<dim, spacedim>>>(*this);
1475  else
1476  policy =
1477  std_cxx14::make_unique<internal::DoFHandlerImplementation::Policy::
1479 
1480  tria_listeners.push_back(this->tria->signals.pre_refinement.connect(
1482  std::ref(*this))));
1483  tria_listeners.push_back(this->tria->signals.post_refinement.connect(
1485  std::ref(*this))));
1486  tria_listeners.push_back(this->tria->signals.create.connect(std::bind(
1488  }
1489 
1490 
1491  template <int dim, int spacedim>
1492  void
1494  {
1495  // release memory
1496  clear_space();
1497  }
1498 
1499 
1500 
1501  template <int dim, int spacedim>
1502  void
1504  const std::vector<types::global_dof_index> &new_numbers)
1505  {
1506  Assert(levels.size() > 0,
1507  ExcMessage(
1508  "You need to distribute DoFs before you can renumber them."));
1509 
1510  AssertDimension(new_numbers.size(), n_locally_owned_dofs());
1511 
1512 #ifdef DEBUG
1513  // assert that the new indices are
1514  // consecutively numbered if we are
1515  // working on a single
1516  // processor. this doesn't need to
1517  // hold in the case of a parallel
1518  // mesh since we map the interval
1519  // [0...n_dofs()) into itself but
1520  // only globally, not on each
1521  // processor
1522  if (n_locally_owned_dofs() == n_dofs())
1523  {
1524  std::vector<types::global_dof_index> tmp(new_numbers);
1525  std::sort(tmp.begin(), tmp.end());
1526  std::vector<types::global_dof_index>::const_iterator p = tmp.begin();
1528  for (; p != tmp.end(); ++p, ++i)
1529  Assert(*p == i, ExcNewNumbersNotConsecutive(i));
1530  }
1531  else
1532  for (types::global_dof_index i = 0; i < new_numbers.size(); ++i)
1533  Assert(new_numbers[i] < n_dofs(),
1534  ExcMessage(
1535  "New DoF index is not less than the total number of dofs."));
1536 #endif
1537 
1538  // uncompress the internal storage scheme of dofs on cells so that
1539  // we can access dofs in turns. uncompress in parallel, starting
1540  // with the most expensive levels (the highest ones)
1541  {
1543  for (int level = levels.size() - 1; level >= 0; --level)
1544  tg += Threads::new_task(
1545  &::internal::hp::DoFLevel::uncompress_data<dim, spacedim>,
1546  *levels[level],
1547  fe_collection);
1548  tg.join_all();
1549  }
1550 
1551  // do the renumbering
1552  number_cache = policy->renumber_dofs(new_numbers);
1553 
1554  // now re-compress the dof indices
1555  {
1557  for (int level = levels.size() - 1; level >= 0; --level)
1558  tg += Threads::new_task(
1559  &::internal::hp::DoFLevel::compress_data<dim, spacedim>,
1560  *levels[level],
1561  fe_collection);
1562  tg.join_all();
1563  }
1564  }
1565 
1566 
1567 
1568  template <int dim, int spacedim>
1569  unsigned int
1571  {
1572  Assert(fe_collection.size() > 0, ExcNoFESelected());
1573  return ::internal::hp::DoFHandlerImplementation::Implementation::
1574  max_couplings_between_dofs(*this);
1575  }
1576 
1577 
1578 
1579  template <int dim, int spacedim>
1580  unsigned int
1582  {
1583  Assert(fe_collection.size() > 0, ExcNoFESelected());
1584 
1585  switch (dim)
1586  {
1587  case 1:
1588  return fe_collection.max_dofs_per_vertex();
1589  case 2:
1590  return (3 * fe_collection.max_dofs_per_vertex() +
1591  2 * fe_collection.max_dofs_per_line());
1592  case 3:
1593  // we need to take refinement of one boundary face into
1594  // consideration here; in fact, this function returns what
1595  // #max_coupling_between_dofs<2> returns
1596  //
1597  // we assume here, that only four faces meet at the boundary;
1598  // this assumption is not justified and needs to be fixed some
1599  // time. fortunately, omitting it for now does no harm since
1600  // the matrix will cry foul if its requirements are not
1601  // satisfied
1602  return (19 * fe_collection.max_dofs_per_vertex() +
1603  28 * fe_collection.max_dofs_per_line() +
1604  8 * fe_collection.max_dofs_per_quad());
1605  default:
1606  Assert(false, ExcNotImplemented());
1607  return 0;
1608  }
1609  }
1610 
1611 
1612 
1613  template <int dim, int spacedim>
1614  void
1616  {
1617  // Create sufficiently many hp::DoFLevels.
1618  while (levels.size() < tria->n_levels())
1619  levels.emplace_back(new ::internal::hp::DoFLevel);
1620 
1621  // then make sure that on each level we have the appropriate size
1622  // of active_fe_indices; preset them to zero, i.e. the default FE
1623  for (unsigned int level = 0; level < levels.size(); ++level)
1624  {
1625  if (levels[level]->active_fe_indices.size() == 0)
1626  levels[level]->active_fe_indices.resize(tria->n_raw_cells(level), 0);
1627  else
1628  {
1629  // Either the active_fe_indices have size zero because
1630  // they were just created, or the correct size. Other
1631  // sizes indicate that something went wrong.
1632  Assert(levels[level]->active_fe_indices.size() ==
1633  tria->n_raw_cells(level),
1634  ExcInternalError());
1635  }
1636 
1637  // it may be that the previous table was compressed; in that
1638  // case, restore the correct active_fe_index. the fact that
1639  // this no longer matches the indices in the table is of no
1640  // importance because the current function is called at a
1641  // point where we have to recreate the dof_indices tables in
1642  // the levels anyway
1643  levels[level]->normalize_active_fe_indices();
1644  }
1645  }
1646 
1647 
1648  template <int dim, int spacedim>
1649  void
1651  {
1653 
1654  // Remember if the cells already have children. That will make the
1655  // transfer of the active_fe_index to the finer levels easier.
1656  Assert(has_children.size() == 0, ExcInternalError());
1657  for (unsigned int i = 0; i < levels.size(); ++i)
1658  {
1659  const unsigned int cells_on_level = tria->n_raw_cells(i);
1660  std::unique_ptr<std::vector<bool>> has_children_level(
1661  new std::vector<bool>(cells_on_level));
1662 
1663  // Check for each cell, if it has children. in 1d, we don't
1664  // store refinement cases, so use the 'children' vector
1665  // instead
1666  if (dim == 1)
1667  std::transform(tria->levels[i]->cells.children.begin(),
1668  tria->levels[i]->cells.children.end(),
1669  has_children_level->begin(),
1670  std::bind(std::not_equal_to<int>(),
1671  std::placeholders::_1,
1672  -1));
1673  else
1674  std::transform(tria->levels[i]->cells.refinement_cases.begin(),
1675  tria->levels[i]->cells.refinement_cases.end(),
1676  has_children_level->begin(),
1677  std::bind(std::not_equal_to<unsigned char>(),
1678  std::placeholders::_1,
1679  static_cast<unsigned char>(
1681 
1682  has_children.emplace_back(std::move(has_children_level));
1683  }
1684  }
1685 
1686 
1687 
1688  template <int dim, int spacedim>
1689  void
1691  {
1692  Assert(has_children.size() == levels.size(), ExcInternalError());
1693 
1694  // Normally only one level is added, but if this Triangulation
1695  // is created by copy_triangulation, it can be more than one level.
1696  while (levels.size() < tria->n_levels())
1697  levels.emplace_back(new ::internal::hp::DoFLevel);
1698 
1699  // Coarsening can lead to the loss of levels. Hence remove them.
1700  while (levels.size() > tria->n_levels())
1701  {
1702  // drop the last element. that also releases the memory pointed to
1703  levels.pop_back();
1704  }
1705 
1706  Assert(levels.size() == tria->n_levels(), ExcInternalError());
1707 
1708  // Resize active_fe_indices vectors. use zero indicator to extend
1709  for (unsigned int i = 0; i < levels.size(); ++i)
1710  levels[i]->active_fe_indices.resize(tria->n_raw_cells(i), 0);
1711 
1712  // if a finite element collection has already been set, then
1713  // actually try to set active_fe_indices for child cells of
1714  // refined cells to the active_fe_index of the mother cell. if no
1715  // finite element collection has been assigned yet, then all
1716  // indicators are zero anyway, and there is no point trying to set
1717  // anything (besides, we would trip over an assertion in
1718  // set_active_fe_index)
1719  if (fe_collection.size() > 0)
1720  {
1721  cell_iterator cell = begin(), endc = end();
1722  for (; cell != endc; ++cell)
1723  {
1724  // Look if the cell got children during refinement by
1725  // checking whether it has children now but didn't have
1726  // children before refinement (the has_children array is
1727  // set in pre-refinement action)
1728  //
1729  // Note: Although one level is added to the DoFHandler
1730  // levels, when the triangulation got one, for the buffer
1731  // has_children this new level is not required, because
1732  // the cells on the finest level never have
1733  // children. Hence cell->has_children () will always
1734  // return false on that level, which would cause shortcut
1735  // evaluation of the following expression. Thus an index
1736  // error in has_children should never occur.
1737  if (cell->has_children() &&
1738  !(*has_children[cell->level()])[cell->index()])
1739  {
1740  // Set active_fe_index in children to the same value
1741  // as in the parent cell. we can't access the
1742  // active_fe_index in the parent cell any more through
1743  // cell->active_fe_index() since that function is not
1744  // allowed for inactive cells, but we can access this
1745  // information from the DoFLevels directly
1746  //
1747  // we don't have to set the active_fe_index for ghost
1748  // cells -- these will be exchanged automatically
1749  // upon distribute_dofs()
1750  for (unsigned int i = 0; i < cell->n_children(); ++i)
1751  if (cell->child(i)->is_locally_owned())
1752  cell->child(i)->set_active_fe_index(
1753  levels[cell->level()]->active_fe_index(cell->index()));
1754  }
1755  }
1756  }
1757 
1758  // Free buffer objects
1759  has_children.clear();
1760  }
1761 
1762 
1763  template <int dim, int spacedim>
1764  template <int structdim>
1766  DoFHandler<dim, spacedim>::get_dof_index(const unsigned int,
1767  const unsigned int,
1768  const unsigned int,
1769  const unsigned int) const
1770  {
1771  Assert(false, ExcNotImplemented());
1773  }
1774 
1775 
1776  template <int dim, int spacedim>
1777  template <int structdim>
1778  void
1779  DoFHandler<dim, spacedim>::set_dof_index(const unsigned int,
1780  const unsigned int,
1781  const unsigned int,
1782  const unsigned int,
1783  const types::global_dof_index) const
1784  {
1785  Assert(false, ExcNotImplemented());
1786  }
1787 
1788 
1789  template <int dim, int spacedim>
1790  void
1792  {
1793  levels.clear();
1794  faces.reset();
1795 
1796  vertex_dofs = std::move(std::vector<types::global_dof_index>());
1797  vertex_dof_offsets = std::move(std::vector<unsigned int>());
1798  }
1799 } // namespace hp
1800 
1801 
1802 
1803 /*-------------- Explicit Instantiations -------------------------------*/
1804 #include "dof_handler.inst"
1805 
1806 
1807 DEAL_II_NAMESPACE_CLOSE
active_cell_iterator end_active(const unsigned int level) const
Definition: dof_handler.cc:979
virtual void distribute_dofs(const hp::FECollection< dim, spacedim > &fe)
static const unsigned int invalid_unsigned_int
Definition: types.h:173
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: dof_handler.cc:943
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1366
static::ExceptionBase & ExcInvalidFEIndex(int arg1, int arg2)
void clear_user_flags()
Definition: tria.cc:11096
unsigned int max_couplings_between_dofs() const
types::global_dof_index n_dofs() const
void renumber_dofs(const std::vector< types::global_dof_index > &new_numbers)
unsigned int offset_type
Definition: dof_level.h:112
unsigned int max_couplings_between_boundary_dofs() const
void set_active_fe_indices(const std::vector< unsigned int > &active_fe_indices)
Task< RT > new_task(const std::function< RT()> &function)
::internal::DoFHandlerImplementation::NumberCache number_cache
Definition: dof_handler.h:954
active_cell_iterator end_active(const unsigned int level) const
std::vector< std::unique_ptr<::internal::hp::DoFLevel > > levels
Definition: dof_handler.h:939
void load_user_flags(std::istream &in)
Definition: tria.cc:11156
std::vector< std::unique_ptr< std::vector< bool > > > has_children
Definition: dof_handler.h:998
void pre_refinement_action()
std::unique_ptr<::internal::hp::DoFIndicesOnFaces< dim > > faces
Definition: dof_handler.h:945
void setup_policy_and_listeners()
static::ExceptionBase & ExcNewNumbersNotConsecutive(types::global_dof_index arg1)
unsigned int size() const
std::vector< boost::signals2::connection > tria_listeners
Definition: dof_handler.h:1004
virtual void clear()
static void reserve_space_cells(DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_handler.cc:247
unsigned long long int global_dof_index
Definition: types.h:72
SmartPointer< const Triangulation< dim, spacedim >, DoFHandler< dim, spacedim > > tria
Definition: dof_handler.h:870
std::unique_ptr<::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > > policy
Definition: dof_handler.h:883
active_cell_iterator begin_active(const unsigned int level=0) const
typename ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:194
IteratorRange< active_cell_iterator > active_cell_iterators_on_level(const unsigned int level) const
static::ExceptionBase & ExcNoFESelected()
static::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:207
const Triangulation< dim, spacedim > & get_triangulation() const
T sum(const T &t, const MPI_Comm &mpi_communicator)
#define Assert(cond, exc)
Definition: exceptions.h:1227
hp::FECollection< dim, spacedim > fe_collection
Definition: dof_handler.h:1090
types::global_dof_index n_dofs() const
static::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::vector< types::global_dof_index > dofs
Definition: dof_faces.h:111
static void reserve_space_release_space(DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_handler.cc:91
virtual ~DoFHandler() override
SmartPointer< const Triangulation< dim, spacedim >, DoFHandler< dim, spacedim > > tria
Definition: dof_handler.h:1083
virtual std::size_t memory_consumption() const
static unsigned int max_couplings_between_dofs(const DoFHandler< 1, spacedim > &dof_handler)
Definition: dof_handler.cc:796
std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces< dim > > faces
Definition: dof_handler.h:1249
cell_iterator end() const
Definition: dof_handler.cc:959
std::vector< unsigned int > dof_offsets
Definition: dof_faces.h:105
static::ExceptionBase & ExcInvalidBoundaryIndicator()
Definition: hp.h:102
std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel< dim > > > levels
Definition: dof_handler.h:1237
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
internal::hp::DoFIndicesOnFacesOrEdges< 2 > quads
Definition: dof_faces.h:320
IteratorRange< cell_iterator > cell_iterators() const
cell_iterator end() const
void initialize(const Triangulation< dim, spacedim > &tria, const hp::FECollection< dim, spacedim > &fe)
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
hp::FECollection< dim, spacedim > fe_collection
Definition: dof_handler.h:875
void clear_space()
unsigned int max_dofs_per_face(const hp::DoFHandler< dim, spacedim > &dh)
static void reserve_space_vertices(DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_handler.cc:129
static::ExceptionBase & ExcNotImplemented()
cell_iterator begin(const unsigned int level=0) const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
Iterator points to a valid object.
void exchange_cell_data_to_ghosts(const MeshType &mesh, const std::function< boost::optional< DataType >(const typename MeshType::active_cell_iterator &)> &pack, const std::function< void(const typename MeshType::active_cell_iterator &, const DataType &)> &unpack)
const types::boundary_id internal_face_boundary_id
Definition: types.h:223
static void communicate_active_fe_indices(::hp::DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_handler.cc:911
static void reserve_space(DoFHandler< 1, spacedim > &dof_handler)
Definition: dof_handler.cc:628
const types::global_dof_index invalid_dof_index
Definition: types.h:188
void create_active_fe_table()
std::vector< types::global_dof_index > vertex_dofs
Definition: dof_handler.h:976
types::global_dof_index n_locally_owned_dofs() const
IteratorRange< active_cell_iterator > active_cell_iterators() const
unsigned int boundary_id
Definition: types.h:111
static::ExceptionBase & ExcNoFESelected()
static void reserve_space_faces(DoFHandler< dim, spacedim > &dof_handler)
Definition: dof_handler.cc:342
std::vector< unsigned int > vertex_dof_offsets
Definition: dof_handler.h:990
types::global_dof_index n_boundary_dofs() const
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:1223
static::ExceptionBase & ExcInternalError()
void get_active_fe_indices(std::vector< unsigned int > &active_fe_indices) const