Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
mg_tools.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1999 - 2022 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
20
24
26#include <deal.II/fe/fe.h>
27
29#include <deal.II/grid/tria.h>
31
38
42
44
45#include <algorithm>
46#include <numeric>
47#include <vector>
48
50
51
52namespace MGTools
53{
54 // specializations for 1D
55 template <>
56 void
58 const unsigned int,
59 std::vector<unsigned int> &,
61 {
62 Assert(false, ExcNotImplemented());
63 }
64
65
66
67 template <>
68 void
70 const unsigned int,
71 std::vector<unsigned int> &,
74 {
75 Assert(false, ExcNotImplemented());
76 }
77
78
79
80 template <>
81 void
83 const unsigned int,
84 std::vector<unsigned int> &,
86 {
87 Assert(false, ExcNotImplemented());
88 }
89
90
91 template <>
92 void
94 const unsigned int,
95 std::vector<unsigned int> &,
98 {
99 Assert(false, ExcNotImplemented());
100 }
101
102
103
104 // Template for 2D and 3D. For 1D see specialization above
105 template <int dim, int spacedim>
106 void
108 const unsigned int level,
109 std::vector<unsigned int> & row_lengths,
110 const DoFTools::Coupling flux_coupling)
111 {
112 Assert(row_lengths.size() == dofs.n_dofs(),
113 ExcDimensionMismatch(row_lengths.size(), dofs.n_dofs()));
114
115 // Function starts here by
116 // resetting the counters.
117 std::fill(row_lengths.begin(), row_lengths.end(), 0);
118 // We need the user flags, so we
119 // save them for later restoration
120 std::vector<bool> old_flags;
121 // We need a non-constant
122 // triangulation for the user
123 // flags. Since we restore them in
124 // the end, this cast is safe.
125 Triangulation<dim, spacedim> &user_flags_triangulation =
127 user_flags_triangulation.save_user_flags(old_flags);
128 user_flags_triangulation.clear_user_flags();
129
130 std::vector<types::global_dof_index> cell_indices;
131 std::vector<types::global_dof_index> neighbor_indices;
132
133 // We loop over cells and go from
134 // cells to lower dimensional
135 // objects. This is the only way to
136 // cope with the fact, that an
137 // unknown number of cells may
138 // share an object of dimension
139 // smaller than dim-1.
140 for (const auto &cell : dofs.cell_iterators_on_level(level))
141 {
142 const FiniteElement<dim> &fe = cell->get_fe();
143 cell_indices.resize(fe.n_dofs_per_cell());
144 cell->get_mg_dof_indices(cell_indices);
145 unsigned int i = 0;
146 // First, dofs on
147 // vertices. We assume that
148 // each vertex dof couples
149 // with all dofs on
150 // adjacent grid cells.
151
152 // Adding all dofs of the cells
153 // will add dofs of the faces
154 // of the cell adjacent to the
155 // vertex twice. Therefore, we
156 // subtract these here and add
157 // them in a loop over the
158 // faces below.
159
160 // in 1d, faces and vertices
161 // are identical. Nevertheless,
162 // this will only work if
163 // dofs_per_face is zero and
164 // n_dofs_per_vertex() is
165 // arbitrary, not the other way
166 // round.
167 // TODO: This assumes that the dofs per face on all faces coincide!
168 const unsigned int face_no = 0;
169
170 Assert(fe.reference_cell() == ReferenceCells::get_hypercube<dim>(),
172
173 unsigned int increment =
174 fe.n_dofs_per_cell() - dim * fe.n_dofs_per_face(face_no);
175 while (i < fe.get_first_line_index())
176 row_lengths[cell_indices[i++]] += increment;
177 // From now on, if an object is
178 // a cell, its dofs only couple
179 // inside the cell. Since the
180 // faces are handled below, we
181 // have to subtract ALL faces
182 // in this case.
183
184 // In all other cases we
185 // subtract adjacent faces to be
186 // added in the loop below.
187 increment =
188 (dim > 1) ?
189 fe.n_dofs_per_cell() - (dim - 1) * fe.n_dofs_per_face(face_no) :
190 fe.n_dofs_per_cell() -
192 while (i < fe.get_first_quad_index(face_no))
193 row_lengths[cell_indices[i++]] += increment;
194
195 // Now quads in 2D and 3D
196 increment =
197 (dim > 2) ?
198 fe.n_dofs_per_cell() - (dim - 2) * fe.n_dofs_per_face(face_no) :
199 fe.n_dofs_per_cell() -
201 while (i < fe.get_first_hex_index())
202 row_lengths[cell_indices[i++]] += increment;
203 // Finally, cells in 3D
205 fe.n_dofs_per_face(face_no);
206 while (i < fe.n_dofs_per_cell())
207 row_lengths[cell_indices[i++]] += increment;
208
209 // At this point, we have
210 // counted all dofs
211 // contributiong from cells
212 // coupled topologically to the
213 // adjacent cells, but we
214 // subtracted some faces.
215
216 // Now, let's go by the faces
217 // and add the missing
218 // contribution as well as the
219 // flux contributions.
220 for (const unsigned int iface : GeometryInfo<dim>::face_indices())
221 {
222 bool level_boundary = cell->at_boundary(iface);
224 if (!level_boundary)
225 {
226 neighbor = cell->neighbor(iface);
227 if (static_cast<unsigned int>(neighbor->level()) != level)
228 level_boundary = true;
229 }
230
231 if (level_boundary)
232 {
233 for (unsigned int local_dof = 0;
234 local_dof < fe.n_dofs_per_cell();
235 ++local_dof)
236 row_lengths[cell_indices[local_dof]] +=
237 fe.n_dofs_per_face(face_no);
238 continue;
239 }
240
241 const FiniteElement<dim> &nfe = neighbor->get_fe();
243 cell->face(iface);
244
245 // Flux couplings are
246 // computed from both sides
247 // for simplicity.
248
249 // The dofs on the common face
250 // will be handled below,
251 // therefore, we subtract them
252 // here.
253 if (flux_coupling != DoFTools::none)
254 {
255 const unsigned int dof_increment =
256 nfe.n_dofs_per_cell() - nfe.n_dofs_per_face(face_no);
257 for (unsigned int local_dof = 0;
258 local_dof < fe.n_dofs_per_cell();
259 ++local_dof)
260 row_lengths[cell_indices[local_dof]] += dof_increment;
261 }
262
263 // Do this only once per
264 // face.
265 if (face->user_flag_set())
266 continue;
267 face->set_user_flag();
268 // At this point, we assume
269 // that each cell added its
270 // dofs minus the face to
271 // the couplings of the
272 // face dofs. Since we
273 // subtracted two faces, we
274 // have to re-add one.
275
276 // If one side of the face
277 // is refined, all the fine
278 // face dofs couple with
279 // the coarse one.
280 neighbor_indices.resize(nfe.n_dofs_per_cell());
281 neighbor->get_mg_dof_indices(neighbor_indices);
282 for (unsigned int local_dof = 0; local_dof < fe.n_dofs_per_cell();
283 ++local_dof)
284 row_lengths[cell_indices[local_dof]] +=
285 nfe.n_dofs_per_face(face_no);
286 for (unsigned int local_dof = 0; local_dof < nfe.n_dofs_per_cell();
287 ++local_dof)
288 row_lengths[neighbor_indices[local_dof]] +=
289 fe.n_dofs_per_face(face_no);
290 }
291 }
292 user_flags_triangulation.load_user_flags(old_flags);
293 }
294
295
296 // This is the template for 2D and 3D. See version for 1D above
297 template <int dim, int spacedim>
298 void
300 const unsigned int level,
301 std::vector<unsigned int> & row_lengths,
302 const Table<2, DoFTools::Coupling> &couplings,
303 const Table<2, DoFTools::Coupling> &flux_couplings)
304 {
305 Assert(row_lengths.size() == dofs.n_dofs(),
306 ExcDimensionMismatch(row_lengths.size(), dofs.n_dofs()));
307
308 // Function starts here by
309 // resetting the counters.
310 std::fill(row_lengths.begin(), row_lengths.end(), 0);
311 // We need the user flags, so we
312 // save them for later restoration
313 std::vector<bool> old_flags;
314 // We need a non-constant
315 // triangulation for the user
316 // flags. Since we restore them in
317 // the end, this cast is safe.
318 Triangulation<dim, spacedim> &user_flags_triangulation =
320 user_flags_triangulation.save_user_flags(old_flags);
321 user_flags_triangulation.clear_user_flags();
322
323 std::vector<types::global_dof_index> cell_indices;
324 std::vector<types::global_dof_index> neighbor_indices;
325
326 // We have to translate the
327 // couplings from components to
328 // blocks, so this works for
329 // nonprimitive elements as well.
330 std::vector<Table<2, DoFTools::Coupling>> couple_cell;
331 std::vector<Table<2, DoFTools::Coupling>> couple_face;
332 DoFTools::convert_couplings_to_blocks(dofs, couplings, couple_cell);
333 DoFTools::convert_couplings_to_blocks(dofs, flux_couplings, couple_face);
334
335 // We loop over cells and go from
336 // cells to lower dimensional
337 // objects. This is the only way to
338 // cope with the fact, that an
339 // unknown number of cells may
340 // share an object of dimension
341 // smaller than dim-1.
342 for (const auto &cell : dofs.cell_iterators_on_level(level))
343 {
344 const FiniteElement<dim> &fe = cell->get_fe();
345 const unsigned int fe_index = cell->active_fe_index();
346
347
348 // TODO: This assumes that the dofs per face on all faces coincide!
349 const unsigned int face_no = 0;
350 Assert(fe.reference_cell() == ReferenceCells::get_hypercube<dim>(),
352
353 Assert(couplings.n_rows() == fe.n_components(),
354 ExcDimensionMismatch(couplings.n_rows(), fe.n_components()));
355 Assert(couplings.n_cols() == fe.n_components(),
356 ExcDimensionMismatch(couplings.n_cols(), fe.n_components()));
357 Assert(flux_couplings.n_rows() == fe.n_components(),
358 ExcDimensionMismatch(flux_couplings.n_rows(),
359 fe.n_components()));
360 Assert(flux_couplings.n_cols() == fe.n_components(),
361 ExcDimensionMismatch(flux_couplings.n_cols(),
362 fe.n_components()));
363
364 cell_indices.resize(fe.n_dofs_per_cell());
365 cell->get_mg_dof_indices(cell_indices);
366 unsigned int i = 0;
367 // First, dofs on
368 // vertices. We assume that
369 // each vertex dof couples
370 // with all dofs on
371 // adjacent grid cells.
372
373 // Adding all dofs of the cells
374 // will add dofs of the faces
375 // of the cell adjacent to the
376 // vertex twice. Therefore, we
377 // subtract these here and add
378 // them in a loop over the
379 // faces below.
380
381 // in 1d, faces and vertices
382 // are identical. Nevertheless,
383 // this will only work if
384 // dofs_per_face is zero and
385 // n_dofs_per_vertex() is
386 // arbitrary, not the other way
387 // round.
388 unsigned int increment;
389 while (i < fe.get_first_line_index())
390 {
391 for (unsigned int base = 0; base < fe.n_base_elements(); ++base)
392 for (unsigned int mult = 0; mult < fe.element_multiplicity(base);
393 ++mult)
394 if (couple_cell[fe_index](fe.system_to_block_index(i).first,
395 fe.first_block_of_base(base) +
396 mult) != DoFTools::none)
397 {
398 increment =
399 fe.base_element(base).n_dofs_per_cell() -
400 dim * fe.base_element(base).n_dofs_per_face(face_no);
401 row_lengths[cell_indices[i]] += increment;
402 }
403 ++i;
404 }
405 // From now on, if an object is
406 // a cell, its dofs only couple
407 // inside the cell. Since the
408 // faces are handled below, we
409 // have to subtract ALL faces
410 // in this case.
411
412 // In all other cases we
413 // subtract adjacent faces to be
414 // added in the loop below.
415 while (i < fe.get_first_quad_index(face_no))
416 {
417 for (unsigned int base = 0; base < fe.n_base_elements(); ++base)
418 for (unsigned int mult = 0; mult < fe.element_multiplicity(base);
419 ++mult)
420 if (couple_cell[fe_index](fe.system_to_block_index(i).first,
421 fe.first_block_of_base(base) +
422 mult) != DoFTools::none)
423 {
424 increment =
425 fe.base_element(base).n_dofs_per_cell() -
426 ((dim > 1) ? (dim - 1) :
428 fe.base_element(base).n_dofs_per_face(face_no);
429 row_lengths[cell_indices[i]] += increment;
430 }
431 ++i;
432 }
433
434 // Now quads in 2D and 3D
435 while (i < fe.get_first_hex_index())
436 {
437 for (unsigned int base = 0; base < fe.n_base_elements(); ++base)
438 for (unsigned int mult = 0; mult < fe.element_multiplicity(base);
439 ++mult)
440 if (couple_cell[fe_index](fe.system_to_block_index(i).first,
441 fe.first_block_of_base(base) +
442 mult) != DoFTools::none)
443 {
444 increment =
445 fe.base_element(base).n_dofs_per_cell() -
446 ((dim > 2) ? (dim - 2) :
448 fe.base_element(base).n_dofs_per_face(face_no);
449 row_lengths[cell_indices[i]] += increment;
450 }
451 ++i;
452 }
453
454 // Finally, cells in 3D
455 while (i < fe.n_dofs_per_cell())
456 {
457 for (unsigned int base = 0; base < fe.n_base_elements(); ++base)
458 for (unsigned int mult = 0; mult < fe.element_multiplicity(base);
459 ++mult)
460 if (couple_cell[fe_index](fe.system_to_block_index(i).first,
461 fe.first_block_of_base(base) +
462 mult) != DoFTools::none)
463 {
464 increment =
465 fe.base_element(base).n_dofs_per_cell() -
467 fe.base_element(base).n_dofs_per_face(face_no);
468 row_lengths[cell_indices[i]] += increment;
469 }
470 ++i;
471 }
472
473 // At this point, we have
474 // counted all dofs
475 // contributiong from cells
476 // coupled topologically to the
477 // adjacent cells, but we
478 // subtracted some faces.
479
480 // Now, let's go by the faces
481 // and add the missing
482 // contribution as well as the
483 // flux contributions.
484 for (const unsigned int iface : GeometryInfo<dim>::face_indices())
485 {
486 bool level_boundary = cell->at_boundary(iface);
488 if (!level_boundary)
489 {
490 neighbor = cell->neighbor(iface);
491 if (static_cast<unsigned int>(neighbor->level()) != level)
492 level_boundary = true;
493 }
494
495 if (level_boundary)
496 {
497 for (unsigned int local_dof = 0;
498 local_dof < fe.n_dofs_per_cell();
499 ++local_dof)
500 row_lengths[cell_indices[local_dof]] +=
501 fe.n_dofs_per_face(face_no);
502 continue;
503 }
504
505 const FiniteElement<dim> &nfe = neighbor->get_fe();
507 cell->face(iface);
508
509 // Flux couplings are
510 // computed from both sides
511 // for simplicity.
512
513 // The dofs on the common face
514 // will be handled below,
515 // therefore, we subtract them
516 // here.
517 for (unsigned int base = 0; base < nfe.n_base_elements(); ++base)
518 for (unsigned int mult = 0; mult < nfe.element_multiplicity(base);
519 ++mult)
520 for (unsigned int local_dof = 0;
521 local_dof < fe.n_dofs_per_cell();
522 ++local_dof)
523 if (couple_face[fe_index](
524 fe.system_to_block_index(local_dof).first,
525 nfe.first_block_of_base(base) + mult) != DoFTools::none)
526 {
527 const unsigned int dof_increment =
528 nfe.base_element(base).n_dofs_per_cell() -
529 nfe.base_element(base).n_dofs_per_face(face_no);
530 row_lengths[cell_indices[local_dof]] += dof_increment;
531 }
532
533 // Do this only once per
534 // face and not on the
535 // hanging faces.
536 if (face->user_flag_set())
537 continue;
538 face->set_user_flag();
539 // At this point, we assume
540 // that each cell added its
541 // dofs minus the face to
542 // the couplings of the
543 // face dofs. Since we
544 // subtracted two faces, we
545 // have to re-add one.
546
547 // If one side of the face
548 // is refined, all the fine
549 // face dofs couple with
550 // the coarse one.
551
552 // Wolfgang, do they couple
553 // with each other by
554 // constraints?
555
556 // This will not work with
557 // different couplings on
558 // different cells.
559 neighbor_indices.resize(nfe.n_dofs_per_cell());
560 neighbor->get_mg_dof_indices(neighbor_indices);
561 for (unsigned int base = 0; base < nfe.n_base_elements(); ++base)
562 for (unsigned int mult = 0; mult < nfe.element_multiplicity(base);
563 ++mult)
564 for (unsigned int local_dof = 0;
565 local_dof < fe.n_dofs_per_cell();
566 ++local_dof)
567 if (couple_cell[fe_index](
568 fe.system_to_component_index(local_dof).first,
569 nfe.first_block_of_base(base) + mult) != DoFTools::none)
570 row_lengths[cell_indices[local_dof]] +=
571 nfe.base_element(base).n_dofs_per_face(face_no);
572 for (unsigned int base = 0; base < fe.n_base_elements(); ++base)
573 for (unsigned int mult = 0; mult < fe.element_multiplicity(base);
574 ++mult)
575 for (unsigned int local_dof = 0;
576 local_dof < nfe.n_dofs_per_cell();
577 ++local_dof)
578 if (couple_cell[fe_index](
579 nfe.system_to_component_index(local_dof).first,
580 fe.first_block_of_base(base) + mult) != DoFTools::none)
581 row_lengths[neighbor_indices[local_dof]] +=
582 fe.base_element(base).n_dofs_per_face(face_no);
583 }
584 }
585 user_flags_triangulation.load_user_flags(old_flags);
586 }
587
588
589
590 template <int dim,
591 int spacedim,
592 typename SparsityPatternType,
593 typename number>
594 void
596 SparsityPatternType & sparsity,
597 const unsigned int level,
598 const AffineConstraints<number> &constraints,
599 const bool keep_constrained_dofs)
600 {
601 const types::global_dof_index n_dofs = dof.n_dofs(level);
602 (void)n_dofs;
603
604 Assert(sparsity.n_rows() == n_dofs,
605 ExcDimensionMismatch(sparsity.n_rows(), n_dofs));
606 Assert(sparsity.n_cols() == n_dofs,
607 ExcDimensionMismatch(sparsity.n_cols(), n_dofs));
608
609 const unsigned int dofs_per_cell = dof.get_fe().n_dofs_per_cell();
610 std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
611 for (const auto &cell : dof.cell_iterators_on_level(level))
612 if (cell->is_locally_owned_on_level())
613 {
614 cell->get_mg_dof_indices(dofs_on_this_cell);
615 constraints.add_entries_local_to_global(dofs_on_this_cell,
616 sparsity,
617 keep_constrained_dofs);
618 }
619 }
620
621
622
623 template <int dim,
624 int spacedim,
625 typename SparsityPatternType,
626 typename number>
627 void
629 SparsityPatternType & sparsity,
630 const unsigned int level,
631 const AffineConstraints<number> &constraints,
632 const bool keep_constrained_dofs)
633 {
634 const types::global_dof_index n_dofs = dof.n_dofs(level);
635 (void)n_dofs;
636
637 Assert(sparsity.n_rows() == n_dofs,
638 ExcDimensionMismatch(sparsity.n_rows(), n_dofs));
639 Assert(sparsity.n_cols() == n_dofs,
640 ExcDimensionMismatch(sparsity.n_cols(), n_dofs));
641
642 const unsigned int dofs_per_cell = dof.get_fe().n_dofs_per_cell();
643 std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
644 std::vector<types::global_dof_index> dofs_on_other_cell(dofs_per_cell);
646 endc = dof.end(level);
647 for (; cell != endc; ++cell)
648 {
649 if (!cell->is_locally_owned_on_level())
650 continue;
651
652 cell->get_mg_dof_indices(dofs_on_this_cell);
653 // make sparsity pattern for this cell
654 constraints.add_entries_local_to_global(dofs_on_this_cell,
655 sparsity,
656 keep_constrained_dofs);
657
658 // Loop over all interior neighbors
659 for (const unsigned int face : GeometryInfo<dim>::face_indices())
660 {
661 bool use_face = false;
662 if ((!cell->at_boundary(face)) &&
663 (static_cast<unsigned int>(cell->neighbor_level(face)) ==
664 level))
665 use_face = true;
666 else if (cell->has_periodic_neighbor(face) &&
667 (static_cast<unsigned int>(
668 cell->periodic_neighbor_level(face)) == level))
669 use_face = true;
670
671 if (use_face)
672 {
674 cell->neighbor_or_periodic_neighbor(face);
675 neighbor->get_mg_dof_indices(dofs_on_other_cell);
676 // only add one direction The other is taken care of by
677 // neighbor (except when the neighbor is not owned by the same
678 // processor)
679 constraints.add_entries_local_to_global(dofs_on_this_cell,
680 dofs_on_other_cell,
681 sparsity,
682 keep_constrained_dofs);
683
684 if (neighbor->is_locally_owned_on_level() == false)
685 {
686 constraints.add_entries_local_to_global(
687 dofs_on_other_cell, sparsity, keep_constrained_dofs);
688
689 constraints.add_entries_local_to_global(
690 dofs_on_other_cell,
691 dofs_on_this_cell,
692 sparsity,
693 keep_constrained_dofs);
694 }
695 }
696 }
697 }
698 }
699
700
701
702 template <int dim, typename SparsityPatternType, int spacedim>
703 void
705 SparsityPatternType & sparsity,
706 const unsigned int level)
707 {
708 Assert((level >= 1) && (level < dof.get_triangulation().n_global_levels()),
710
711 const types::global_dof_index fine_dofs = dof.n_dofs(level);
712 const types::global_dof_index coarse_dofs = dof.n_dofs(level - 1);
713 (void)fine_dofs;
714 (void)coarse_dofs;
715
716 // Matrix maps from fine level to coarse level
717
718 Assert(sparsity.n_rows() == coarse_dofs,
719 ExcDimensionMismatch(sparsity.n_rows(), coarse_dofs));
720 Assert(sparsity.n_cols() == fine_dofs,
721 ExcDimensionMismatch(sparsity.n_cols(), fine_dofs));
722
723 const unsigned int dofs_per_cell = dof.get_fe().n_dofs_per_cell();
724 std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
725 std::vector<types::global_dof_index> dofs_on_other_cell(dofs_per_cell);
727 endc = dof.end(level);
728 for (; cell != endc; ++cell)
729 {
730 if (!cell->is_locally_owned_on_level())
731 continue;
732
733 cell->get_mg_dof_indices(dofs_on_this_cell);
734 // Loop over all interior neighbors
735 for (const unsigned int face : GeometryInfo<dim>::face_indices())
736 {
737 // Neighbor is coarser
738 bool use_face = false;
739 if ((!cell->at_boundary(face)) &&
740 (static_cast<unsigned int>(cell->neighbor_level(face)) !=
741 level))
742 use_face = true;
743 else if (cell->has_periodic_neighbor(face) &&
744 (static_cast<unsigned int>(
745 cell->periodic_neighbor_level(face)) != level))
746 use_face = true;
747
748 if (use_face)
749 {
751 cell->neighbor_or_periodic_neighbor(face);
752 neighbor->get_mg_dof_indices(dofs_on_other_cell);
753
754 for (unsigned int i = 0; i < dofs_per_cell; ++i)
755 {
756 for (unsigned int j = 0; j < dofs_per_cell; ++j)
757 {
758 sparsity.add(dofs_on_other_cell[i],
759 dofs_on_this_cell[j]);
760 sparsity.add(dofs_on_other_cell[j],
761 dofs_on_this_cell[i]);
762 }
763 }
764 }
765 }
766 }
767 }
768
769
770
771 template <int dim, typename SparsityPatternType, int spacedim>
772 void
774 SparsityPatternType & sparsity,
775 const unsigned int level,
776 const Table<2, DoFTools::Coupling> &int_mask,
777 const Table<2, DoFTools::Coupling> &flux_mask)
778 {
779 const FiniteElement<dim> & fe = dof.get_fe();
780 const types::global_dof_index n_dofs = dof.n_dofs(level);
781 const unsigned int n_comp = fe.n_components();
782 (void)n_dofs;
783 (void)n_comp;
784
785 Assert(sparsity.n_rows() == n_dofs,
786 ExcDimensionMismatch(sparsity.n_rows(), n_dofs));
787 Assert(sparsity.n_cols() == n_dofs,
788 ExcDimensionMismatch(sparsity.n_cols(), n_dofs));
789 Assert(int_mask.n_rows() == n_comp,
790 ExcDimensionMismatch(int_mask.n_rows(), n_comp));
791 Assert(int_mask.n_cols() == n_comp,
792 ExcDimensionMismatch(int_mask.n_cols(), n_comp));
793 Assert(flux_mask.n_rows() == n_comp,
794 ExcDimensionMismatch(flux_mask.n_rows(), n_comp));
795 Assert(flux_mask.n_cols() == n_comp,
796 ExcDimensionMismatch(flux_mask.n_cols(), n_comp));
797
798 const unsigned int total_dofs = fe.n_dofs_per_cell();
799 std::vector<types::global_dof_index> dofs_on_this_cell(total_dofs);
800 std::vector<types::global_dof_index> dofs_on_other_cell(total_dofs);
801 Table<2, bool> support_on_face(total_dofs,
803
805 endc = dof.end(level);
806
808 int_dof_mask =
810 flux_dof_mask =
812
813 for (unsigned int i = 0; i < total_dofs; ++i)
814 for (auto f : GeometryInfo<dim>::face_indices())
815 support_on_face(i, f) = fe.has_support_on_face(i, f);
816
817 // Clear user flags because we will
818 // need them. But first we save
819 // them and make sure that we
820 // restore them later such that at
821 // the end of this function the
822 // Triangulation will be in the
823 // same state as it was at the
824 // beginning of this function.
825 std::vector<bool> user_flags;
826 dof.get_triangulation().save_user_flags(user_flags);
828 .clear_user_flags();
829
830 for (; cell != endc; ++cell)
831 {
832 if (!cell->is_locally_owned_on_level())
833 continue;
834
835 cell->get_mg_dof_indices(dofs_on_this_cell);
836 // make sparsity pattern for this cell
837 for (unsigned int i = 0; i < total_dofs; ++i)
838 for (unsigned int j = 0; j < total_dofs; ++j)
839 if (int_dof_mask[i][j] != DoFTools::none)
840 sparsity.add(dofs_on_this_cell[i], dofs_on_this_cell[j]);
841
842 // Loop over all interior neighbors
843 for (const unsigned int face : GeometryInfo<dim>::face_indices())
844 {
846 cell->face(face);
847 if (cell_face->user_flag_set())
848 continue;
849
850 if (cell->at_boundary(face) && !cell->has_periodic_neighbor(face))
851 {
852 for (unsigned int i = 0; i < total_dofs; ++i)
853 {
854 const bool i_non_zero_i = support_on_face(i, face);
855 for (unsigned int j = 0; j < total_dofs; ++j)
856 {
857 const bool j_non_zero_i = support_on_face(j, face);
858
859 if (flux_dof_mask(i, j) == DoFTools::always)
860 sparsity.add(dofs_on_this_cell[i],
861 dofs_on_this_cell[j]);
862 if (flux_dof_mask(i, j) == DoFTools::nonzero &&
863 i_non_zero_i && j_non_zero_i)
864 sparsity.add(dofs_on_this_cell[i],
865 dofs_on_this_cell[j]);
866 }
867 }
868 }
869 else
870 {
872 cell->neighbor_or_periodic_neighbor(face);
873
874 if (neighbor->level() < cell->level())
875 continue;
876
877 unsigned int neighbor_face =
878 cell->has_periodic_neighbor(face) ?
879 cell->periodic_neighbor_of_periodic_neighbor(face) :
880 cell->neighbor_of_neighbor(face);
881
882 neighbor->get_mg_dof_indices(dofs_on_other_cell);
883 for (unsigned int i = 0; i < total_dofs; ++i)
884 {
885 const bool i_non_zero_i = support_on_face(i, face);
886 const bool i_non_zero_e = support_on_face(i, neighbor_face);
887 for (unsigned int j = 0; j < total_dofs; ++j)
888 {
889 const bool j_non_zero_i = support_on_face(j, face);
890 const bool j_non_zero_e =
891 support_on_face(j, neighbor_face);
892 if (flux_dof_mask(i, j) == DoFTools::always)
893 {
894 sparsity.add(dofs_on_this_cell[i],
895 dofs_on_other_cell[j]);
896 sparsity.add(dofs_on_other_cell[i],
897 dofs_on_this_cell[j]);
898 sparsity.add(dofs_on_this_cell[i],
899 dofs_on_this_cell[j]);
900 sparsity.add(dofs_on_other_cell[i],
901 dofs_on_other_cell[j]);
902 }
903 if (flux_dof_mask(i, j) == DoFTools::nonzero)
904 {
905 if (i_non_zero_i && j_non_zero_e)
906 sparsity.add(dofs_on_this_cell[i],
907 dofs_on_other_cell[j]);
908 if (i_non_zero_e && j_non_zero_i)
909 sparsity.add(dofs_on_other_cell[i],
910 dofs_on_this_cell[j]);
911 if (i_non_zero_i && j_non_zero_i)
912 sparsity.add(dofs_on_this_cell[i],
913 dofs_on_this_cell[j]);
914 if (i_non_zero_e && j_non_zero_e)
915 sparsity.add(dofs_on_other_cell[i],
916 dofs_on_other_cell[j]);
917 }
918
919 if (flux_dof_mask(j, i) == DoFTools::always)
920 {
921 sparsity.add(dofs_on_this_cell[j],
922 dofs_on_other_cell[i]);
923 sparsity.add(dofs_on_other_cell[j],
924 dofs_on_this_cell[i]);
925 sparsity.add(dofs_on_this_cell[j],
926 dofs_on_this_cell[i]);
927 sparsity.add(dofs_on_other_cell[j],
928 dofs_on_other_cell[i]);
929 }
930 if (flux_dof_mask(j, i) == DoFTools::nonzero)
931 {
932 if (j_non_zero_i && i_non_zero_e)
933 sparsity.add(dofs_on_this_cell[j],
934 dofs_on_other_cell[i]);
935 if (j_non_zero_e && i_non_zero_i)
936 sparsity.add(dofs_on_other_cell[j],
937 dofs_on_this_cell[i]);
938 if (j_non_zero_i && i_non_zero_i)
939 sparsity.add(dofs_on_this_cell[j],
940 dofs_on_this_cell[i]);
941 if (j_non_zero_e && i_non_zero_e)
942 sparsity.add(dofs_on_other_cell[j],
943 dofs_on_other_cell[i]);
944 }
945 }
946 }
947 neighbor->face(neighbor_face)->set_user_flag();
948 }
949 }
950 }
951
952 // finally restore the user flags
954 .load_user_flags(user_flags);
955 }
956
957
958
959 template <int dim, typename SparsityPatternType, int spacedim>
960 void
962 SparsityPatternType & sparsity,
963 const unsigned int level,
964 const Table<2, DoFTools::Coupling> &flux_mask)
965 {
966 const FiniteElement<dim> &fe = dof.get_fe();
967 const unsigned int n_comp = fe.n_components();
968 (void)n_comp;
969
970 Assert((level >= 1) && (level < dof.get_triangulation().n_global_levels()),
972
973 const types::global_dof_index fine_dofs = dof.n_dofs(level);
974 const types::global_dof_index coarse_dofs = dof.n_dofs(level - 1);
975 (void)fine_dofs;
976 (void)coarse_dofs;
977
978 // Matrix maps from fine level to coarse level
979
980 Assert(sparsity.n_rows() == coarse_dofs,
981 ExcDimensionMismatch(sparsity.n_rows(), coarse_dofs));
982 Assert(sparsity.n_cols() == fine_dofs,
983 ExcDimensionMismatch(sparsity.n_cols(), fine_dofs));
984 Assert(flux_mask.n_rows() == n_comp,
985 ExcDimensionMismatch(flux_mask.n_rows(), n_comp));
986 Assert(flux_mask.n_cols() == n_comp,
987 ExcDimensionMismatch(flux_mask.n_cols(), n_comp));
988
989 const unsigned int dofs_per_cell = dof.get_fe().n_dofs_per_cell();
990 std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
991 std::vector<types::global_dof_index> dofs_on_other_cell(dofs_per_cell);
992 Table<2, bool> support_on_face(dofs_per_cell,
994
996 endc = dof.end(level);
997
998 const Table<2, DoFTools::Coupling> flux_dof_mask =
1000
1001 for (unsigned int i = 0; i < dofs_per_cell; ++i)
1002 for (auto f : GeometryInfo<dim>::face_indices())
1003 support_on_face(i, f) = fe.has_support_on_face(i, f);
1004
1005 for (; cell != endc; ++cell)
1006 {
1007 if (!cell->is_locally_owned_on_level())
1008 continue;
1009
1010 cell->get_mg_dof_indices(dofs_on_this_cell);
1011 // Loop over all interior neighbors
1012 for (const unsigned int face : GeometryInfo<dim>::face_indices())
1013 {
1014 // Neighbor is coarser
1015 bool use_face = false;
1016 if ((!cell->at_boundary(face)) &&
1017 (static_cast<unsigned int>(cell->neighbor_level(face)) !=
1018 level))
1019 use_face = true;
1020 else if (cell->has_periodic_neighbor(face) &&
1021 (static_cast<unsigned int>(
1022 cell->periodic_neighbor_level(face)) != level))
1023 use_face = true;
1024
1025 if (use_face)
1026 {
1028 cell->neighbor_or_periodic_neighbor(face);
1029 neighbor->get_mg_dof_indices(dofs_on_other_cell);
1030
1031 for (unsigned int i = 0; i < dofs_per_cell; ++i)
1032 {
1033 for (unsigned int j = 0; j < dofs_per_cell; ++j)
1034 {
1035 if (flux_dof_mask(i, j) != DoFTools::none)
1036 {
1037 sparsity.add(dofs_on_other_cell[i],
1038 dofs_on_this_cell[j]);
1039 sparsity.add(dofs_on_other_cell[j],
1040 dofs_on_this_cell[i]);
1041 }
1042 }
1043 }
1044 }
1045 }
1046 }
1047 }
1048
1049
1050
1051 template <int dim, int spacedim, typename SparsityPatternType>
1052 void
1054 const MGConstrainedDoFs &mg_constrained_dofs,
1055 SparsityPatternType & sparsity,
1056 const unsigned int level)
1057 {
1058 const types::global_dof_index n_dofs = dof.n_dofs(level);
1059 (void)n_dofs;
1060
1061 Assert(sparsity.n_rows() == n_dofs,
1062 ExcDimensionMismatch(sparsity.n_rows(), n_dofs));
1063 Assert(sparsity.n_cols() == n_dofs,
1064 ExcDimensionMismatch(sparsity.n_cols(), n_dofs));
1065
1066 const unsigned int dofs_per_cell = dof.get_fe().n_dofs_per_cell();
1067 std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
1069 endc = dof.end(level);
1070 for (; cell != endc; ++cell)
1071 if (cell->is_locally_owned_on_level())
1072 {
1073 cell->get_mg_dof_indices(dofs_on_this_cell);
1074 for (unsigned int i = 0; i < dofs_per_cell; ++i)
1075 for (unsigned int j = 0; j < dofs_per_cell; ++j)
1076 if (mg_constrained_dofs.is_interface_matrix_entry(
1077 level, dofs_on_this_cell[i], dofs_on_this_cell[j]))
1078 sparsity.add(dofs_on_this_cell[i], dofs_on_this_cell[j]);
1079 }
1080 }
1081
1082
1083
1084 template <int dim, int spacedim>
1085 void
1087 const DoFHandler<dim, spacedim> & dof_handler,
1088 std::vector<std::vector<types::global_dof_index>> &result,
1089 bool only_once,
1090 std::vector<unsigned int> target_component)
1091 {
1092 const FiniteElement<dim> &fe = dof_handler.get_fe();
1093 const unsigned int n_components = fe.n_components();
1094 const unsigned int nlevels =
1095 dof_handler.get_triangulation().n_global_levels();
1096
1097 Assert(result.size() == nlevels,
1098 ExcDimensionMismatch(result.size(), nlevels));
1099
1100 if (target_component.size() == 0)
1101 {
1102 target_component.resize(n_components);
1103 for (unsigned int i = 0; i < n_components; ++i)
1104 target_component[i] = i;
1105 }
1106
1107 Assert(target_component.size() == n_components,
1108 ExcDimensionMismatch(target_component.size(), n_components));
1109
1110 for (unsigned int l = 0; l < nlevels; ++l)
1111 {
1112 result[l].resize(n_components);
1113 std::fill(result[l].begin(), result[l].end(), 0U);
1114
1115 // special case for only one
1116 // component. treat this first
1117 // since it does not require any
1118 // computations
1119 if (n_components == 1)
1120 {
1121 result[l][0] = dof_handler.n_dofs(l);
1122 }
1123 else
1124 {
1125 // otherwise determine the number
1126 // of dofs in each component
1127 // separately. do so in parallel
1128 std::vector<std::vector<bool>> dofs_in_component(
1129 n_components, std::vector<bool>(dof_handler.n_dofs(l), false));
1130 std::vector<ComponentMask> component_select(n_components);
1132 for (unsigned int i = 0; i < n_components; ++i)
1133 {
1134 void (*fun_ptr)(const unsigned int level,
1136 const ComponentMask &,
1137 std::vector<bool> &) =
1138 &DoFTools::extract_level_dofs<dim, spacedim>;
1139
1140 std::vector<bool> tmp(n_components, false);
1141 tmp[i] = true;
1142 component_select[i] = ComponentMask(tmp);
1143
1144 tasks += Threads::new_task(fun_ptr,
1145 l,
1146 dof_handler,
1147 component_select[i],
1148 dofs_in_component[i]);
1149 }
1150 tasks.join_all();
1151
1152 // next count what we got
1153 unsigned int component = 0;
1154 for (unsigned int b = 0; b < fe.n_base_elements(); ++b)
1155 {
1156 const FiniteElement<dim> &base = fe.base_element(b);
1157 // Dimension of base element
1158 unsigned int d = base.n_components();
1159
1160 for (unsigned int m = 0; m < fe.element_multiplicity(b); ++m)
1161 {
1162 for (unsigned int dd = 0; dd < d; ++dd)
1163 {
1164 if (base.is_primitive() || (!only_once || dd == 0))
1165 result[l][target_component[component]] +=
1166 std::count(dofs_in_component[component].begin(),
1167 dofs_in_component[component].end(),
1168 true);
1169 ++component;
1170 }
1171 }
1172 }
1173 // finally sanity check
1174 Assert(!dof_handler.get_fe().is_primitive() ||
1175 std::accumulate(result[l].begin(),
1176 result[l].end(),
1178 dof_handler.n_dofs(l),
1180 }
1181 }
1182 }
1183
1184
1185
1186 template <int dim, int spacedim>
1187 void
1189 const DoFHandler<dim, spacedim> & dof_handler,
1190 std::vector<std::vector<types::global_dof_index>> &dofs_per_block,
1191 std::vector<unsigned int> target_block)
1192 {
1193 const FiniteElement<dim, spacedim> &fe = dof_handler.get_fe();
1194 const unsigned int n_blocks = fe.n_blocks();
1195 const unsigned int n_levels =
1196 dof_handler.get_triangulation().n_global_levels();
1197
1198 AssertDimension(dofs_per_block.size(), n_levels);
1199
1200 for (unsigned int l = 0; l < n_levels; ++l)
1201 std::fill(dofs_per_block[l].begin(), dofs_per_block[l].end(), 0U);
1202 // If the empty vector was given as
1203 // default argument, set up this
1204 // vector as identity.
1205 if (target_block.size() == 0)
1206 {
1207 target_block.resize(n_blocks);
1208 for (unsigned int i = 0; i < n_blocks; ++i)
1209 target_block[i] = i;
1210 }
1211 Assert(target_block.size() == n_blocks,
1212 ExcDimensionMismatch(target_block.size(), n_blocks));
1213
1214 const unsigned int max_block =
1215 *std::max_element(target_block.begin(), target_block.end());
1216 const unsigned int n_target_blocks = max_block + 1;
1217 (void)n_target_blocks;
1218
1219 for (unsigned int l = 0; l < n_levels; ++l)
1220 AssertDimension(dofs_per_block[l].size(), n_target_blocks);
1221
1222 // special case for only one
1223 // block. treat this first
1224 // since it does not require any
1225 // computations
1226 if (n_blocks == 1)
1227 {
1228 for (unsigned int l = 0; l < n_levels; ++l)
1229 dofs_per_block[l][0] = dof_handler.n_dofs(l);
1230 return;
1231 }
1232 // otherwise determine the number
1233 // of dofs in each block
1234 // separately. do so in parallel
1235 for (unsigned int l = 0; l < n_levels; ++l)
1236 {
1237 std::vector<std::vector<bool>> dofs_in_block(
1238 n_blocks, std::vector<bool>(dof_handler.n_dofs(l), false));
1239 std::vector<BlockMask> block_select(n_blocks);
1241 for (unsigned int i = 0; i < n_blocks; ++i)
1242 {
1243 void (*fun_ptr)(const unsigned int level,
1245 const BlockMask &,
1246 std::vector<bool> &) =
1247 &DoFTools::extract_level_dofs<dim, spacedim>;
1248
1249 std::vector<bool> tmp(n_blocks, false);
1250 tmp[i] = true;
1251 block_select[i] = tmp;
1252
1253 tasks += Threads::new_task(
1254 fun_ptr, l, dof_handler, block_select[i], dofs_in_block[i]);
1255 }
1256 tasks.join_all();
1257
1258 // next count what we got
1259 for (unsigned int block = 0; block < fe.n_blocks(); ++block)
1260 dofs_per_block[l][target_block[block]] +=
1261 std::count(dofs_in_block[block].begin(),
1262 dofs_in_block[block].end(),
1263 true);
1264 }
1265 }
1266
1267
1268
1269 template <int dim, int spacedim>
1270 void
1272 const DoFHandler<dim, spacedim> &dof,
1273 const std::map<types::boundary_id, const Function<spacedim> *>
1274 & function_map,
1275 std::vector<std::set<types::global_dof_index>> &boundary_indices,
1276 const ComponentMask & component_mask)
1277 {
1278 Assert(boundary_indices.size() == dof.get_triangulation().n_global_levels(),
1279 ExcDimensionMismatch(boundary_indices.size(),
1281
1282 std::set<types::boundary_id> boundary_ids;
1283 for (const auto &boundary_function : function_map)
1284 boundary_ids.insert(boundary_function.first);
1285
1286 std::vector<IndexSet> boundary_indexset;
1287 make_boundary_list(dof, boundary_ids, boundary_indexset, component_mask);
1288 for (unsigned int i = 0; i < dof.get_triangulation().n_global_levels(); ++i)
1289 boundary_indices[i].insert(boundary_indexset[i].begin(),
1290 boundary_indexset[i].end());
1291 }
1292
1293
1294 template <int dim, int spacedim>
1295 void
1297 const std::map<types::boundary_id,
1298 const Function<spacedim> *> &function_map,
1299 std::vector<IndexSet> &boundary_indices,
1300 const ComponentMask & component_mask)
1301 {
1302 Assert(boundary_indices.size() == dof.get_triangulation().n_global_levels(),
1303 ExcDimensionMismatch(boundary_indices.size(),
1305
1306 std::set<types::boundary_id> boundary_ids;
1307 for (const auto &boundary_function : function_map)
1308 boundary_ids.insert(boundary_function.first);
1309
1310 make_boundary_list(dof, boundary_ids, boundary_indices, component_mask);
1311 }
1312
1313
1314
1315 template <int dim, int spacedim>
1316 void
1318 const std::set<types::boundary_id> &boundary_ids,
1319 std::vector<IndexSet> & boundary_indices,
1320 const ComponentMask & component_mask)
1321 {
1322 boundary_indices.resize(dof.get_triangulation().n_global_levels());
1323
1324 // if for whatever reason we were passed an empty set, return immediately
1325 if (boundary_ids.size() == 0)
1326 return;
1327
1328 for (unsigned int i = 0; i < dof.get_triangulation().n_global_levels(); ++i)
1329 if (boundary_indices[i].size() == 0)
1330 boundary_indices[i] = IndexSet(dof.n_dofs(i));
1331
1332 const unsigned int n_components = dof.get_fe_collection().n_components();
1333 const bool fe_is_system = (n_components != 1);
1334
1335 std::vector<types::global_dof_index> local_dofs;
1336 local_dofs.reserve(dof.get_fe_collection().max_dofs_per_face());
1337 std::fill(local_dofs.begin(), local_dofs.end(), numbers::invalid_dof_index);
1338
1339 std::vector<std::vector<types::global_dof_index>> dofs_by_level(
1340 dof.get_triangulation().n_levels());
1341
1342 // First, deal with the simpler case when we have to identify all boundary
1343 // dofs
1344 if (component_mask.n_selected_components(n_components) == n_components)
1345 {
1346 for (const auto &cell : dof.cell_iterators())
1347 {
1348 if (cell->is_artificial_on_level())
1349 continue;
1350 const FiniteElement<dim> &fe = cell->get_fe();
1351 const unsigned int level = cell->level();
1352
1353 for (const unsigned int face_no : GeometryInfo<dim>::face_indices())
1354 if (cell->at_boundary(face_no) == true)
1355 {
1356 const typename DoFHandler<dim, spacedim>::face_iterator face =
1357 cell->face(face_no);
1358 const types::boundary_id bi = face->boundary_id();
1359 // Face is listed in boundary map
1360 if (boundary_ids.find(bi) != boundary_ids.end())
1361 {
1362 local_dofs.resize(fe.n_dofs_per_face(face_no));
1363 face->get_mg_dof_indices(level, local_dofs);
1364 dofs_by_level[level].insert(dofs_by_level[level].end(),
1365 local_dofs.begin(),
1366 local_dofs.end());
1367 }
1368 }
1369 }
1370 }
1371 else
1372 {
1373 Assert(component_mask.n_selected_components(n_components) > 0,
1374 ExcMessage(
1375 "It's probably worthwhile to select at least one component."));
1376
1377 for (const auto &cell : dof.cell_iterators())
1378 if (!cell->is_artificial_on_level())
1379 for (const unsigned int face_no : GeometryInfo<dim>::face_indices())
1380 {
1381 if (cell->at_boundary(face_no) == false)
1382 continue;
1383
1384 const FiniteElement<dim> &fe = cell->get_fe();
1385 const unsigned int level = cell->level();
1386
1388 cell->face(face_no);
1389 const types::boundary_id boundary_component =
1390 face->boundary_id();
1391 if (boundary_ids.find(boundary_component) != boundary_ids.end())
1392 // we want to constrain this boundary
1393 {
1394 for (unsigned int i = 0;
1395 i < cell->get_fe().n_dofs_per_cell();
1396 ++i)
1397 {
1398 const ComponentMask &nonzero_component_array =
1399 cell->get_fe().get_nonzero_components(i);
1400 // if we want to constrain one of the nonzero
1401 // components, we have to constrain all of them
1402
1403 bool selected = false;
1404 for (unsigned int c = 0; c < n_components; ++c)
1405 if (nonzero_component_array[c] == true &&
1406 component_mask[c] == true)
1407 {
1408 selected = true;
1409 break;
1410 }
1411 if (selected)
1412 for (unsigned int c = 0; c < n_components; ++c)
1413 Assert(
1414 nonzero_component_array[c] == false ||
1415 component_mask[c] == true,
1416 ExcMessage(
1417 "You are using a non-primitive FiniteElement "
1418 "and try to constrain just some of its components!"));
1419 }
1420
1421 // get indices, physical location and boundary values of
1422 // dofs on this face
1423 local_dofs.resize(fe.n_dofs_per_face(face_no));
1424 face->get_mg_dof_indices(level, local_dofs);
1425 if (fe_is_system)
1426 {
1427 for (unsigned int i = 0; i < local_dofs.size(); ++i)
1428 {
1429 unsigned int component =
1431 if (fe.is_primitive())
1432 component =
1433 fe.face_system_to_component_index(i, face_no)
1434 .first;
1435 else
1436 {
1437 // Just pick the first of the components
1438 // We already know that either all or none
1439 // of the components are selected
1440 const ComponentMask &nonzero_component_array =
1441 cell->get_fe().get_nonzero_components(i);
1442 for (unsigned int c = 0; c < n_components; ++c)
1443 if (nonzero_component_array[c] == true)
1444 {
1445 component = c;
1446 break;
1447 }
1448 }
1451 if (component_mask[component] == true)
1452 dofs_by_level[level].push_back(local_dofs[i]);
1453 }
1454 }
1455 else
1456 dofs_by_level[level].insert(dofs_by_level[level].end(),
1457 local_dofs.begin(),
1458 local_dofs.end());
1459 }
1460 }
1461 }
1462 for (unsigned int level = 0; level < dof.get_triangulation().n_levels();
1463 ++level)
1464 {
1465 std::sort(dofs_by_level[level].begin(), dofs_by_level[level].end());
1466 boundary_indices[level].add_indices(
1467 dofs_by_level[level].begin(),
1468 std::unique(dofs_by_level[level].begin(),
1469 dofs_by_level[level].end()));
1470 }
1471 }
1472
1473
1474
1475 template <int dim, int spacedim>
1476 void
1478 std::vector<IndexSet> & interface_dofs)
1479 {
1480 Assert(interface_dofs.size() ==
1481 mg_dof_handler.get_triangulation().n_global_levels(),
1483 interface_dofs.size(),
1484 mg_dof_handler.get_triangulation().n_global_levels()));
1485
1486 std::vector<std::vector<types::global_dof_index>> tmp_interface_dofs(
1487 interface_dofs.size());
1488
1489 const FiniteElement<dim, spacedim> &fe = mg_dof_handler.get_fe();
1490
1491 const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
1492
1493 std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1494
1495 std::vector<bool> cell_dofs(dofs_per_cell, false);
1496
1497 typename DoFHandler<dim>::cell_iterator cell = mg_dof_handler.begin(),
1498 endc = mg_dof_handler.end();
1499
1500 for (; cell != endc; ++cell)
1501 {
1502 // Do not look at artificial level cells (in a serial computation we
1503 // need to ignore the level_subdomain_id() because it is never set).
1504 if (cell->is_artificial_on_level())
1505 continue;
1506
1507 bool has_coarser_neighbor = false;
1508
1509 std::fill(cell_dofs.begin(), cell_dofs.end(), false);
1510
1511 for (const unsigned int face_nr : GeometryInfo<dim>::face_indices())
1512 {
1513 const typename DoFHandler<dim, spacedim>::face_iterator face =
1514 cell->face(face_nr);
1515 if (!face->at_boundary() || cell->has_periodic_neighbor(face_nr))
1516 {
1517 // interior face
1518 const typename DoFHandler<dim>::cell_iterator neighbor =
1519 cell->neighbor_or_periodic_neighbor(face_nr);
1520
1521 // only process cell pairs if one or both of them are owned by
1522 // me (ignore if running in serial)
1523 if (neighbor->is_artificial_on_level())
1524 continue;
1525
1526 // Do refinement face from the coarse side
1527 if (neighbor->level() < cell->level())
1528 {
1529 for (unsigned int j = 0; j < fe.n_dofs_per_face(face_nr);
1530 ++j)
1531 cell_dofs[fe.face_to_cell_index(j, face_nr)] = true;
1532
1533 has_coarser_neighbor = true;
1534 }
1535 }
1536 }
1537
1538 if (has_coarser_neighbor == false)
1539 continue;
1540
1541 const unsigned int level = cell->level();
1542 cell->get_mg_dof_indices(local_dof_indices);
1543
1544 for (unsigned int i = 0; i < dofs_per_cell; ++i)
1545 {
1546 if (cell_dofs[i])
1547 tmp_interface_dofs[level].push_back(local_dof_indices[i]);
1548 }
1549 }
1550
1551 for (unsigned int l = 0;
1552 l < mg_dof_handler.get_triangulation().n_global_levels();
1553 ++l)
1554 {
1555 interface_dofs[l].clear();
1556 std::sort(tmp_interface_dofs[l].begin(), tmp_interface_dofs[l].end());
1557 interface_dofs[l].add_indices(tmp_interface_dofs[l].begin(),
1558 std::unique(tmp_interface_dofs[l].begin(),
1559 tmp_interface_dofs[l].end()));
1560 interface_dofs[l].compress();
1561 }
1562 }
1563
1564
1565
1566 template <int dim, int spacedim>
1567 unsigned int
1569 {
1570 // Find minimum level for an active cell in
1571 // this locally owned subdomain
1572 // Note: with the way active cells are traversed,
1573 // the first locally owned cell we find will have
1574 // the lowest level in the particular subdomain.
1575 unsigned int min_level = tria.n_global_levels();
1577 cell = tria.begin_active(),
1578 endc = tria.end();
1579 for (; cell != endc; ++cell)
1580 if (cell->is_locally_owned())
1581 {
1582 min_level = cell->level();
1583 break;
1584 }
1585
1586 unsigned int global_min = min_level;
1587 // If necessary, communicate to find minimum
1588 // level for an active cell over all subdomains
1590 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1591 &tria))
1592 global_min = Utilities::MPI::min(min_level, tr->get_communicator());
1593
1594 AssertIndexRange(global_min, tria.n_global_levels());
1595
1596 return global_min;
1597 }
1598
1599
1600
1601 namespace internal
1602 {
1603 double
1605 const std::vector<types::global_dof_index> &n_cells_on_levels,
1606 const MPI_Comm comm)
1607 {
1608 std::vector<types::global_dof_index> n_cells_on_leves_max(
1609 n_cells_on_levels.size());
1610 std::vector<types::global_dof_index> n_cells_on_leves_sum(
1611 n_cells_on_levels.size());
1612
1613 Utilities::MPI::max(n_cells_on_levels, comm, n_cells_on_leves_max);
1614 Utilities::MPI::sum(n_cells_on_levels, comm, n_cells_on_leves_sum);
1615
1616 const unsigned int n_proc = Utilities::MPI::n_mpi_processes(comm);
1617
1618 const double ideal_work = std::accumulate(n_cells_on_leves_sum.begin(),
1619 n_cells_on_leves_sum.end(),
1620 0) /
1621 static_cast<double>(n_proc);
1622 const double workload_imbalance =
1623 std::accumulate(n_cells_on_leves_max.begin(),
1624 n_cells_on_leves_max.end(),
1625 0) /
1626 ideal_work;
1627
1628 return workload_imbalance;
1629 }
1630
1631
1632
1633 double
1635 const std::vector<
1636 std::pair<types::global_dof_index, types::global_dof_index>> &cells,
1637 const MPI_Comm comm)
1638 {
1639 std::vector<types::global_dof_index> cells_local(cells.size());
1640 std::vector<types::global_dof_index> cells_remote(cells.size());
1641
1642 for (unsigned int i = 0; i < cells.size(); ++i)
1643 {
1644 cells_local[i] = cells[i].first;
1645 cells_remote[i] = cells[i].second;
1646 }
1647
1648 std::vector<types::global_dof_index> cells_local_sum(cells_local.size());
1649 Utilities::MPI::sum(cells_local, comm, cells_local_sum);
1650
1651 std::vector<types::global_dof_index> cells_remote_sum(
1652 cells_remote.size());
1653 Utilities::MPI::sum(cells_remote, comm, cells_remote_sum);
1654
1655 const auto n_cells_local =
1656 std::accumulate(cells_local_sum.begin(), cells_local_sum.end(), 0);
1657 const auto n_cells_remote =
1658 std::accumulate(cells_remote_sum.begin(), cells_remote_sum.end(), 0);
1659
1660 return static_cast<double>(n_cells_local) /
1661 (n_cells_local + n_cells_remote);
1662 }
1663 } // namespace internal
1664
1665
1666
1667 template <int dim, int spacedim>
1668 std::vector<types::global_dof_index>
1670 {
1672 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1673 &tria))
1674 Assert(
1675 tr->is_multilevel_hierarchy_constructed(),
1676 ExcMessage(
1677 "We can only compute the workload imbalance if the multilevel hierarchy has been constructed!"));
1678
1679 const unsigned int n_global_levels = tria.n_global_levels();
1680
1681 std::vector<types::global_dof_index> n_cells_on_levels(n_global_levels);
1682
1683 for (unsigned int lvl = 0; lvl < n_global_levels; ++lvl)
1684 for (const auto &cell : tria.cell_iterators_on_level(lvl))
1685 if (cell->is_locally_owned_on_level())
1686 ++n_cells_on_levels[lvl];
1687
1688 return n_cells_on_levels;
1689 }
1690
1691
1692
1693 template <int dim, int spacedim>
1694 std::vector<types::global_dof_index>
1696 const std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
1697 &trias)
1698 {
1699 const unsigned int n_global_levels = trias.size();
1700
1701 std::vector<types::global_dof_index> n_cells_on_levels(n_global_levels);
1702
1703 for (unsigned int lvl = 0; lvl < n_global_levels; ++lvl)
1704 for (const auto &cell : trias[lvl]->active_cell_iterators())
1705 if (cell->is_locally_owned())
1706 ++n_cells_on_levels[lvl];
1707
1708 return n_cells_on_levels;
1709 }
1710
1711
1712
1713 template <int dim, int spacedim>
1714 double
1716 {
1719 }
1720
1721
1722
1723 template <int dim, int spacedim>
1724 double
1726 const std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
1727 &trias)
1728 {
1730 trias.back()->get_communicator());
1731 }
1732
1733
1734
1735 template <int dim, int spacedim>
1736 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>
1738 {
1739 const unsigned int n_global_levels = tria.n_global_levels();
1740
1741 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>
1742 cells(n_global_levels);
1743
1744 const MPI_Comm communicator = tria.get_communicator();
1745
1746 const unsigned int my_rank = Utilities::MPI::this_mpi_process(communicator);
1747
1748 for (unsigned int lvl = 0; lvl < n_global_levels - 1; ++lvl)
1749 for (const auto &cell : tria.cell_iterators_on_level(lvl))
1750 if (cell->is_locally_owned_on_level() && cell->has_children())
1751 for (unsigned int i = 0; i < GeometryInfo<dim>::max_children_per_cell;
1752 ++i)
1753 {
1754 const auto level_subdomain_id =
1755 cell->child(i)->level_subdomain_id();
1756 if (level_subdomain_id == my_rank)
1757 ++cells[lvl + 1].first;
1758 else if (level_subdomain_id != numbers::invalid_unsigned_int)
1759 ++cells[lvl + 1].second;
1760 else
1762 }
1763
1764 return cells;
1765 }
1766
1767
1768
1769 template <int dim, int spacedim>
1770 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>
1772 const std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
1773 &trias)
1774 {
1775 const unsigned int n_global_levels = trias.size();
1776
1777 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>
1778 cells(n_global_levels);
1779
1780 const MPI_Comm communicator = trias.back()->get_communicator();
1781
1782 const unsigned int my_rank = Utilities::MPI::this_mpi_process(communicator);
1783
1784 for (unsigned int lvl = 0; lvl < n_global_levels - 1; ++lvl)
1785 {
1786 const auto &tria_coarse = *trias[lvl];
1787 const auto &tria_fine = *trias[lvl + 1];
1788
1789 const unsigned int n_coarse_cells = tria_fine.n_global_coarse_cells();
1790 const unsigned int n_global_levels = tria_fine.n_global_levels();
1791
1792 const ::internal::CellIDTranslator<dim> cell_id_translator(
1793 n_coarse_cells, n_global_levels);
1794
1795 IndexSet is_fine_owned(cell_id_translator.size());
1796 IndexSet is_fine_required(cell_id_translator.size());
1797
1798 for (const auto &cell : tria_fine.active_cell_iterators())
1799 if (!cell->is_artificial() && cell->is_locally_owned())
1800 is_fine_owned.add_index(cell_id_translator.translate(cell));
1801
1802 for (const auto &cell : tria_coarse.active_cell_iterators())
1803 if (!cell->is_artificial() && cell->is_locally_owned())
1804 {
1805 if (cell->level() + 1u == tria_fine.n_global_levels())
1806 continue;
1807
1808 for (unsigned int i = 0;
1809 i < GeometryInfo<dim>::max_children_per_cell;
1810 ++i)
1811 is_fine_required.add_index(
1812 cell_id_translator.translate(cell, i));
1813 }
1814
1815 std::vector<unsigned int> is_fine_required_ranks(
1816 is_fine_required.n_elements());
1817
1819 process(is_fine_owned,
1820 is_fine_required,
1821 communicator,
1822 is_fine_required_ranks,
1823 false);
1824
1826 std::vector<
1827 std::pair<types::global_cell_index, types::global_cell_index>>,
1828 std::vector<unsigned int>>
1829 consensus_algorithm;
1830 consensus_algorithm.run(process, communicator);
1831
1832 for (unsigned i = 0; i < is_fine_required.n_elements(); ++i)
1833 if (is_fine_required_ranks[i] == my_rank)
1834 ++cells[lvl + 1].first;
1835 else if (is_fine_required_ranks[i] != numbers::invalid_unsigned_int)
1836 ++cells[lvl + 1].second;
1837 }
1838
1839 return cells;
1840 }
1841
1842
1843
1844 template <int dim, int spacedim>
1845 double
1847 {
1850 }
1851
1852
1853
1854 template <int dim, int spacedim>
1855 double
1857 const std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
1858 &trias)
1859 {
1862 trias.back()->get_communicator());
1863 }
1864
1865} // namespace MGTools
1866
1867
1868// explicit instantiations
1869#include "mg_tools.inst"
1870
void add_entries_local_to_global(const std::vector< size_type > &local_dof_indices, SparsityPatternType &sparsity_pattern, const bool keep_constrained_entries=true, const Table< 2, bool > &dof_mask=Table< 2, bool >()) const
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
cell_iterator end() const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
const Triangulation< dim, spacedim > & get_triangulation() const
types::global_dof_index n_dofs() const
cell_iterator begin(const unsigned int level=0) const
void clear()
unsigned int get_first_line_index() const
unsigned int n_dofs_per_cell() const
unsigned int get_first_quad_index(const unsigned int quad_no=0) const
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
unsigned int n_blocks() const
unsigned int n_components() const
ReferenceCell reference_cell() const
unsigned int get_first_hex_index() const
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
bool is_primitive() const
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
unsigned int element_multiplicity(const unsigned int index) const
unsigned int n_base_elements() const
std::pair< unsigned int, unsigned int > face_system_to_component_index(const unsigned int index, const unsigned int face_no=0) const
types::global_dof_index first_block_of_base(const unsigned int b) const
size_type n_elements() const
Definition: index_set.h:1834
void add_index(const size_type index)
Definition: index_set.h:1655
bool is_interface_matrix_entry(const unsigned int level, const types::global_dof_index i, const types::global_dof_index j) const
virtual MPI_Comm get_communicator() const
void save_user_flags(std::ostream &out) const
unsigned int n_levels() const
cell_iterator end() const
void clear_user_flags()
virtual unsigned int n_global_levels() const
void load_user_flags(std::istream &in)
active_cell_iterator begin_active(const unsigned int level=0) const
virtual std::vector< unsigned int > run(const std::vector< unsigned int > &targets, const std::function< RequestType(const unsigned int)> &create_request, const std::function< AnswerType(const unsigned int, const RequestType &)> &answer_request, const std::function< void(const unsigned int, const AnswerType &)> &process_answer, const MPI_Comm &comm) override
unsigned int max_dofs_per_face() const
unsigned int n_components() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
unsigned int level
Definition: grid_out.cc:4606
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
IteratorRange< cell_iterator > cell_iterators() const
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Definition: exceptions.h:1473
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1667
#define AssertIndexRange(index, range)
Definition: exceptions.h:1732
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1583
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:466
typename ActiveSelector::face_iterator face_iterator
Definition: dof_handler.h:484
Task< RT > new_task(const std::function< RT()> &function)
void convert_couplings_to_blocks(const DoFHandler< dim, spacedim > &dof_handler, const Table< 2, Coupling > &table_by_component, std::vector< Table< 2, Coupling > > &tables_by_block)
Definition: dof_tools.cc:2380
Table< 2, Coupling > dof_couplings_from_component_couplings(const FiniteElement< dim, spacedim > &fe, const Table< 2, Coupling > &component_couplings)
double workload_imbalance(const std::vector< types::global_dof_index > &n_cells_on_levels, const MPI_Comm comm)
Definition: mg_tools.cc:1604
double vertical_communication_efficiency(const std::vector< std::pair< types::global_dof_index, types::global_dof_index > > &cells, const MPI_Comm comm)
Definition: mg_tools.cc:1634
void count_dofs_per_component(const DoFHandler< dim, spacedim > &mg_dof, std::vector< std::vector< types::global_dof_index > > &result, const bool only_once=false, std::vector< unsigned int > target_component={})
Definition: mg_tools.cc:1086
std::vector< std::pair< types::global_dof_index, types::global_dof_index > > local_vertical_communication_cost(const Triangulation< dim, spacedim > &tria)
Definition: mg_tools.cc:1737
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity, const unsigned int level, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const bool keep_constrained_dofs=true)
Definition: mg_tools.cc:595
std::vector< types::global_dof_index > local_workload(const Triangulation< dim, spacedim > &tria)
Definition: mg_tools.cc:1669
void make_boundary_list(const DoFHandler< dim, spacedim > &mg_dof, const std::map< types::boundary_id, const Function< spacedim > * > &function_map, std::vector< std::set< types::global_dof_index > > &boundary_indices, const ComponentMask &component_mask=ComponentMask())
Definition: mg_tools.cc:1271
void compute_row_length_vector(const DoFHandler< dim, spacedim > &dofs, const unsigned int level, std::vector< unsigned int > &row_lengths, const DoFTools::Coupling flux_couplings=DoFTools::none)
Definition: mg_tools.cc:107
void extract_inner_interface_dofs(const DoFHandler< dim, spacedim > &mg_dof_handler, std::vector< IndexSet > &interface_dofs)
Definition: mg_tools.cc:1477
double workload_imbalance(const Triangulation< dim, spacedim > &tria)
Definition: mg_tools.cc:1715
void count_dofs_per_block(const DoFHandler< dim, spacedim > &dof_handler, std::vector< std::vector< types::global_dof_index > > &dofs_per_block, std::vector< unsigned int > target_block={})
Definition: mg_tools.cc:1188
void make_interface_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, const MGConstrainedDoFs &mg_constrained_dofs, SparsityPatternType &sparsity, const unsigned int level)
Definition: mg_tools.cc:1053
void make_flux_sparsity_pattern_edge(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity, const unsigned int level)
Definition: mg_tools.cc:704
double vertical_communication_efficiency(const Triangulation< dim, spacedim > &tria)
Definition: mg_tools.cc:1846
void make_flux_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity, const unsigned int level, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const bool keep_constrained_dofs=true)
Definition: mg_tools.cc:628
unsigned int max_level_for_coarse_mesh(const Triangulation< dim, spacedim > &tria)
Definition: mg_tools.cc:1568
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:151
T min(const T &t, const MPI_Comm &mpi_communicator)
T sum(const T &t, const MPI_Comm &mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:140
T max(const T &t, const MPI_Comm &mpi_communicator)
static const unsigned int invalid_unsigned_int
Definition: types.h:201
const types::global_dof_index invalid_dof_index
Definition: types.h:216
const MPI_Comm & comm
const ::Triangulation< dim, spacedim > & tria