Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
dof_handler.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2020 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 
19 
22 
24 #include <deal.II/dofs/dof_faces.h>
28 
29 #include <deal.II/fe/fe.h>
30 
31 #include <deal.II/grid/tria.h>
35 
36 #include <algorithm>
37 #include <set>
38 #include <unordered_set>
39 
41 
42 
43 template <int dim, int spacedim>
44 const unsigned int DoFHandler<dim, spacedim>::dimension;
45 
46 template <int dim, int spacedim>
48 
49 template <int dim, int spacedim>
51 
52 
53 // reference the invalid_dof_index variable explicitly to work around
54 // a bug in the icc8 compiler
55 namespace internal
56 {
57  template <int dim, int spacedim>
60  {
62  }
63 } // namespace internal
64 
65 
66 
67 namespace internal
68 {
69  template <int dim, int spacedim>
70  std::string
71  policy_to_string(const ::internal::DoFHandlerImplementation::Policy::
72  PolicyBase<dim, spacedim> &policy)
73  {
74  std::string policy_name;
75  if (dynamic_cast<const typename ::internal::DoFHandlerImplementation::
76  Policy::Sequential<::DoFHandler<dim, spacedim>> *>(
77  &policy) ||
78  dynamic_cast<
79  const typename ::internal::DoFHandlerImplementation::Policy::
80  Sequential<::hp::DoFHandler<dim, spacedim>> *>(&policy))
81  policy_name = "Policy::Sequential<";
82  else if (dynamic_cast<
83  const typename ::internal::DoFHandlerImplementation::
84  Policy::ParallelDistributed<::DoFHandler<dim, spacedim>>
85  *>(&policy) ||
86  dynamic_cast<
87  const typename ::internal::DoFHandlerImplementation::
88  Policy::ParallelDistributed<
90  policy_name = "Policy::ParallelDistributed<";
91  else if (dynamic_cast<
92  const typename ::internal::DoFHandlerImplementation::
93  Policy::ParallelShared<::DoFHandler<dim, spacedim>> *>(
94  &policy) ||
95  dynamic_cast<const typename ::internal::
98  &policy))
99  policy_name = "Policy::ParallelShared<";
100  else
101  AssertThrow(false, ExcNotImplemented());
102  policy_name += Utilities::int_to_string(dim) + "," +
103  Utilities::int_to_string(spacedim) + ">";
104  return policy_name;
105  }
106 
107 
108  namespace DoFHandlerImplementation
109  {
110  // access class
111  // ::DoFHandler instead of
112  // namespace internal::DoFHandler
113  using ::DoFHandler;
114 
115 
121  {
126  template <int spacedim>
127  static unsigned int
129  {
130  return std::min(static_cast<types::global_dof_index>(
131  3 * dof_handler.get_fe().dofs_per_vertex +
132  2 * dof_handler.get_fe().dofs_per_line),
133  dof_handler.n_dofs());
134  }
135 
136 
137 
138  template <int spacedim>
139  static unsigned int
141  {
142  // get these numbers by drawing pictures
143  // and counting...
144  // example:
145  // | | |
146  // --x-----x--x--X--
147  // | | | |
148  // | x--x--x
149  // | | | |
150  // --x--x--*--x--x--
151  // | | | |
152  // x--x--x |
153  // | | | |
154  // --X--x--x-----x--
155  // | | |
156  // x = vertices connected with center vertex *;
157  // = total of 19
158  // (the X vertices are connected with * if
159  // the vertices adjacent to X are hanging
160  // nodes)
161  // count lines -> 28 (don't forget to count
162  // mother and children separately!)
163  types::global_dof_index max_couplings;
164  switch (dof_handler.tria->max_adjacent_cells())
165  {
166  case 4:
167  max_couplings = 19 * dof_handler.get_fe().dofs_per_vertex +
168  28 * dof_handler.get_fe().dofs_per_line +
169  8 * dof_handler.get_fe().dofs_per_quad;
170  break;
171  case 5:
172  max_couplings = 21 * dof_handler.get_fe().dofs_per_vertex +
173  31 * dof_handler.get_fe().dofs_per_line +
174  9 * dof_handler.get_fe().dofs_per_quad;
175  break;
176  case 6:
177  max_couplings = 28 * dof_handler.get_fe().dofs_per_vertex +
178  42 * dof_handler.get_fe().dofs_per_line +
179  12 * dof_handler.get_fe().dofs_per_quad;
180  break;
181  case 7:
182  max_couplings = 30 * dof_handler.get_fe().dofs_per_vertex +
183  45 * dof_handler.get_fe().dofs_per_line +
184  13 * dof_handler.get_fe().dofs_per_quad;
185  break;
186  case 8:
187  max_couplings = 37 * dof_handler.get_fe().dofs_per_vertex +
188  56 * dof_handler.get_fe().dofs_per_line +
189  16 * dof_handler.get_fe().dofs_per_quad;
190  break;
191 
192  // the following numbers are not based on actual counting but by
193  // extrapolating the number sequences from the previous ones (for
194  // example, for dofs_per_vertex, the sequence above is 19, 21, 28,
195  // 30, 37, and is continued as follows):
196  case 9:
197  max_couplings = 39 * dof_handler.get_fe().dofs_per_vertex +
198  59 * dof_handler.get_fe().dofs_per_line +
199  17 * dof_handler.get_fe().dofs_per_quad;
200  break;
201  case 10:
202  max_couplings = 46 * dof_handler.get_fe().dofs_per_vertex +
203  70 * dof_handler.get_fe().dofs_per_line +
204  20 * dof_handler.get_fe().dofs_per_quad;
205  break;
206  case 11:
207  max_couplings = 48 * dof_handler.get_fe().dofs_per_vertex +
208  73 * dof_handler.get_fe().dofs_per_line +
209  21 * dof_handler.get_fe().dofs_per_quad;
210  break;
211  case 12:
212  max_couplings = 55 * dof_handler.get_fe().dofs_per_vertex +
213  84 * dof_handler.get_fe().dofs_per_line +
214  24 * dof_handler.get_fe().dofs_per_quad;
215  break;
216  case 13:
217  max_couplings = 57 * dof_handler.get_fe().dofs_per_vertex +
218  87 * dof_handler.get_fe().dofs_per_line +
219  25 * dof_handler.get_fe().dofs_per_quad;
220  break;
221  case 14:
222  max_couplings = 63 * dof_handler.get_fe().dofs_per_vertex +
223  98 * dof_handler.get_fe().dofs_per_line +
224  28 * dof_handler.get_fe().dofs_per_quad;
225  break;
226  case 15:
227  max_couplings = 65 * dof_handler.get_fe().dofs_per_vertex +
228  103 * dof_handler.get_fe().dofs_per_line +
229  29 * dof_handler.get_fe().dofs_per_quad;
230  break;
231  case 16:
232  max_couplings = 72 * dof_handler.get_fe().dofs_per_vertex +
233  114 * dof_handler.get_fe().dofs_per_line +
234  32 * dof_handler.get_fe().dofs_per_quad;
235  break;
236 
237  default:
238  Assert(false, ExcNotImplemented());
239  max_couplings = 0;
240  }
241  return std::min(max_couplings, dof_handler.n_dofs());
242  }
243 
244 
245  template <int spacedim>
246  static unsigned int
248  {
249  // TODO:[?] Invent significantly better estimates than the ones in this
250  // function
251 
252  // doing the same thing here is a
253  // rather complicated thing, compared
254  // to the 2d case, since it is hard
255  // to draw pictures with several
256  // refined hexahedra :-) so I
257  // presently only give a coarse
258  // estimate for the case that at most
259  // 8 hexes meet at each vertex
260  //
261  // can anyone give better estimate
262  // here?
263  const unsigned int max_adjacent_cells =
264  dof_handler.tria->max_adjacent_cells();
265 
266  types::global_dof_index max_couplings;
267  if (max_adjacent_cells <= 8)
268  max_couplings = 7 * 7 * 7 * dof_handler.get_fe().dofs_per_vertex +
269  7 * 6 * 7 * 3 * dof_handler.get_fe().dofs_per_line +
270  9 * 4 * 7 * 3 * dof_handler.get_fe().dofs_per_quad +
271  27 * dof_handler.get_fe().dofs_per_hex;
272  else
273  {
274  Assert(false, ExcNotImplemented());
275  max_couplings = 0;
276  }
277 
278  return std::min(max_couplings, dof_handler.n_dofs());
279  }
280 
281 
291  template <int spacedim>
292  static void reserve_space(DoFHandler<1, spacedim> &dof_handler)
293  {
294  dof_handler.vertex_dofs.resize(dof_handler.tria->n_vertices() *
295  dof_handler.get_fe().dofs_per_vertex,
297 
298  for (unsigned int i = 0; i < dof_handler.tria->n_levels(); ++i)
299  {
300  dof_handler.levels.emplace_back(
302 
303  dof_handler.levels.back()->dof_object.dofs.resize(
304  dof_handler.tria->n_raw_cells(i) *
305  dof_handler.get_fe().dofs_per_line,
307 
308  dof_handler.levels.back()->cell_dof_indices_cache.resize(
309  dof_handler.tria->n_raw_cells(i) *
310  dof_handler.get_fe().dofs_per_cell,
312  }
313  }
314 
315 
316  template <int spacedim>
317  static void reserve_space(DoFHandler<2, spacedim> &dof_handler)
318  {
319  dof_handler.vertex_dofs.resize(dof_handler.tria->n_vertices() *
320  dof_handler.get_fe().dofs_per_vertex,
322 
323  for (unsigned int i = 0; i < dof_handler.tria->n_levels(); ++i)
324  {
325  dof_handler.levels.emplace_back(
327 
328  dof_handler.levels.back()->dof_object.dofs.resize(
329  dof_handler.tria->n_raw_cells(i) *
330  dof_handler.get_fe().dofs_per_quad,
332 
333  dof_handler.levels.back()->cell_dof_indices_cache.resize(
334  dof_handler.tria->n_raw_cells(i) *
335  dof_handler.get_fe().dofs_per_cell,
337  }
338 
339  dof_handler.faces = std_cxx14::make_unique<
341  // avoid access to n_raw_lines when there are no cells
342  if (dof_handler.tria->n_cells() > 0)
343  {
344  dof_handler.faces->lines.dofs.resize(
345  dof_handler.tria->n_raw_lines() *
346  dof_handler.get_fe().dofs_per_line,
348  }
349  }
350 
351 
352  template <int spacedim>
353  static void reserve_space(DoFHandler<3, spacedim> &dof_handler)
354  {
355  dof_handler.vertex_dofs.resize(dof_handler.tria->n_vertices() *
356  dof_handler.get_fe().dofs_per_vertex,
358 
359  for (unsigned int i = 0; i < dof_handler.tria->n_levels(); ++i)
360  {
361  dof_handler.levels.emplace_back(
363 
364  dof_handler.levels.back()->dof_object.dofs.resize(
365  dof_handler.tria->n_raw_cells(i) *
366  dof_handler.get_fe().dofs_per_hex,
368 
369  dof_handler.levels.back()->cell_dof_indices_cache.resize(
370  dof_handler.tria->n_raw_cells(i) *
371  dof_handler.get_fe().dofs_per_cell,
373  }
374  dof_handler.faces = std_cxx14::make_unique<
376 
377  // avoid access to n_raw_lines when there are no cells
378  if (dof_handler.tria->n_cells() > 0)
379  {
380  dof_handler.faces->lines.dofs.resize(
381  dof_handler.tria->n_raw_lines() *
382  dof_handler.get_fe().dofs_per_line,
384  dof_handler.faces->quads.dofs.resize(
385  dof_handler.tria->n_raw_quads() *
386  dof_handler.get_fe().dofs_per_quad,
388  }
389  }
390 
391  template <int spacedim>
392  static void reserve_space_mg(DoFHandler<1, spacedim> &dof_handler)
393  {
394  Assert(dof_handler.get_triangulation().n_levels() > 0,
395  ExcMessage("Invalid triangulation"));
396  dof_handler.clear_mg_space();
397 
398  const ::Triangulation<1, spacedim> &tria =
399  dof_handler.get_triangulation();
400  const unsigned int dofs_per_line = dof_handler.get_fe().dofs_per_line;
401  const unsigned int n_levels = tria.n_levels();
402 
403  for (unsigned int i = 0; i < n_levels; ++i)
404  {
405  dof_handler.mg_levels.emplace_back(
407  dof_handler.mg_levels.back()->dof_object.dofs =
408  std::vector<types::global_dof_index>(tria.n_raw_lines(i) *
409  dofs_per_line,
411  }
412 
413  const unsigned int n_vertices = tria.n_vertices();
414 
415  dof_handler.mg_vertex_dofs.resize(n_vertices);
416 
417  std::vector<unsigned int> max_level(n_vertices, 0);
418  std::vector<unsigned int> min_level(n_vertices, n_levels);
419 
420  for (typename ::Triangulation<1, spacedim>::cell_iterator cell =
421  tria.begin();
422  cell != tria.end();
423  ++cell)
424  {
425  const unsigned int level = cell->level();
426 
427  for (const unsigned int vertex : GeometryInfo<1>::vertex_indices())
428  {
429  const unsigned int vertex_index = cell->vertex_index(vertex);
430 
431  if (min_level[vertex_index] > level)
432  min_level[vertex_index] = level;
433 
434  if (max_level[vertex_index] < level)
435  max_level[vertex_index] = level;
436  }
437  }
438 
439  for (unsigned int vertex = 0; vertex < n_vertices; ++vertex)
440  if (tria.vertex_used(vertex))
441  {
442  Assert(min_level[vertex] < n_levels, ExcInternalError());
443  Assert(max_level[vertex] >= min_level[vertex],
444  ExcInternalError());
445  dof_handler.mg_vertex_dofs[vertex].init(
446  min_level[vertex],
447  max_level[vertex],
448  dof_handler.get_fe().dofs_per_vertex);
449  }
450 
451  else
452  {
453  Assert(min_level[vertex] == n_levels, ExcInternalError());
454  Assert(max_level[vertex] == 0, ExcInternalError());
455  dof_handler.mg_vertex_dofs[vertex].init(1, 0, 0);
456  }
457  }
458 
459  template <int spacedim>
460  static void reserve_space_mg(DoFHandler<2, spacedim> &dof_handler)
461  {
462  Assert(dof_handler.get_triangulation().n_levels() > 0,
463  ExcMessage("Invalid triangulation"));
464  dof_handler.clear_mg_space();
465 
466  const ::FiniteElement<2, spacedim> &fe = dof_handler.get_fe();
467  const ::Triangulation<2, spacedim> &tria =
468  dof_handler.get_triangulation();
469  const unsigned int n_levels = tria.n_levels();
470 
471  for (unsigned int i = 0; i < n_levels; ++i)
472  {
473  dof_handler.mg_levels.emplace_back(
474  std_cxx14::make_unique<
476  dof_handler.mg_levels.back()->dof_object.dofs =
477  std::vector<types::global_dof_index>(tria.n_raw_quads(i) *
478  fe.dofs_per_quad,
480  }
481 
482  dof_handler.mg_faces = std_cxx14::make_unique<
484  dof_handler.mg_faces->lines.dofs = std::vector<types::global_dof_index>(
485  tria.n_raw_lines() * fe.dofs_per_line, numbers::invalid_dof_index);
486 
487  const unsigned int n_vertices = tria.n_vertices();
488 
489  dof_handler.mg_vertex_dofs.resize(n_vertices);
490 
491  std::vector<unsigned int> max_level(n_vertices, 0);
492  std::vector<unsigned int> min_level(n_vertices, n_levels);
493 
494  for (typename ::Triangulation<2, spacedim>::cell_iterator cell =
495  tria.begin();
496  cell != tria.end();
497  ++cell)
498  {
499  const unsigned int level = cell->level();
500 
501  for (const unsigned int vertex : GeometryInfo<2>::vertex_indices())
502  {
503  const unsigned int vertex_index = cell->vertex_index(vertex);
504 
505  if (min_level[vertex_index] > level)
506  min_level[vertex_index] = level;
507 
508  if (max_level[vertex_index] < level)
509  max_level[vertex_index] = level;
510  }
511  }
512 
513  for (unsigned int vertex = 0; vertex < n_vertices; ++vertex)
514  if (tria.vertex_used(vertex))
515  {
516  Assert(min_level[vertex] < n_levels, ExcInternalError());
517  Assert(max_level[vertex] >= min_level[vertex],
518  ExcInternalError());
519  dof_handler.mg_vertex_dofs[vertex].init(min_level[vertex],
520  max_level[vertex],
521  fe.dofs_per_vertex);
522  }
523 
524  else
525  {
526  Assert(min_level[vertex] == n_levels, ExcInternalError());
527  Assert(max_level[vertex] == 0, ExcInternalError());
528  dof_handler.mg_vertex_dofs[vertex].init(1, 0, 0);
529  }
530  }
531 
532  template <int spacedim>
533  static void reserve_space_mg(DoFHandler<3, spacedim> &dof_handler)
534  {
535  Assert(dof_handler.get_triangulation().n_levels() > 0,
536  ExcMessage("Invalid triangulation"));
537  dof_handler.clear_mg_space();
538 
539  const ::FiniteElement<3, spacedim> &fe = dof_handler.get_fe();
540  const ::Triangulation<3, spacedim> &tria =
541  dof_handler.get_triangulation();
542  const unsigned int n_levels = tria.n_levels();
543 
544  for (unsigned int i = 0; i < n_levels; ++i)
545  {
546  dof_handler.mg_levels.emplace_back(
547  std_cxx14::make_unique<
549  dof_handler.mg_levels.back()->dof_object.dofs =
550  std::vector<types::global_dof_index>(tria.n_raw_hexs(i) *
551  fe.dofs_per_hex,
553  }
554 
555  dof_handler.mg_faces = std_cxx14::make_unique<
557  dof_handler.mg_faces->lines.dofs = std::vector<types::global_dof_index>(
558  tria.n_raw_lines() * fe.dofs_per_line, numbers::invalid_dof_index);
559  dof_handler.mg_faces->quads.dofs = std::vector<types::global_dof_index>(
560  tria.n_raw_quads() * fe.dofs_per_quad, numbers::invalid_dof_index);
561 
562  const unsigned int n_vertices = tria.n_vertices();
563 
564  dof_handler.mg_vertex_dofs.resize(n_vertices);
565 
566  std::vector<unsigned int> max_level(n_vertices, 0);
567  std::vector<unsigned int> min_level(n_vertices, n_levels);
568 
569  for (typename ::Triangulation<3, spacedim>::cell_iterator cell =
570  tria.begin();
571  cell != tria.end();
572  ++cell)
573  {
574  const unsigned int level = cell->level();
575 
576  for (const unsigned int vertex : GeometryInfo<3>::vertex_indices())
577  {
578  const unsigned int vertex_index = cell->vertex_index(vertex);
579 
580  if (min_level[vertex_index] > level)
581  min_level[vertex_index] = level;
582 
583  if (max_level[vertex_index] < level)
584  max_level[vertex_index] = level;
585  }
586  }
587 
588  for (unsigned int vertex = 0; vertex < n_vertices; ++vertex)
589  if (tria.vertex_used(vertex))
590  {
591  Assert(min_level[vertex] < n_levels, ExcInternalError());
592  Assert(max_level[vertex] >= min_level[vertex],
593  ExcInternalError());
594  dof_handler.mg_vertex_dofs[vertex].init(min_level[vertex],
595  max_level[vertex],
596  fe.dofs_per_vertex);
597  }
598 
599  else
600  {
601  Assert(min_level[vertex] == n_levels, ExcInternalError());
602  Assert(max_level[vertex] == 0, ExcInternalError());
603  dof_handler.mg_vertex_dofs[vertex].init(1, 0, 0);
604  }
605  }
606 
607 
608 
609  template <int spacedim>
612  const DoFHandler<1, spacedim> &dof_handler,
614  &mg_level,
616  &,
617  const unsigned int obj_index,
618  const unsigned int fe_index,
619  const unsigned int local_index,
620  const std::integral_constant<int, 1>)
621  {
622  return mg_level->dof_object.get_dof_index(dof_handler,
623  obj_index,
624  fe_index,
625  local_index);
626  }
627 
628  template <int spacedim>
631  const DoFHandler<2, spacedim> &dof_handler,
633  &,
635  & mg_faces,
636  const unsigned int obj_index,
637  const unsigned int fe_index,
638  const unsigned int local_index,
639  const std::integral_constant<int, 1>)
640  {
641  return mg_faces->lines.get_dof_index(dof_handler,
642  obj_index,
643  fe_index,
644  local_index);
645  }
646 
647  template <int spacedim>
650  const DoFHandler<2, spacedim> &dof_handler,
652  &mg_level,
654  &,
655  const unsigned int obj_index,
656  const unsigned int fe_index,
657  const unsigned int local_index,
658  const std::integral_constant<int, 2>)
659  {
660  return mg_level->dof_object.get_dof_index(dof_handler,
661  obj_index,
662  fe_index,
663  local_index);
664  }
665 
666  template <int spacedim>
669  const DoFHandler<3, spacedim> &dof_handler,
671  &,
673  & mg_faces,
674  const unsigned int obj_index,
675  const unsigned int fe_index,
676  const unsigned int local_index,
677  const std::integral_constant<int, 1>)
678  {
679  return mg_faces->lines.get_dof_index(dof_handler,
680  obj_index,
681  fe_index,
682  local_index);
683  }
684 
685  template <int spacedim>
688  const DoFHandler<3, spacedim> &dof_handler,
690  &,
692  & mg_faces,
693  const unsigned int obj_index,
694  const unsigned int fe_index,
695  const unsigned int local_index,
696  const std::integral_constant<int, 2>)
697  {
698  return mg_faces->quads.get_dof_index(dof_handler,
699  obj_index,
700  fe_index,
701  local_index);
702  }
703 
704  template <int spacedim>
707  const DoFHandler<3, spacedim> &dof_handler,
709  &mg_level,
711  &,
712  const unsigned int obj_index,
713  const unsigned int fe_index,
714  const unsigned int local_index,
715  const std::integral_constant<int, 3>)
716  {
717  return mg_level->dof_object.get_dof_index(dof_handler,
718  obj_index,
719  fe_index,
720  local_index);
721  }
722 
723  template <int spacedim>
724  static void
726  const DoFHandler<1, spacedim> &dof_handler,
728  &mg_level,
730  &,
731  const unsigned int obj_index,
732  const unsigned int fe_index,
733  const unsigned int local_index,
735  const std::integral_constant<int, 1>)
736  {
737  mg_level->dof_object.set_dof_index(
738  dof_handler, obj_index, fe_index, local_index, global_index);
739  }
740 
741  template <int spacedim>
742  static void
744  const DoFHandler<2, spacedim> &dof_handler,
746  &,
748  & mg_faces,
749  const unsigned int obj_index,
750  const unsigned int fe_index,
751  const unsigned int local_index,
753  const std::integral_constant<int, 1>)
754  {
755  mg_faces->lines.set_dof_index(
756  dof_handler, obj_index, fe_index, local_index, global_index);
757  }
758 
759  template <int spacedim>
760  static void
762  const DoFHandler<2, spacedim> &dof_handler,
764  &mg_level,
766  &,
767  const unsigned int obj_index,
768  const unsigned int fe_index,
769  const unsigned int local_index,
771  const std::integral_constant<int, 2>)
772  {
773  mg_level->dof_object.set_dof_index(
774  dof_handler, obj_index, fe_index, local_index, global_index);
775  }
776 
777  template <int spacedim>
778  static void
780  const DoFHandler<3, spacedim> &dof_handler,
782  &,
784  & mg_faces,
785  const unsigned int obj_index,
786  const unsigned int fe_index,
787  const unsigned int local_index,
789  const std::integral_constant<int, 1>)
790  {
791  mg_faces->lines.set_dof_index(
792  dof_handler, obj_index, fe_index, local_index, global_index);
793  }
794 
795  template <int spacedim>
796  static void
798  const DoFHandler<3, spacedim> &dof_handler,
800  &,
802  & mg_faces,
803  const unsigned int obj_index,
804  const unsigned int fe_index,
805  const unsigned int local_index,
807  const std::integral_constant<int, 2>)
808  {
809  mg_faces->quads.set_dof_index(
810  dof_handler, obj_index, fe_index, local_index, global_index);
811  }
812 
813  template <int spacedim>
814  static void
816  const DoFHandler<3, spacedim> &dof_handler,
818  &mg_level,
820  &,
821  const unsigned int obj_index,
822  const unsigned int fe_index,
823  const unsigned int local_index,
825  const std::integral_constant<int, 3>)
826  {
827  mg_level->dof_object.set_dof_index(
828  dof_handler, obj_index, fe_index, local_index, global_index);
829  }
830  };
831  } // namespace DoFHandlerImplementation
832 } // namespace internal
833 
834 
835 
836 template <int dim, int spacedim>
838  : tria(&tria, typeid(*this).name())
839  , faces(nullptr)
840  , mg_faces(nullptr)
841 {
842  // decide whether we need a sequential or a parallel distributed policy
843  if (dynamic_cast<const parallel::shared::Triangulation<dim, spacedim> *>(
844  &tria) != nullptr)
845  policy =
846  std_cxx14::make_unique<internal::DoFHandlerImplementation::Policy::
847  ParallelShared<DoFHandler<dim, spacedim>>>(
848  *this);
849  else if (dynamic_cast<
851  &tria) == nullptr)
852  policy =
853  std_cxx14::make_unique<internal::DoFHandlerImplementation::Policy::
854  Sequential<DoFHandler<dim, spacedim>>>(*this);
855  else
856  policy =
857  std_cxx14::make_unique<internal::DoFHandlerImplementation::Policy::
858  ParallelDistributed<DoFHandler<dim, spacedim>>>(
859  *this);
860 }
861 
862 
863 template <int dim, int spacedim>
865  : tria(nullptr, typeid(*this).name())
866 {}
867 
868 
869 template <int dim, int spacedim>
871 {
872  // release allocated memory
873  // virtual functions called in constructors and destructors never use the
874  // override in a derived class
875  // for clarity be explicit on which function is called
877 
878  // also release the policy. this needs to happen before the
879  // current object disappears because the policy objects
880  // store references to the DoFhandler object they work on
881  policy.reset();
882 }
883 
884 
885 template <int dim, int spacedim>
886 void
889 {
890  tria = &t;
891  faces = nullptr;
892  number_cache.n_global_dofs = 0;
893 
894  // decide whether we need a sequential or a parallel distributed policy
895  if (dynamic_cast<const parallel::shared::Triangulation<dim, spacedim> *>(
896  &t) != nullptr)
897  policy =
898  std_cxx14::make_unique<internal::DoFHandlerImplementation::Policy::
899  ParallelShared<DoFHandler<dim, spacedim>>>(
900  *this);
901  else if (dynamic_cast<
903  &t) != nullptr)
904  policy =
905  std_cxx14::make_unique<internal::DoFHandlerImplementation::Policy::
906  ParallelDistributed<DoFHandler<dim, spacedim>>>(
907  *this);
908  else
909  policy =
910  std_cxx14::make_unique<internal::DoFHandlerImplementation::Policy::
911  Sequential<DoFHandler<dim, spacedim>>>(*this);
912 
913  distribute_dofs(fe);
914 }
915 
916 
917 
918 /*------------------------ Cell iterator functions ------------------------*/
919 
920 template <int dim, int spacedim>
922 DoFHandler<dim, spacedim>::begin(const unsigned int level) const
923 {
925  this->get_triangulation().begin(level);
926  if (cell == this->get_triangulation().end(level))
927  return end(level);
928  return cell_iterator(*cell, this);
929 }
930 
931 
932 
933 template <int dim, int spacedim>
936 {
937  // level is checked in begin
939  if (i.state() != IteratorState::valid)
940  return i;
941  while (i->has_children())
942  if ((++i).state() != IteratorState::valid)
943  return i;
944  return i;
945 }
946 
947 
948 
949 template <int dim, int spacedim>
952 {
953  return cell_iterator(&this->get_triangulation(), -1, -1, this);
954 }
955 
956 
957 template <int dim, int spacedim>
959 DoFHandler<dim, spacedim>::end(const unsigned int level) const
960 {
962  this->get_triangulation().end(level);
963  if (cell.state() != IteratorState::valid)
964  return end();
965  return cell_iterator(*cell, this);
966 }
967 
968 
969 template <int dim, int spacedim>
972 {
974  this->get_triangulation().end_active(level);
975  if (cell.state() != IteratorState::valid)
976  return active_cell_iterator(end());
977  return active_cell_iterator(*cell, this);
978 }
979 
980 
981 
982 template <int dim, int spacedim>
985 {
986  // Assert(this->has_level_dofs(), ExcMessage("You can only iterate over mg "
987  // "levels if mg dofs got distributed."));
989  this->get_triangulation().begin(level);
990  if (cell == this->get_triangulation().end(level))
991  return end_mg(level);
992  return level_cell_iterator(*cell, this);
993 }
994 
995 
996 template <int dim, int spacedim>
998 DoFHandler<dim, spacedim>::end_mg(const unsigned int level) const
999 {
1000  // Assert(this->has_level_dofs(), ExcMessage("You can only iterate over mg "
1001  // "levels if mg dofs got distributed."));
1003  this->get_triangulation().end(level);
1004  if (cell.state() != IteratorState::valid)
1005  return end();
1006  return level_cell_iterator(*cell, this);
1007 }
1008 
1009 
1010 template <int dim, int spacedim>
1013 {
1014  return level_cell_iterator(&this->get_triangulation(), -1, -1, this);
1015 }
1016 
1017 
1018 
1019 template <int dim, int spacedim>
1022 {
1024  begin(), end());
1025 }
1026 
1027 
1028 template <int dim, int spacedim>
1031 {
1032  return IteratorRange<
1033  typename DoFHandler<dim, spacedim>::active_cell_iterator>(begin_active(),
1034  end());
1035 }
1036 
1037 
1038 
1039 template <int dim, int spacedim>
1042 {
1044  begin_mg(), end_mg());
1045 }
1046 
1047 
1048 
1049 template <int dim, int spacedim>
1052  const unsigned int level) const
1053 {
1055  begin(level), end(level));
1056 }
1057 
1058 
1059 
1060 template <int dim, int spacedim>
1063  const unsigned int level) const
1064 {
1065  return IteratorRange<
1067  begin_active(level), end_active(level));
1068 }
1069 
1070 
1071 
1072 template <int dim, int spacedim>
1075  const unsigned int level) const
1076 {
1078  begin_mg(level), end_mg(level));
1079 }
1080 
1081 
1082 
1083 //---------------------------------------------------------------------------
1084 
1085 
1086 
1087 template <int dim, int spacedim>
1090 {
1091  const FiniteElement<dim, spacedim> &fe = get_fe();
1092  const unsigned int dofs_per_face = fe.n_dofs_per_face();
1093  std::vector<::types::global_dof_index> dofs_on_face(dofs_per_face);
1094 
1095  const IndexSet &owned_dofs = locally_owned_dofs();
1096 
1097  std::unordered_set<unsigned int> boundary_dof_indices;
1098  // Loop over all faces of all cells and see whether they are at a boundary.
1099  // note (i) that we visit interior faces twice (which we don't care about) but
1100  // exterior faces only once as is appropriate, and (ii) that we need not take
1101  // special care of single lines (using @p{cell->has_boundary_lines}), since we
1102  // do not support boundaries of dimension dim-2, and so every boundary line is
1103  // also part of a boundary face.
1104  for (const auto &cell : this->active_cell_iterators())
1105  {
1106  if (cell->is_locally_owned() && cell->at_boundary())
1107  {
1108  for (const unsigned int iface :
1110  {
1111  const auto face = cell->face(iface);
1112  if (face->at_boundary())
1113  {
1114  cell->face(iface)->get_dof_indices(dofs_on_face);
1115  for (unsigned int idof_face = 0; idof_face < dofs_per_face;
1116  ++idof_face)
1117  {
1118  const unsigned int global_idof_index =
1119  dofs_on_face[idof_face];
1120 
1121  // If dof is not already in our hashmap, then insert it
1122  // into our set of surface indices. However, do not add if
1123  // it is not locally owned, it will get picked up by
1124  // another processor
1125  if (owned_dofs.is_element(global_idof_index))
1126  {
1127  boundary_dof_indices.insert(global_idof_index);
1128  }
1129  }
1130  }
1131  }
1132  }
1133  }
1134  return boundary_dof_indices.size();
1135 }
1136 
1137 
1138 
1139 template <int dim, int spacedim>
1142  const std::set<types::boundary_id> &boundary_ids) const
1143 {
1144  Assert(boundary_ids.find(numbers::internal_face_boundary_id) ==
1145  boundary_ids.end(),
1147 
1148  const FiniteElement<dim, spacedim> &fe = get_fe();
1149  const unsigned int dofs_per_face = fe.n_dofs_per_face();
1150  std::vector<::types::global_dof_index> dofs_on_face(dofs_per_face);
1151 
1152  const IndexSet &owned_dofs = locally_owned_dofs();
1153 
1154  std::unordered_set<unsigned int> boundary_dof_indices;
1155  // Loop over all faces of all cells and see whether they are at a boundary.
1156  // note (i) that we visit interior faces twice (which we don't care about) but
1157  // exterior faces only once as is appropriate, and (ii) that we need not take
1158  // special care of single lines (using @p{cell->has_boundary_lines}), since we
1159  // do not support boundaries of dimension dim-2, and so every boundary line is
1160  // also part of a boundary face.
1161  for (const auto &cell : this->active_cell_iterators())
1162  {
1163  if (cell->is_locally_owned() && cell->at_boundary())
1164  {
1165  for (const unsigned int iface :
1167  {
1168  const auto face = cell->face(iface);
1169  const unsigned int boundary_id = face->boundary_id();
1170 
1171  if (face->at_boundary() &&
1172  std::find(boundary_ids.begin(),
1173  boundary_ids.end(),
1174  boundary_id) != boundary_ids.end())
1175  {
1176  cell->face(iface)->get_dof_indices(dofs_on_face);
1177  for (unsigned int idof_face = 0; idof_face < dofs_per_face;
1178  ++idof_face)
1179  {
1180  const unsigned int global_idof_index =
1181  dofs_on_face[idof_face];
1182 
1183  // If dof is not already in our hashmap, then insert it
1184  // into our set of surface indices. However, do not add if
1185  // it is not locally owned, it will get picked up by
1186  // another processor
1187  if (owned_dofs.is_element(global_idof_index))
1188  {
1189  boundary_dof_indices.insert(global_idof_index);
1190  }
1191  }
1192  }
1193  }
1194  }
1195  }
1196  return boundary_dof_indices.size();
1197 }
1198 
1199 
1200 
1201 template <int dim, int spacedim>
1202 std::size_t
1204 {
1205  std::size_t mem =
1207  MemoryConsumption::memory_consumption(fe_collection) +
1208  MemoryConsumption::memory_consumption(block_info_object) +
1211  MemoryConsumption::memory_consumption(faces) + sizeof(number_cache) +
1212  MemoryConsumption::memory_consumption(mg_number_cache) +
1214  for (unsigned int i = 0; i < levels.size(); ++i)
1215  mem += MemoryConsumption::memory_consumption(*levels[i]);
1216 
1217  for (unsigned int level = 0; level < mg_levels.size(); ++level)
1218  mem += mg_levels[level]->memory_consumption();
1219 
1220  if (mg_faces != nullptr)
1221  mem += MemoryConsumption::memory_consumption(*mg_faces);
1222 
1223  for (unsigned int i = 0; i < mg_vertex_dofs.size(); ++i)
1224  mem += sizeof(MGVertexDoFs) + (1 + mg_vertex_dofs[i].get_finest_level() -
1225  mg_vertex_dofs[i].get_coarsest_level()) *
1226  sizeof(types::global_dof_index);
1227 
1228  return mem;
1229 }
1230 
1231 
1232 
1233 template <int dim, int spacedim>
1234 void
1236 {
1237  // Only recreate the FECollection if we don't already store
1238  // the exact same FiniteElement object.
1239  if (fe_collection.size() == 0 || fe_collection[0] != ff)
1240  fe_collection = hp::FECollection<dim, spacedim>(ff);
1241 }
1242 
1243 
1244 
1245 template <int dim, int spacedim>
1246 void
1248  const FiniteElement<dim, spacedim> &ff)
1249 {
1250  Assert(
1251  tria != nullptr,
1252  ExcMessage(
1253  "You need to set the Triangulation in the DoFHandler using initialize() or "
1254  "in the constructor before you can distribute DoFs."));
1255  Assert(tria->n_levels() > 0,
1256  ExcMessage("The Triangulation you are using is empty!"));
1257 
1258  // first, assign the finite_element
1259  set_fe(ff);
1260 
1261  // delete all levels and set them
1262  // up newly. note that we still
1263  // have to allocate space for all
1264  // degrees of freedom on this mesh
1265  // (including ghost and cells that
1266  // are entirely stored on different
1267  // processors), though we may not
1268  // assign numbers to some of them
1269  // (i.e. they will remain at
1270  // invalid_dof_index). We need to
1271  // allocate the space because we
1272  // will want to be able to query
1273  // the dof_indices on each cell,
1274  // and simply be told that we don't
1275  // know them on some cell (i.e. get
1276  // back invalid_dof_index)
1277  clear_space();
1279 
1280  // hand things off to the policy
1281  number_cache = policy->distribute_dofs();
1282 
1283  // initialize the block info object
1284  // only if this is a sequential
1285  // triangulation. it doesn't work
1286  // correctly yet if it is parallel
1288  *>(&*tria) == nullptr)
1289  block_info_object.initialize(*this, false, true);
1290 }
1291 
1292 
1293 
1294 template <int dim, int spacedim>
1295 void
1297 {
1298  Assert(
1299  levels.size() > 0,
1300  ExcMessage(
1301  "Distribute active DoFs using distribute_dofs() before calling distribute_mg_dofs()."));
1302 
1303  Assert(
1304  ((tria->get_mesh_smoothing() &
1307  ExcMessage(
1308  "The mesh smoothing requirement 'limit_level_difference_at_vertices' has to be set for using multigrid!"));
1309 
1310  clear_mg_space();
1311 
1313  mg_number_cache = policy->distribute_mg_dofs();
1314 
1315  // initialize the block info object
1316  // only if this is a sequential
1317  // triangulation. it doesn't work
1318  // correctly yet if it is parallel
1319  if (dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1320  &*tria) == nullptr)
1321  block_info_object.initialize(*this, true, false);
1322 }
1323 
1324 
1325 
1326 template <int dim, int spacedim>
1327 void
1329 {
1330  mg_levels.clear();
1331  mg_faces.reset();
1332 
1333  std::vector<MGVertexDoFs> tmp;
1334 
1335  std::swap(mg_vertex_dofs, tmp);
1336 
1337  mg_number_cache.clear();
1338 }
1339 
1340 
1341 template <int dim, int spacedim>
1342 void
1344 {
1345  block_info_object.initialize_local(*this);
1346 }
1347 
1348 
1349 
1350 template <int dim, int spacedim>
1351 void
1353 {
1354  // release memory
1355  clear_space();
1356  clear_mg_space();
1357 }
1358 
1359 
1360 
1361 template <int dim, int spacedim>
1362 void
1364  const std::vector<types::global_dof_index> &new_numbers)
1365 {
1366  Assert(levels.size() > 0,
1367  ExcMessage(
1368  "You need to distribute DoFs before you can renumber them."));
1369 
1370 #ifdef DEBUG
1371  if (dynamic_cast<const parallel::shared::Triangulation<dim, spacedim> *>(
1372  &*tria) != nullptr)
1373  {
1374  Assert(new_numbers.size() == n_dofs() ||
1375  new_numbers.size() == n_locally_owned_dofs(),
1376  ExcMessage("Incorrect size of the input array."));
1377  }
1378  else if (dynamic_cast<
1380  &*tria) != nullptr)
1381  {
1382  AssertDimension(new_numbers.size(), n_locally_owned_dofs());
1383  }
1384  else
1385  {
1386  AssertDimension(new_numbers.size(), n_dofs());
1387  }
1388 
1389  // assert that the new indices are
1390  // consecutively numbered if we are
1391  // working on a single
1392  // processor. this doesn't need to
1393  // hold in the case of a parallel
1394  // mesh since we map the interval
1395  // [0...n_dofs()) into itself but
1396  // only globally, not on each
1397  // processor
1398  if (n_locally_owned_dofs() == n_dofs())
1399  {
1400  std::vector<types::global_dof_index> tmp(new_numbers);
1401  std::sort(tmp.begin(), tmp.end());
1402  std::vector<types::global_dof_index>::const_iterator p = tmp.begin();
1404  for (; p != tmp.end(); ++p, ++i)
1405  Assert(*p == i, ExcNewNumbersNotConsecutive(i));
1406  }
1407  else
1408  for (const auto new_number : new_numbers)
1409  Assert(new_number < n_dofs(),
1410  ExcMessage(
1411  "New DoF index is not less than the total number of dofs."));
1412 #endif
1413 
1414  number_cache = policy->renumber_dofs(new_numbers);
1415 }
1416 
1417 
1418 template <int dim, int spacedim>
1419 void
1421  const unsigned int level,
1422  const std::vector<types::global_dof_index> &new_numbers)
1423 {
1424  Assert(
1425  mg_levels.size() > 0 && levels.size() > 0,
1426  ExcMessage(
1427  "You need to distribute active and level DoFs before you can renumber level DoFs."));
1428  AssertIndexRange(level, get_triangulation().n_global_levels());
1429  AssertDimension(new_numbers.size(),
1430  locally_owned_mg_dofs(level).n_elements());
1431 
1432 #ifdef DEBUG
1433  // assert that the new indices are consecutively numbered if we are working
1434  // on a single processor. this doesn't need to hold in the case of a
1435  // parallel mesh since we map the interval [0...n_dofs(level)) into itself
1436  // but only globally, not on each processor
1437  if (n_locally_owned_dofs() == n_dofs())
1438  {
1439  std::vector<types::global_dof_index> tmp(new_numbers);
1440  std::sort(tmp.begin(), tmp.end());
1441  std::vector<types::global_dof_index>::const_iterator p = tmp.begin();
1443  for (; p != tmp.end(); ++p, ++i)
1444  Assert(*p == i, ExcNewNumbersNotConsecutive(i));
1445  }
1446  else
1447  for (const auto new_number : new_numbers)
1448  Assert(new_number < n_dofs(level),
1449  ExcMessage(
1450  "New DoF index is not less than the total number of dofs."));
1451 #endif
1452 
1453  mg_number_cache[level] = policy->renumber_mg_dofs(level, new_numbers);
1454 }
1455 
1456 
1457 
1458 template <int dim, int spacedim>
1459 unsigned int
1461 {
1464 }
1465 
1466 
1467 
1468 template <int dim, int spacedim>
1469 unsigned int
1471 {
1472  switch (dim)
1473  {
1474  case 1:
1475  return get_fe().dofs_per_vertex;
1476  case 2:
1477  return (3 * get_fe().dofs_per_vertex + 2 * get_fe().dofs_per_line);
1478  case 3:
1479  // we need to take refinement of
1480  // one boundary face into
1481  // consideration here; in fact,
1482  // this function returns what
1483  // #max_coupling_between_dofs<2>
1484  // returns
1485  //
1486  // we assume here, that only four
1487  // faces meet at the boundary;
1488  // this assumption is not
1489  // justified and needs to be
1490  // fixed some time. fortunately,
1491  // omitting it for now does no
1492  // harm since the matrix will cry
1493  // foul if its requirements are
1494  // not satisfied
1495  return (19 * get_fe().dofs_per_vertex + 28 * get_fe().dofs_per_line +
1496  8 * get_fe().dofs_per_quad);
1497  default:
1498  Assert(false, ExcNotImplemented());
1500  }
1501 }
1502 
1503 
1504 
1505 template <int dim, int spacedim>
1506 void
1508 {
1509  levels.clear();
1510  faces.reset();
1511 
1512  std::vector<types::global_dof_index> tmp;
1513  std::swap(vertex_dofs, tmp);
1514 
1515  number_cache.clear();
1516 }
1517 
1518 
1519 
1520 template <int dim, int spacedim>
1521 template <int structdim>
1523 DoFHandler<dim, spacedim>::get_dof_index(const unsigned int obj_level,
1524  const unsigned int obj_index,
1525  const unsigned int fe_index,
1526  const unsigned int local_index) const
1527 {
1529  *this,
1530  this->mg_levels[obj_level],
1531  this->mg_faces,
1532  obj_index,
1533  fe_index,
1534  local_index,
1535  std::integral_constant<int, structdim>());
1536 }
1537 
1538 
1539 
1540 template <int dim, int spacedim>
1541 template <int structdim>
1542 void
1544  const unsigned int obj_level,
1545  const unsigned int obj_index,
1546  const unsigned int fe_index,
1547  const unsigned int local_index,
1549 {
1551  *this,
1552  this->mg_levels[obj_level],
1553  this->mg_faces,
1554  obj_index,
1555  fe_index,
1556  local_index,
1557  global_index,
1558  std::integral_constant<int, structdim>());
1559 }
1560 
1561 
1562 
1563 template <int dim, int spacedim>
1565  : coarsest_level(numbers::invalid_unsigned_int)
1566  , finest_level(0)
1567 {}
1568 
1569 
1570 
1571 template <int dim, int spacedim>
1572 void
1574  const unsigned int cl,
1575  const unsigned int fl,
1576  const unsigned int dofs_per_vertex)
1577 {
1578  coarsest_level = cl;
1579  finest_level = fl;
1580 
1581  if (coarsest_level <= finest_level)
1582  {
1583  const unsigned int n_levels = finest_level - coarsest_level + 1;
1584  const unsigned int n_indices = n_levels * dofs_per_vertex;
1585 
1586  indices = std_cxx14::make_unique<types::global_dof_index[]>(n_indices);
1587  std::fill(indices.get(),
1588  indices.get() + n_indices,
1590  }
1591  else
1592  indices.reset();
1593 }
1594 
1595 
1596 
1597 template <int dim, int spacedim>
1598 unsigned int
1600 {
1601  return coarsest_level;
1602 }
1603 
1604 
1605 
1606 template <int dim, int spacedim>
1607 unsigned int
1609 {
1610  return finest_level;
1611 }
1612 
1613 
1614 /*-------------- Explicit Instantiations -------------------------------*/
1615 #include "dof_handler.inst"
1616 
1617 
internal::DoFHandlerImplementation::Implementation::set_dof_index
static void set_dof_index(const DoFHandler< 2, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< 2 >> &, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< 2 >> &mg_faces, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index, const std::integral_constant< int, 1 >)
Definition: dof_handler.cc:743
internal::DoFHandlerImplementation::Policy
Definition: dof_handler_policy.h:51
IndexSet
Definition: index_set.h:74
internal::DoFHandlerImplementation::Implementation::reserve_space_mg
static void reserve_space_mg(DoFHandler< 1, spacedim > &dof_handler)
Definition: dof_handler.cc:392
internal::DoFHandlerImplementation::Implementation::get_dof_index
static types::global_dof_index get_dof_index(const DoFHandler< 2, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< 2 >> &mg_level, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< 2 >> &, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const std::integral_constant< int, 2 >)
Definition: dof_handler.cc:649
numbers::invalid_dof_index
const types::global_dof_index invalid_dof_index
Definition: types.h:206
parallel::shared::Triangulation
Definition: shared_tria.h:103
DoFHandler::mg_cell_iterators_on_level
IteratorRange< level_cell_iterator > mg_cell_iterators_on_level(const unsigned int level) const
Definition: dof_handler.cc:1074
DoFHandler::faces
std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces< dim > > faces
Definition: dof_handler.h:1349
DoFHandler::cell_iterator
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:357
DoFHandler::active_cell_iterators_on_level
IteratorRange< active_cell_iterator > active_cell_iterators_on_level(const unsigned int level) const
Definition: dof_handler.cc:1062
parallel::DistributedTriangulationBase
Definition: tria_base.h:352
tria_accessor.h
DoFHandler::active_cell_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:329
DoFHandler::levels
std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel< dim > > > levels
Definition: dof_handler.h:1337
internal::DoFHandlerImplementation::DoFLevel
Definition: dof_levels.h:71
Utilities::int_to_string
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:474
DoFHandler::begin_mg
level_cell_iterator begin_mg(const unsigned int level=0) const
Definition: dof_handler.cc:984
DoFHandler::level_cell_iterator
typename LevelSelector::cell_iterator level_cell_iterator
Definition: dof_handler.h:393
internal::DoFHandlerImplementation::Implementation::get_dof_index
static types::global_dof_index get_dof_index(const DoFHandler< 2, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< 2 >> &, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< 2 >> &mg_faces, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const std::integral_constant< int, 1 >)
Definition: dof_handler.cc:630
DoFHandler::policy
std::unique_ptr<::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > > policy
Definition: dof_handler.h:1198
tria_levels.h
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
Triangulation
Definition: tria.h:1109
tria.h
DoFHandler::end_active
active_cell_iterator end_active(const unsigned int level) const
Definition: dof_handler.cc:971
DoFHandler::mg_faces
std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces< dim > > mg_faces
Definition: dof_handler.h:1352
memory_consumption.h
DoFHandler::MGVertexDoFs::get_finest_level
unsigned int get_finest_level() const
Definition: dof_handler.cc:1608
tria_iterator.h
GeometryInfo
Definition: geometry_info.h:1224
DoFHandler::initialize_local_block_info
void initialize_local_block_info()
Definition: dof_handler.cc:1343
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
internal::DoFHandlerImplementation::Implementation::set_dof_index
static void set_dof_index(const DoFHandler< 3, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< 3 >> &, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< 3 >> &mg_faces, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index, const std::integral_constant< int, 1 >)
Definition: dof_handler.cc:779
FiniteElementData::n_dofs_per_face
unsigned int n_dofs_per_face() const
IteratorState::valid
@ valid
Iterator points to a valid object.
Definition: tria_iterator_base.h:38
internal::DoFHandlerImplementation::Implementation::reserve_space
static void reserve_space(DoFHandler< 3, spacedim > &dof_handler)
Definition: dof_handler.cc:353
types::boundary_id
unsigned int boundary_id
Definition: types.h:129
DoFHandler::renumber_dofs
void renumber_dofs(const std::vector< types::global_dof_index > &new_numbers)
Definition: dof_handler.cc:1363
hp::DoFHandler
Definition: dof_handler.h:203
DoFHandler::DoFHandler
DoFHandler()
Definition: dof_handler.cc:864
DoFHandler::clear_space
void clear_space()
Definition: dof_handler.cc:1507
internal::DoFHandlerImplementation::Implementation::set_dof_index
static void set_dof_index(const DoFHandler< 2, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< 2 >> &mg_level, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< 2 >> &, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index, const std::integral_constant< int, 2 >)
Definition: dof_handler.cc:761
DoFHandler
Definition: dof_handler.h:205
DoFHandler::MGVertexDoFs::get_coarsest_level
unsigned int get_coarsest_level() const
Definition: dof_handler.cc:1599
DoFHandler::max_couplings_between_boundary_dofs
unsigned int max_couplings_between_boundary_dofs() const
Definition: dof_handler.cc:1470
DoFHandler::distribute_dofs
virtual void distribute_dofs(const FiniteElement< dim, spacedim > &fe)
Definition: dof_handler.cc:1247
internal::DoFHandlerImplementation::Implementation::max_couplings_between_dofs
static unsigned int max_couplings_between_dofs(const DoFHandler< 2, spacedim > &dof_handler)
Definition: dof_handler.cc:140
internal::DoFHandlerImplementation::Implementation::get_dof_index
static types::global_dof_index get_dof_index(const DoFHandler< 3, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< 3 >> &mg_level, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< 3 >> &, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const std::integral_constant< int, 3 >)
Definition: dof_handler.cc:706
level
unsigned int level
Definition: grid_out.cc:4355
internal::DoFHandlerImplementation
Definition: dof_faces.h:38
DoFHandler::active_cell_iterators
IteratorRange< active_cell_iterator > active_cell_iterators() const
Definition: dof_handler.cc:1030
internal::DoFHandlerImplementation::DoFFaces< 1 >
Definition: dof_faces.h:87
DoFHandler::begin
cell_iterator begin(const unsigned int level=0) const
Definition: dof_handler.cc:922
DoFHandler::mg_cell_iterators
IteratorRange< level_cell_iterator > mg_cell_iterators() const
Definition: dof_handler.cc:1041
DoFHandler::~DoFHandler
virtual ~DoFHandler() override
Definition: dof_handler.cc:870
DoFHandler::get_dof_index
types::global_dof_index get_dof_index(const unsigned int obj_level, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index) const
Definition: dof_handler.cc:1523
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
TrilinosWrappers::internal::begin
VectorType::value_type * begin(VectorType &V)
Definition: trilinos_sparse_matrix.cc:51
geometry_info.h
IteratorRange
Definition: iterator_range.h:129
DoFHandler::mg_vertex_dofs
std::vector< MGVertexDoFs > mg_vertex_dofs
Definition: dof_handler.h:1329
DoFHandler::get_triangulation
const Triangulation< dim, spacedim > & get_triangulation() const
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
MemoryConsumption::memory_consumption
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
Definition: memory_consumption.h:268
internal::DoFHandlerImplementation::Implementation::get_dof_index
static types::global_dof_index get_dof_index(const DoFHandler< 1, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< 1 >> &mg_level, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< 1 >> &, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const std::integral_constant< int, 1 >)
Definition: dof_handler.cc:611
DoFHandler::vertex_dofs
std::vector< types::global_dof_index > vertex_dofs
Definition: dof_handler.h:1323
fe.h
shared_tria.h
IndexSet::is_element
bool is_element(const size_type index) const
Definition: index_set.h:1766
internal::DoFHandlerImplementation::Implementation::get_dof_index
static types::global_dof_index get_dof_index(const DoFHandler< 3, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< 3 >> &, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< 3 >> &mg_faces, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const std::integral_constant< int, 2 >)
Definition: dof_handler.cc:687
TrilinosWrappers::internal::end
VectorType::value_type * end(VectorType &V)
Definition: trilinos_sparse_matrix.cc:65
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
MemorySpace::swap
void swap(MemorySpaceData< Number, MemorySpace > &, MemorySpaceData< Number, MemorySpace > &)
Definition: memory_space.h:103
DoFHandler::end_mg
level_cell_iterator end_mg() const
Definition: dof_handler.cc:1012
numbers::internal_face_boundary_id
const types::boundary_id internal_face_boundary_id
Definition: types.h:250
numbers
Definition: numbers.h:207
internal::DoFHandlerImplementation::Implementation::reserve_space_mg
static void reserve_space_mg(DoFHandler< 3, spacedim > &dof_handler)
Definition: dof_handler.cc:533
DoFHandler::cell_iterators
IteratorRange< cell_iterator > cell_iterators() const
Definition: dof_handler.cc:1021
FiniteElement
Definition: fe.h:648
internal::DoFHandlerImplementation::Implementation::get_dof_index
static types::global_dof_index get_dof_index(const DoFHandler< 3, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< 3 >> &, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< 3 >> &mg_faces, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const std::integral_constant< int, 1 >)
Definition: dof_handler.cc:668
DoFHandler::end
cell_iterator end() const
Definition: dof_handler.cc:951
DoFHandler::max_couplings_between_dofs
unsigned int max_couplings_between_dofs() const
Definition: dof_handler.cc:1460
DoFHandler::cell_iterators_on_level
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
Definition: dof_handler.cc:1051
TriaRawIterator::state
IteratorState::IteratorStates state() const
Definition: tria_iterator.h:1051
internal::dummy
const types::global_dof_index * dummy()
Definition: dof_handler.cc:59
unsigned int
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
tria.h
DoFHandler::MGVertexDoFs::init
void init(const unsigned int coarsest_level, const unsigned int finest_level, const unsigned int dofs_per_vertex)
Definition: dof_handler.cc:1573
DoFTools::ExcInvalidBoundaryIndicator
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
internal::DoFHandlerImplementation::DoFFaces< 3 >
Definition: dof_faces.h:144
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
DoFHandler::memory_consumption
virtual std::size_t memory_consumption() const
Definition: dof_handler.cc:1203
DoFHandler::clear
virtual void clear()
Definition: dof_handler.cc:1352
dof_handler.h
dof_accessor.h
internal::DoFHandlerImplementation::DoFFaces< 2 >
Definition: dof_faces.h:113
dof_faces.h
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
internal::policy_to_string
std::string policy_to_string(const ::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > &policy)
Definition: dof_handler.cc:71
Utilities::MPI::min
T min(const T &t, const MPI_Comm &mpi_communicator)
internal::DoFHandlerImplementation::Implementation::set_dof_index
static void set_dof_index(const DoFHandler< 3, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< 3 >> &, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< 3 >> &mg_faces, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index, const std::integral_constant< int, 2 >)
Definition: dof_handler.cc:797
internal::DoFHandlerImplementation::Implementation
Definition: dof_handler.cc:120
internal::DoFHandlerImplementation::Implementation::reserve_space_mg
static void reserve_space_mg(DoFHandler< 2, spacedim > &dof_handler)
Definition: dof_handler.cc:460
parallel::TriangulationBase
Definition: tria_base.h:78
internal::DoFHandlerImplementation::Implementation::max_couplings_between_dofs
static unsigned int max_couplings_between_dofs(const DoFHandler< 1, spacedim > &dof_handler)
Definition: dof_handler.cc:128
DoFHandler::set_dof_index
void set_dof_index(const unsigned int obj_level, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index) const
Definition: dof_handler.cc:1543
dof_handler_policy.h
DoFHandler::distribute_mg_dofs
virtual void distribute_mg_dofs()
Definition: dof_handler.cc:1296
internal
Definition: aligned_vector.h:369
memory.h
internal::DoFHandlerImplementation::Implementation::set_dof_index
static void set_dof_index(const DoFHandler< 3, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< 3 >> &mg_level, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< 3 >> &, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index, const std::integral_constant< int, 3 >)
Definition: dof_handler.cc:815
DoFHandler::n_boundary_dofs
types::global_dof_index n_boundary_dofs() const
Definition: dof_handler.cc:1089
DoFHandler::clear_mg_space
void clear_mg_space()
Definition: dof_handler.cc:1328
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
dof_levels.h
internal::DoFHandlerImplementation::Implementation::reserve_space
static void reserve_space(DoFHandler< 2, spacedim > &dof_handler)
Definition: dof_handler.cc:317
internal::DoFHandlerImplementation::Policy::ParallelShared
Definition: dof_handler_policy.h:161
internal::DoFHandlerImplementation::Implementation::set_dof_index
static void set_dof_index(const DoFHandler< 1, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< 1 >> &mg_level, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< 1 >> &, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index, const std::integral_constant< int, 1 >)
Definition: dof_handler.cc:725
DoFHandler::tria
SmartPointer< const Triangulation< dim, spacedim >, DoFHandler< dim, spacedim > > tria
Definition: dof_handler.h:1183
internal::DoFHandlerImplementation::Implementation::reserve_space
static void reserve_space(DoFHandler< 1, spacedim > &dof_handler)
Definition: dof_handler.cc:292
DoFHandler::begin_active
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: dof_handler.cc:935
internal::DoFHandlerImplementation::Implementation::max_couplings_between_dofs
static unsigned int max_couplings_between_dofs(const DoFHandler< 3, spacedim > &dof_handler)
Definition: dof_handler.cc:247
TriaIterator
Definition: tria_iterator.h:578
TrilinosWrappers::global_index
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
Definition: trilinos_index_access.h:82
DoFHandler::initialize
void initialize(const Triangulation< dim, spacedim > &tria, const FiniteElement< dim, spacedim > &fe)
Definition: dof_handler.cc:887
DoFHandler::mg_levels
std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel< dim > > > mg_levels
Definition: dof_handler.h:1341
DoFHandler::get_fe
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
DoFHandler::MGVertexDoFs::MGVertexDoFs
MGVertexDoFs()
Definition: dof_handler.cc:1564
DoFHandler::set_fe
virtual void set_fe(const FiniteElement< dim, spacedim > &fe)
Definition: dof_handler.cc:1235
hp::FECollection
Definition: fe_collection.h:54
DoFHandler::n_dofs
types::global_dof_index n_dofs() const