deal.II version GIT relicensing-2206-gaa53ff9447 2024-12-02 09:10:00+00:00
\(\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
tria.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2010 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
18#include <deal.II/base/point.h>
20
23
25#include <deal.II/grid/tria.h>
28
29#include <algorithm>
30#include <fstream>
31#include <iostream>
32#include <limits>
33#include <numeric>
34
35
37
38
39namespace internal
40{
41 namespace parallel
42 {
43 namespace distributed
44 {
45 namespace TriangulationImplementation
46 {
54 template <int dim, int spacedim>
55 void
58 {
59 auto pack =
61 &cell) -> std::uint8_t {
62 if (cell->refine_flag_set())
63 return 1;
64 if (cell->coarsen_flag_set())
65 return 2;
66 return 0;
67 };
68
69 auto unpack =
71 &cell,
72 const std::uint8_t &flag) -> void {
73 cell->clear_coarsen_flag();
74 cell->clear_refine_flag();
75 if (flag == 1)
76 cell->set_refine_flag();
77 else if (flag == 2)
78 cell->set_coarsen_flag();
79 };
80
81 GridTools::exchange_cell_data_to_ghosts<std::uint8_t>(tria,
82 pack,
83 unpack);
84 }
85 } // namespace TriangulationImplementation
86 } // namespace distributed
87 } // namespace parallel
88} // namespace internal
89
90
91
92#ifdef DEAL_II_WITH_P4EST
93
94namespace
95{
96 template <int dim, int spacedim>
97 void
98 get_vertex_to_cell_mappings(
100 std::vector<unsigned int> &vertex_touch_count,
101 std::vector<std::list<
103 unsigned int>>> &vertex_to_cell)
104 {
105 vertex_touch_count.resize(triangulation.n_vertices());
106 vertex_to_cell.resize(triangulation.n_vertices());
107
108 for (const auto &cell : triangulation.active_cell_iterators())
109 for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
110 {
111 ++vertex_touch_count[cell->vertex_index(v)];
112 vertex_to_cell[cell->vertex_index(v)].emplace_back(cell, v);
113 }
114 }
115
116
117
118 template <int dim, int spacedim>
119 void
120 get_edge_to_cell_mappings(
122 std::vector<unsigned int> &edge_touch_count,
123 std::vector<std::list<
125 unsigned int>>> &edge_to_cell)
126 {
127 Assert(triangulation.n_levels() == 1, ExcInternalError());
128
129 edge_touch_count.resize(triangulation.n_active_lines());
130 edge_to_cell.resize(triangulation.n_active_lines());
131
132 for (const auto &cell : triangulation.active_cell_iterators())
133 for (unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++l)
134 {
135 ++edge_touch_count[cell->line(l)->index()];
136 edge_to_cell[cell->line(l)->index()].emplace_back(cell, l);
137 }
138 }
139
140
141
146 template <int dim, int spacedim>
147 void
148 set_vertex_and_cell_info(
150 const std::vector<unsigned int> &vertex_touch_count,
151 const std::vector<std::list<
153 unsigned int>>> &vertex_to_cell,
154 const std::vector<types::global_dof_index>
155 &coarse_cell_to_p4est_tree_permutation,
156 const bool set_vertex_info,
157 typename internal::p4est::types<dim>::connectivity *connectivity)
158 {
159 // copy the vertices into the connectivity structure. the triangulation
160 // exports the array of vertices, but some of the entries are sometimes
161 // unused; this shouldn't be the case for a newly created triangulation,
162 // but make sure
163 //
164 // note that p4est stores coordinates as a triplet of values even in 2d
165 Assert(triangulation.get_used_vertices().size() ==
166 triangulation.get_vertices().size(),
168 Assert(std::find(triangulation.get_used_vertices().begin(),
169 triangulation.get_used_vertices().end(),
170 false) == triangulation.get_used_vertices().end(),
172 if (set_vertex_info == true)
173 for (unsigned int v = 0; v < triangulation.n_vertices(); ++v)
174 {
175 connectivity->vertices[3 * v] = triangulation.get_vertices()[v][0];
176 connectivity->vertices[3 * v + 1] =
177 triangulation.get_vertices()[v][1];
178 connectivity->vertices[3 * v + 2] =
179 (spacedim == 2 ? 0 : triangulation.get_vertices()[v][2]);
180 }
181
182 // next store the tree_to_vertex indices (each tree is here only a single
183 // cell in the coarse mesh). p4est requires vertex numbering in clockwise
184 // orientation
185 //
186 // while we're at it, also copy the neighborship information between cells
188 cell = triangulation.begin_active(),
189 endc = triangulation.end();
190 for (; cell != endc; ++cell)
191 {
192 const unsigned int index =
193 coarse_cell_to_p4est_tree_permutation[cell->index()];
194
195 for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
196 {
197 if (set_vertex_info == true)
198 connectivity
200 v] = cell->vertex_index(v);
201 connectivity
203 v] = cell->vertex_index(v);
204 }
205
206 // neighborship information. if a cell is at a boundary, then enter
207 // the index of the cell itself here
208 for (auto f : GeometryInfo<dim>::face_indices())
209 if (cell->face(f)->at_boundary() == false)
210 connectivity
211 ->tree_to_tree[index * GeometryInfo<dim>::faces_per_cell + f] =
212 coarse_cell_to_p4est_tree_permutation[cell->neighbor(f)->index()];
213 else
214 connectivity
215 ->tree_to_tree[index * GeometryInfo<dim>::faces_per_cell + f] =
216 coarse_cell_to_p4est_tree_permutation[cell->index()];
217
218 // fill tree_to_face, which is essentially neighbor_to_neighbor;
219 // however, we have to remap the resulting face number as well
220 for (auto f : GeometryInfo<dim>::face_indices())
221 if (cell->face(f)->at_boundary() == false)
222 {
223 switch (dim)
224 {
225 case 2:
226 {
227 connectivity->tree_to_face
229 cell->neighbor_of_neighbor(f);
230 break;
231 }
232
233 case 3:
234 {
235 /*
236 * The values for tree_to_face are in 0..23 where ttf % 6
237 * gives the face number and ttf / 4 the face orientation
238 * code. The orientation is determined as follows. Let
239 * my_face and other_face be the two face numbers of the
240 * connecting trees in 0..5. Then the first face vertex
241 * of the lower of my_face and other_face connects to a
242 * face vertex numbered 0..3 in the higher of my_face and
243 * other_face. The face orientation is defined as this
244 * number. If my_face == other_face, treating either of
245 * both faces as the lower one leads to the same result.
246 */
247
248 connectivity->tree_to_face[index * 6 + f] =
249 cell->neighbor_of_neighbor(f);
250
251 unsigned int face_idx_list[2] = {
252 f, cell->neighbor_of_neighbor(f)};
254 cell_list[2] = {cell, cell->neighbor(f)};
255 unsigned int smaller_idx = 0;
256
257 if (f > cell->neighbor_of_neighbor(f))
258 smaller_idx = 1;
259
260 unsigned int larger_idx = (smaller_idx + 1) % 2;
261 // smaller = *_list[smaller_idx]
262 // larger = *_list[larger_idx]
263
264 unsigned int v = 0;
265
266 // global vertex index of vertex 0 on face of cell with
267 // smaller local face index
268 unsigned int g_idx = cell_list[smaller_idx]->vertex_index(
270 face_idx_list[smaller_idx],
271 0,
272 cell_list[smaller_idx]->face_orientation(
273 face_idx_list[smaller_idx]),
274 cell_list[smaller_idx]->face_flip(
275 face_idx_list[smaller_idx]),
276 cell_list[smaller_idx]->face_rotation(
277 face_idx_list[smaller_idx])));
278
279 // loop over vertices on face from other cell and compare
280 // global vertex numbers
281 for (unsigned int i = 0;
282 i < GeometryInfo<dim>::vertices_per_face;
283 ++i)
284 {
285 unsigned int idx =
286 cell_list[larger_idx]->vertex_index(
288 face_idx_list[larger_idx], i));
289
290 if (idx == g_idx)
291 {
292 v = i;
293 break;
294 }
295 }
296
297 connectivity->tree_to_face[index * 6 + f] += 6 * v;
298 break;
299 }
300
301 default:
303 }
304 }
305 else
306 connectivity
307 ->tree_to_face[index * GeometryInfo<dim>::faces_per_cell + f] = f;
308 }
309
310 // now fill the vertex information
311 connectivity->ctt_offset[0] = 0;
312 std::partial_sum(vertex_touch_count.begin(),
313 vertex_touch_count.end(),
314 &connectivity->ctt_offset[1]);
315
316 [[maybe_unused]] const typename internal::p4est::types<dim>::locidx
317 num_vtt = std::accumulate(vertex_touch_count.begin(),
318 vertex_touch_count.end(),
319 0u);
320 Assert(connectivity->ctt_offset[triangulation.n_vertices()] == num_vtt,
322
323 for (unsigned int v = 0; v < triangulation.n_vertices(); ++v)
324 {
325 Assert(vertex_to_cell[v].size() == vertex_touch_count[v],
327
328 typename std::list<
329 std::pair<typename Triangulation<dim, spacedim>::active_cell_iterator,
330 unsigned int>>::const_iterator p =
331 vertex_to_cell[v].begin();
332 for (unsigned int c = 0; c < vertex_touch_count[v]; ++c, ++p)
333 {
334 connectivity->corner_to_tree[connectivity->ctt_offset[v] + c] =
335 coarse_cell_to_p4est_tree_permutation[p->first->index()];
336 connectivity->corner_to_corner[connectivity->ctt_offset[v] + c] =
337 p->second;
338 }
339 }
340 }
341
342
343
344 template <int dim, int spacedim>
345 bool
347 const typename internal::p4est::types<dim>::forest *parallel_forest,
348 const typename internal::p4est::types<dim>::topidx coarse_grid_cell)
349 {
350 Assert(coarse_grid_cell < parallel_forest->connectivity->num_trees,
352 return ((coarse_grid_cell >= parallel_forest->first_local_tree) &&
353 (coarse_grid_cell <= parallel_forest->last_local_tree));
354 }
355
356
357 template <int dim, int spacedim>
358 void
359 delete_all_children_and_self(
361 {
362 if (cell->has_children())
363 for (unsigned int c = 0; c < cell->n_children(); ++c)
364 delete_all_children_and_self<dim, spacedim>(cell->child(c));
365 else
366 cell->set_coarsen_flag();
367 }
368
369
370
371 template <int dim, int spacedim>
372 void
373 delete_all_children(
375 {
376 if (cell->has_children())
377 for (unsigned int c = 0; c < cell->n_children(); ++c)
378 delete_all_children_and_self<dim, spacedim>(cell->child(c));
379 }
380
381
382 template <int dim, int spacedim>
383 void
384 determine_level_subdomain_id_recursively(
385 const typename internal::p4est::types<dim>::tree &tree,
386 const typename internal::p4est::types<dim>::locidx &tree_index,
387 const typename Triangulation<dim, spacedim>::cell_iterator &dealii_cell,
388 const typename internal::p4est::types<dim>::quadrant &p4est_cell,
390 const types::subdomain_id my_subdomain,
391 const std::vector<std::vector<bool>> &marked_vertices)
392 {
393 if (dealii_cell->level_subdomain_id() == numbers::artificial_subdomain_id)
394 {
395 // important: only assign the level_subdomain_id if it is a ghost cell
396 // even though we could fill in all.
397 bool used = false;
398 for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
399 {
400 if (marked_vertices[dealii_cell->level()]
401 [dealii_cell->vertex_index(v)])
402 {
403 used = true;
404 break;
405 }
406 }
407
408 // Special case: if this cell is active we might be a ghost neighbor
409 // to a locally owned cell across a vertex that is finer.
410 // Example (M= my, O=dealii_cell, owned by somebody else):
411 // *------*
412 // | |
413 // | O |
414 // | |
415 // *---*---*------*
416 // | M | M |
417 // *---*---*
418 // | | M |
419 // *---*---*
420 if (!used && dealii_cell->is_active() &&
421 dealii_cell->is_artificial() == false &&
422 dealii_cell->level() + 1 < static_cast<int>(marked_vertices.size()))
423 {
424 for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
425 {
426 if (marked_vertices[dealii_cell->level() + 1]
427 [dealii_cell->vertex_index(v)])
428 {
429 used = true;
430 break;
431 }
432 }
433 }
434
435 // Like above, but now the other way around
436 if (!used && dealii_cell->is_active() &&
437 dealii_cell->is_artificial() == false && dealii_cell->level() > 0)
438 {
439 for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
440 {
441 if (marked_vertices[dealii_cell->level() - 1]
442 [dealii_cell->vertex_index(v)])
443 {
444 used = true;
445 break;
446 }
447 }
448 }
449
450 if (used)
451 {
453 &forest, tree_index, &p4est_cell, my_subdomain);
454 Assert((owner != -2) && (owner != -1),
455 ExcMessage("p4est should know the owner."));
456 dealii_cell->set_level_subdomain_id(owner);
457 }
458 }
459
460 if (dealii_cell->has_children())
461 {
464 for (unsigned int c = 0; c < GeometryInfo<dim>::max_children_per_cell;
465 ++c)
466 switch (dim)
467 {
468 case 2:
469 P4EST_QUADRANT_INIT(&p4est_child[c]);
470 break;
471 case 3:
472 P8EST_QUADRANT_INIT(&p4est_child[c]);
473 break;
474 default:
476 }
477
478
480 p4est_child);
481
482 for (unsigned int c = 0; c < GeometryInfo<dim>::max_children_per_cell;
483 ++c)
484 {
485 determine_level_subdomain_id_recursively<dim, spacedim>(
486 tree,
487 tree_index,
488 dealii_cell->child(c),
489 p4est_child[c],
490 forest,
491 my_subdomain,
492 marked_vertices);
493 }
494 }
495 }
496
497
498 template <int dim, int spacedim>
499 void
500 match_tree_recursively(
501 const typename internal::p4est::types<dim>::tree &tree,
502 const typename Triangulation<dim, spacedim>::cell_iterator &dealii_cell,
503 const typename internal::p4est::types<dim>::quadrant &p4est_cell,
504 const typename internal::p4est::types<dim>::forest &forest,
505 const types::subdomain_id my_subdomain)
506 {
507 // check if this cell exists in the local p4est cell
508 if (sc_array_bsearch(const_cast<sc_array_t *>(&tree.quadrants),
509 &p4est_cell,
511 -1)
512 {
513 // yes, cell found in local part of p4est
514 delete_all_children<dim, spacedim>(dealii_cell);
515 if (dealii_cell->is_active())
516 dealii_cell->set_subdomain_id(my_subdomain);
517 }
518 else
519 {
520 // no, cell not found in local part of p4est. this means that the
521 // local part is more refined than the current cell. if this cell has
522 // no children of its own, we need to refine it, and if it does
523 // already have children then loop over all children and see if they
524 // are locally available as well
525 if (dealii_cell->is_active())
526 dealii_cell->set_refine_flag();
527 else
528 {
531 for (unsigned int c = 0;
532 c < GeometryInfo<dim>::max_children_per_cell;
533 ++c)
534 switch (dim)
535 {
536 case 2:
537 P4EST_QUADRANT_INIT(&p4est_child[c]);
538 break;
539 case 3:
540 P8EST_QUADRANT_INIT(&p4est_child[c]);
541 break;
542 default:
544 }
545
546
548 p4est_child);
549
550 for (unsigned int c = 0;
551 c < GeometryInfo<dim>::max_children_per_cell;
552 ++c)
554 const_cast<typename internal::p4est::types<dim>::tree *>(
555 &tree),
556 &p4est_child[c]) == false)
557 {
558 // no, this child is locally not available in the p4est.
559 // delete all its children but, because this may not be
560 // successful, make sure to mark all children recursively
561 // as not local.
562 delete_all_children<dim, spacedim>(dealii_cell->child(c));
563 dealii_cell->child(c)->recursively_set_subdomain_id(
565 }
566 else
567 {
568 // at least some part of the tree rooted in this child is
569 // locally available
570 match_tree_recursively<dim, spacedim>(tree,
571 dealii_cell->child(c),
572 p4est_child[c],
573 forest,
574 my_subdomain);
575 }
576 }
577 }
578 }
579
580
581 template <int dim, int spacedim>
582 void
583 match_quadrant(
584 const ::Triangulation<dim, spacedim> *tria,
585 unsigned int dealii_index,
586 const typename internal::p4est::types<dim>::quadrant &ghost_quadrant,
587 types::subdomain_id ghost_owner)
588 {
589 const int l = ghost_quadrant.level;
590
591 for (int i = 0; i < l; ++i)
592 {
594 i,
595 dealii_index);
596 if (cell->is_active())
597 {
598 cell->clear_coarsen_flag();
599 cell->set_refine_flag();
600 return;
601 }
602
603 const int child_id =
605 i + 1);
606 dealii_index = cell->child_index(child_id);
607 }
608
610 l,
611 dealii_index);
612 if (cell->has_children())
613 delete_all_children<dim, spacedim>(cell);
614 else
615 {
616 cell->clear_coarsen_flag();
617 cell->set_subdomain_id(ghost_owner);
618 }
619 }
620
621 template <int dim>
622 class PartitionSearch
623 {
624 public:
625 PartitionSearch()
626 {
627 Assert(dim > 1, ExcNotImplemented());
628 }
629
630 PartitionSearch(const PartitionSearch<dim> &other) = delete;
631
632 PartitionSearch<dim> &
633 operator=(const PartitionSearch<dim> &other) = delete;
634
635 public:
645 static int
646 local_quadrant_fn(typename internal::p4est::types<dim>::forest *forest,
647 typename internal::p4est::types<dim>::topidx which_tree,
648 typename internal::p4est::types<dim>::quadrant *quadrant,
649 int rank_begin,
650 int rank_end,
651 void *point);
652
665 static int
666 local_point_fn(typename internal::p4est::types<dim>::forest *forest,
667 typename internal::p4est::types<dim>::topidx which_tree,
668 typename internal::p4est::types<dim>::quadrant *quadrant,
669 int rank_begin,
670 int rank_end,
671 void *point);
672
673 private:
678 class QuadrantData
679 {
680 public:
681 QuadrantData();
682
683 void
684 set_cell_vertices(
686 typename internal::p4est::types<dim>::topidx which_tree,
687 typename internal::p4est::types<dim>::quadrant *quadrant,
689 quad_length_on_level);
690
691 void
692 initialize_mapping();
693
695 map_real_to_unit_cell(const Point<dim> &p) const;
696
697 bool
698 is_in_this_quadrant(const Point<dim> &p) const;
699
700 private:
701 std::vector<Point<dim>> cell_vertices;
702
707 FullMatrix<double> quadrant_mapping_matrix;
708
709 bool are_vertices_initialized;
710
711 bool is_reference_mapping_initialized;
712 };
713
717 QuadrantData quadrant_data;
718 }; // class PartitionSearch
719
720
721
722 template <int dim>
723 int
724 PartitionSearch<dim>::local_quadrant_fn(
726 typename internal::p4est::types<dim>::topidx which_tree,
727 typename internal::p4est::types<dim>::quadrant *quadrant,
728 int /* rank_begin */,
729 int /* rank_end */,
730 void * /* this is always nullptr */ point)
731 {
732 // point must be nullptr here
733 Assert(point == nullptr, ::ExcInternalError());
734
735 // we need the user pointer
736 // note that this is not available since function is static
737 PartitionSearch<dim> *this_object =
738 reinterpret_cast<PartitionSearch<dim> *>(forest->user_pointer);
739
740 // Avoid p4est macros, instead do bitshifts manually with fixed size types
742 quad_length_on_level =
743 1 << (static_cast<typename internal::p4est::types<dim>::quadrant_coord>(
744 (dim == 2 ? P4EST_MAXLEVEL : P8EST_MAXLEVEL)) -
746 quadrant->level));
747
748 this_object->quadrant_data.set_cell_vertices(forest,
749 which_tree,
750 quadrant,
751 quad_length_on_level);
752
753 // from cell vertices we can initialize the mapping
754 this_object->quadrant_data.initialize_mapping();
755
756 // always return true since we must decide by point
757 return /* true */ 1;
758 }
759
760
761
762 template <int dim>
763 int
764 PartitionSearch<dim>::local_point_fn(
766 typename internal::p4est::types<dim>::topidx /* which_tree */,
767 typename internal::p4est::types<dim>::quadrant * /* quadrant */,
768 int rank_begin,
769 int rank_end,
770 void *point)
771 {
772 // point must NOT be be nullptr here
773 Assert(point != nullptr, ::ExcInternalError());
774
775 // we need the user pointer
776 // note that this is not available since function is static
777 PartitionSearch<dim> *this_object =
778 reinterpret_cast<PartitionSearch<dim> *>(forest->user_pointer);
779
780 // point with rank as double pointer
781 double *this_point_dptr = static_cast<double *>(point);
782
783 Point<dim> this_point =
784 (dim == 2 ? Point<dim>(this_point_dptr[0], this_point_dptr[1]) :
785 Point<dim>(this_point_dptr[0],
786 this_point_dptr[1],
787 this_point_dptr[2]));
788
789 // use reference mapping to decide whether this point is in this quadrant
790 const bool is_in_this_quadrant =
791 this_object->quadrant_data.is_in_this_quadrant(this_point);
792
793
794
795 if (!is_in_this_quadrant)
796 {
797 // no need to search further, stop recursion
798 return /* false */ 0;
799 }
800
801
802
803 // From here we have a candidate
804 if (rank_begin < rank_end)
805 {
806 // continue recursion
807 return /* true */ 1;
808 }
809
810 // Now, we know that the point is found (rank_begin==rank_end) and we have
811 // the MPI rank, so no need to search further.
812 this_point_dptr[dim] = static_cast<double>(rank_begin);
813
814 // stop recursion.
815 return /* false */ 0;
816 }
817
818
819
820 template <int dim>
821 bool
822 PartitionSearch<dim>::QuadrantData::is_in_this_quadrant(
823 const Point<dim> &p) const
824 {
825 const Point<dim> p_ref = map_real_to_unit_cell(p);
826
828 }
829
830
831
832 template <int dim>
834 PartitionSearch<dim>::QuadrantData::map_real_to_unit_cell(
835 const Point<dim> &p) const
836 {
837 Assert(is_reference_mapping_initialized,
839 "Cell vertices and mapping coefficients must be fully "
840 "initialized before transforming a point to the unit cell."));
841
842 Point<dim> p_out;
843
844 if (dim == 2)
845 {
846 for (unsigned int alpha = 0;
847 alpha < GeometryInfo<dim>::vertices_per_cell;
848 ++alpha)
849 {
850 const Point<dim> &p_ref =
852
853 p_out += (quadrant_mapping_matrix(alpha, 0) +
854 quadrant_mapping_matrix(alpha, 1) * p(0) +
855 quadrant_mapping_matrix(alpha, 2) * p(1) +
856 quadrant_mapping_matrix(alpha, 3) * p(0) * p(1)) *
857 p_ref;
858 }
859 }
860 else
861 {
862 for (unsigned int alpha = 0;
863 alpha < GeometryInfo<dim>::vertices_per_cell;
864 ++alpha)
865 {
866 const Point<dim> &p_ref =
868
869 p_out += (quadrant_mapping_matrix(alpha, 0) +
870 quadrant_mapping_matrix(alpha, 1) * p(0) +
871 quadrant_mapping_matrix(alpha, 2) * p(1) +
872 quadrant_mapping_matrix(alpha, 3) * p(2) +
873 quadrant_mapping_matrix(alpha, 4) * p(0) * p(1) +
874 quadrant_mapping_matrix(alpha, 5) * p(1) * p(2) +
875 quadrant_mapping_matrix(alpha, 6) * p(0) * p(2) +
876 quadrant_mapping_matrix(alpha, 7) * p(0) * p(1) * p(2)) *
877 p_ref;
878 }
879 }
880
881 return p_out;
882 }
883
884
885 template <int dim>
886 PartitionSearch<dim>::QuadrantData::QuadrantData()
887 : cell_vertices(GeometryInfo<dim>::vertices_per_cell)
888 , quadrant_mapping_matrix(GeometryInfo<dim>::vertices_per_cell,
889 GeometryInfo<dim>::vertices_per_cell)
890 , are_vertices_initialized(false)
891 , is_reference_mapping_initialized(false)
892 {}
893
894
895
896 template <int dim>
897 void
898 PartitionSearch<dim>::QuadrantData::initialize_mapping()
899 {
900 Assert(
901 are_vertices_initialized,
903 "Cell vertices must be initialized before the cell mapping can be filled."));
904
907
908 if (dim == 2)
909 {
910 for (unsigned int alpha = 0;
911 alpha < GeometryInfo<dim>::vertices_per_cell;
912 ++alpha)
913 {
914 // point matrix to be inverted
915 point_matrix(0, alpha) = 1;
916 point_matrix(1, alpha) = cell_vertices[alpha](0);
917 point_matrix(2, alpha) = cell_vertices[alpha](1);
918 point_matrix(3, alpha) =
919 cell_vertices[alpha](0) * cell_vertices[alpha](1);
920 }
921
922 /*
923 * Rows of quadrant_mapping_matrix are the coefficients of the basis
924 * on the physical cell
925 */
926 quadrant_mapping_matrix.invert(point_matrix);
927 }
928 else
929 {
930 for (unsigned int alpha = 0;
931 alpha < GeometryInfo<dim>::vertices_per_cell;
932 ++alpha)
933 {
934 // point matrix to be inverted
935 point_matrix(0, alpha) = 1;
936 point_matrix(1, alpha) = cell_vertices[alpha](0);
937 point_matrix(2, alpha) = cell_vertices[alpha](1);
938 point_matrix(3, alpha) = cell_vertices[alpha](2);
939 point_matrix(4, alpha) =
940 cell_vertices[alpha](0) * cell_vertices[alpha](1);
941 point_matrix(5, alpha) =
942 cell_vertices[alpha](1) * cell_vertices[alpha](2);
943 point_matrix(6, alpha) =
944 cell_vertices[alpha](0) * cell_vertices[alpha](2);
945 point_matrix(7, alpha) = cell_vertices[alpha](0) *
946 cell_vertices[alpha](1) *
947 cell_vertices[alpha](2);
948 }
949
950 /*
951 * Rows of quadrant_mapping_matrix are the coefficients of the basis
952 * on the physical cell
953 */
954 quadrant_mapping_matrix.invert(point_matrix);
955 }
956
957 is_reference_mapping_initialized = true;
958 }
959
960
961
962 template <>
963 void
964 PartitionSearch<2>::QuadrantData::set_cell_vertices(
965 typename internal::p4est::types<2>::forest *forest,
966 typename internal::p4est::types<2>::topidx which_tree,
967 typename internal::p4est::types<2>::quadrant *quadrant,
969 quad_length_on_level)
970 {
971 constexpr unsigned int dim = 2;
972
973 // p4est for some reason always needs double vxyz[3] as last argument to
974 // quadrant_coord_to_vertex
975 double corner_point[dim + 1] = {0};
976
977 // A lambda to avoid code duplication.
978 const auto copy_vertex = [&](unsigned int vertex_index) -> void {
979 // copy into local struct
980 for (unsigned int d = 0; d < dim; ++d)
981 {
982 cell_vertices[vertex_index](d) = corner_point[d];
983 // reset
984 corner_point[d] = 0;
985 }
986 };
987
988 // Fill points of QuadrantData in lexicographic order
989 /*
990 * Corner #0
991 */
992 unsigned int vertex_index = 0;
994 forest->connectivity, which_tree, quadrant->x, quadrant->y, corner_point);
995
996 // copy into local struct
997 copy_vertex(vertex_index);
998
999 /*
1000 * Corner #1
1001 */
1002 vertex_index = 1;
1004 forest->connectivity,
1005 which_tree,
1006 quadrant->x + quad_length_on_level,
1007 quadrant->y,
1008 corner_point);
1009
1010 // copy into local struct
1011 copy_vertex(vertex_index);
1012
1013 /*
1014 * Corner #2
1015 */
1016 vertex_index = 2;
1018 forest->connectivity,
1019 which_tree,
1020 quadrant->x,
1021 quadrant->y + quad_length_on_level,
1022 corner_point);
1023
1024 // copy into local struct
1025 copy_vertex(vertex_index);
1026
1027 /*
1028 * Corner #3
1029 */
1030 vertex_index = 3;
1032 forest->connectivity,
1033 which_tree,
1034 quadrant->x + quad_length_on_level,
1035 quadrant->y + quad_length_on_level,
1036 corner_point);
1037
1038 // copy into local struct
1039 copy_vertex(vertex_index);
1040
1041 are_vertices_initialized = true;
1042 }
1043
1044
1045
1046 template <>
1047 void
1048 PartitionSearch<3>::QuadrantData::set_cell_vertices(
1049 typename internal::p4est::types<3>::forest *forest,
1050 typename internal::p4est::types<3>::topidx which_tree,
1051 typename internal::p4est::types<3>::quadrant *quadrant,
1053 quad_length_on_level)
1054 {
1055 constexpr unsigned int dim = 3;
1056
1057 double corner_point[dim] = {0};
1058
1059 // A lambda to avoid code duplication.
1060 auto copy_vertex = [&](unsigned int vertex_index) -> void {
1061 // copy into local struct
1062 for (unsigned int d = 0; d < dim; ++d)
1063 {
1064 cell_vertices[vertex_index](d) = corner_point[d];
1065 // reset
1066 corner_point[d] = 0;
1067 }
1068 };
1069
1070 // Fill points of QuadrantData in lexicographic order
1071 /*
1072 * Corner #0
1073 */
1074 unsigned int vertex_index = 0;
1076 forest->connectivity,
1077 which_tree,
1078 quadrant->x,
1079 quadrant->y,
1080 quadrant->z,
1081 corner_point);
1082
1083 // copy into local struct
1084 copy_vertex(vertex_index);
1085
1086
1087 /*
1088 * Corner #1
1089 */
1090 vertex_index = 1;
1092 forest->connectivity,
1093 which_tree,
1094 quadrant->x + quad_length_on_level,
1095 quadrant->y,
1096 quadrant->z,
1097 corner_point);
1098
1099 // copy into local struct
1100 copy_vertex(vertex_index);
1101
1102 /*
1103 * Corner #2
1104 */
1105 vertex_index = 2;
1107 forest->connectivity,
1108 which_tree,
1109 quadrant->x,
1110 quadrant->y + quad_length_on_level,
1111 quadrant->z,
1112 corner_point);
1113
1114 // copy into local struct
1115 copy_vertex(vertex_index);
1116
1117 /*
1118 * Corner #3
1119 */
1120 vertex_index = 3;
1122 forest->connectivity,
1123 which_tree,
1124 quadrant->x + quad_length_on_level,
1125 quadrant->y + quad_length_on_level,
1126 quadrant->z,
1127 corner_point);
1128
1129 // copy into local struct
1130 copy_vertex(vertex_index);
1131
1132 /*
1133 * Corner #4
1134 */
1135 vertex_index = 4;
1137 forest->connectivity,
1138 which_tree,
1139 quadrant->x,
1140 quadrant->y,
1141 quadrant->z + quad_length_on_level,
1142 corner_point);
1143
1144 // copy into local struct
1145 copy_vertex(vertex_index);
1146
1147 /*
1148 * Corner #5
1149 */
1150 vertex_index = 5;
1152 forest->connectivity,
1153 which_tree,
1154 quadrant->x + quad_length_on_level,
1155 quadrant->y,
1156 quadrant->z + quad_length_on_level,
1157 corner_point);
1158
1159 // copy into local struct
1160 copy_vertex(vertex_index);
1161
1162 /*
1163 * Corner #6
1164 */
1165 vertex_index = 6;
1167 forest->connectivity,
1168 which_tree,
1169 quadrant->x,
1170 quadrant->y + quad_length_on_level,
1171 quadrant->z + quad_length_on_level,
1172 corner_point);
1173
1174 // copy into local struct
1175 copy_vertex(vertex_index);
1176
1177 /*
1178 * Corner #7
1179 */
1180 vertex_index = 7;
1182 forest->connectivity,
1183 which_tree,
1184 quadrant->x + quad_length_on_level,
1185 quadrant->y + quad_length_on_level,
1186 quadrant->z + quad_length_on_level,
1187 corner_point);
1188
1189 // copy into local struct
1190 copy_vertex(vertex_index);
1191
1192
1193 are_vertices_initialized = true;
1194 }
1195
1196
1197
1203 template <int dim, int spacedim>
1204 class RefineAndCoarsenList
1205 {
1206 public:
1207 RefineAndCoarsenList(const Triangulation<dim, spacedim> &triangulation,
1208 const std::vector<types::global_dof_index>
1209 &p4est_tree_to_coarse_cell_permutation,
1210 const types::subdomain_id my_subdomain);
1211
1220 static int
1221 refine_callback(
1222 typename internal::p4est::types<dim>::forest *forest,
1223 typename internal::p4est::types<dim>::topidx coarse_cell_index,
1224 typename internal::p4est::types<dim>::quadrant *quadrant);
1225
1230 static int
1231 coarsen_callback(
1232 typename internal::p4est::types<dim>::forest *forest,
1233 typename internal::p4est::types<dim>::topidx coarse_cell_index,
1234 typename internal::p4est::types<dim>::quadrant *children[]);
1235
1236 bool
1237 pointers_are_at_end() const;
1238
1239 private:
1240 std::vector<typename internal::p4est::types<dim>::quadrant> refine_list;
1241 typename std::vector<typename internal::p4est::types<dim>::quadrant>::
1242 const_iterator current_refine_pointer;
1243
1244 std::vector<typename internal::p4est::types<dim>::quadrant> coarsen_list;
1245 typename std::vector<typename internal::p4est::types<dim>::quadrant>::
1246 const_iterator current_coarsen_pointer;
1247
1248 void
1249 build_lists(
1251 const typename internal::p4est::types<dim>::quadrant &p4est_cell,
1252 const types::subdomain_id myid);
1253 };
1254
1255
1256
1257 template <int dim, int spacedim>
1258 bool
1259 RefineAndCoarsenList<dim, spacedim>::pointers_are_at_end() const
1260 {
1261 return ((current_refine_pointer == refine_list.end()) &&
1262 (current_coarsen_pointer == coarsen_list.end()));
1263 }
1264
1265
1266
1267 template <int dim, int spacedim>
1268 RefineAndCoarsenList<dim, spacedim>::RefineAndCoarsenList(
1270 const std::vector<types::global_dof_index>
1271 &p4est_tree_to_coarse_cell_permutation,
1272 const types::subdomain_id my_subdomain)
1273 {
1274 // count how many flags are set and allocate that much memory
1275 unsigned int n_refine_flags = 0, n_coarsen_flags = 0;
1276 for (const auto &cell : triangulation.active_cell_iterators())
1277 {
1278 // skip cells that are not local
1279 if (cell->subdomain_id() != my_subdomain)
1280 continue;
1281
1282 if (cell->refine_flag_set())
1283 ++n_refine_flags;
1284 else if (cell->coarsen_flag_set())
1285 ++n_coarsen_flags;
1286 }
1287
1288 refine_list.reserve(n_refine_flags);
1289 coarsen_list.reserve(n_coarsen_flags);
1290
1291
1292 // now build the lists of cells that are flagged. note that p4est will
1293 // traverse its cells in the order in which trees appear in the
1294 // forest. this order is not the same as the order of coarse cells in the
1295 // deal.II Triangulation because we have translated everything by the
1296 // coarse_cell_to_p4est_tree_permutation permutation. in order to make
1297 // sure that the output array is already in the correct order, traverse
1298 // our coarse cells in the same order in which p4est will:
1299 for (unsigned int c = 0; c < triangulation.n_cells(0); ++c)
1300 {
1301 unsigned int coarse_cell_index =
1302 p4est_tree_to_coarse_cell_permutation[c];
1303
1305 &triangulation, 0, coarse_cell_index);
1306
1307 typename internal::p4est::types<dim>::quadrant p4est_cell;
1309 /*level=*/0,
1310 /*index=*/0);
1311 p4est_cell.p.which_tree = c;
1312 build_lists(cell, p4est_cell, my_subdomain);
1313 }
1314
1315
1316 Assert(refine_list.size() == n_refine_flags, ExcInternalError());
1317 Assert(coarsen_list.size() == n_coarsen_flags, ExcInternalError());
1318
1319 // make sure that our ordering in fact worked
1320 for (unsigned int i = 1; i < refine_list.size(); ++i)
1321 Assert(refine_list[i].p.which_tree >= refine_list[i - 1].p.which_tree,
1323 for (unsigned int i = 1; i < coarsen_list.size(); ++i)
1324 Assert(coarsen_list[i].p.which_tree >= coarsen_list[i - 1].p.which_tree,
1326
1327 current_refine_pointer = refine_list.begin();
1328 current_coarsen_pointer = coarsen_list.begin();
1329 }
1330
1331
1332
1333 template <int dim, int spacedim>
1334 void
1335 RefineAndCoarsenList<dim, spacedim>::build_lists(
1337 const typename internal::p4est::types<dim>::quadrant &p4est_cell,
1338 const types::subdomain_id my_subdomain)
1339 {
1340 if (cell->is_active())
1341 {
1342 if (cell->subdomain_id() == my_subdomain)
1343 {
1344 if (cell->refine_flag_set())
1345 refine_list.push_back(p4est_cell);
1346 else if (cell->coarsen_flag_set())
1347 coarsen_list.push_back(p4est_cell);
1348 }
1349 }
1350 else
1351 {
1354 for (unsigned int c = 0; c < GeometryInfo<dim>::max_children_per_cell;
1355 ++c)
1356 switch (dim)
1357 {
1358 case 2:
1359 P4EST_QUADRANT_INIT(&p4est_child[c]);
1360 break;
1361 case 3:
1362 P8EST_QUADRANT_INIT(&p4est_child[c]);
1363 break;
1364 default:
1366 }
1368 p4est_child);
1369 for (unsigned int c = 0; c < GeometryInfo<dim>::max_children_per_cell;
1370 ++c)
1371 {
1372 p4est_child[c].p.which_tree = p4est_cell.p.which_tree;
1373 build_lists(cell->child(c), p4est_child[c], my_subdomain);
1374 }
1375 }
1376 }
1377
1378
1379 template <int dim, int spacedim>
1380 int
1381 RefineAndCoarsenList<dim, spacedim>::refine_callback(
1382 typename internal::p4est::types<dim>::forest *forest,
1383 typename internal::p4est::types<dim>::topidx coarse_cell_index,
1384 typename internal::p4est::types<dim>::quadrant *quadrant)
1385 {
1386 RefineAndCoarsenList<dim, spacedim> *this_object =
1387 reinterpret_cast<RefineAndCoarsenList<dim, spacedim> *>(
1388 forest->user_pointer);
1389
1390 // if there are no more cells in our list the current cell can't be
1391 // flagged for refinement
1392 if (this_object->current_refine_pointer == this_object->refine_list.end())
1393 return 0;
1394
1395 Assert(coarse_cell_index <=
1396 this_object->current_refine_pointer->p.which_tree,
1398
1399 // if p4est hasn't yet reached the tree of the next flagged cell the
1400 // current cell can't be flagged for refinement
1401 if (coarse_cell_index < this_object->current_refine_pointer->p.which_tree)
1402 return 0;
1403
1404 // now we're in the right tree in the forest
1405 Assert(coarse_cell_index <=
1406 this_object->current_refine_pointer->p.which_tree,
1408
1409 // make sure that the p4est loop over cells hasn't gotten ahead of our own
1410 // pointer
1412 quadrant, &*this_object->current_refine_pointer) <= 0,
1414
1415 // now, if the p4est cell is one in the list, it is supposed to be refined
1417 quadrant, &*this_object->current_refine_pointer))
1418 {
1419 ++this_object->current_refine_pointer;
1420 return 1;
1421 }
1422
1423 // p4est cell is not in list
1424 return 0;
1425 }
1426
1427
1428
1429 template <int dim, int spacedim>
1430 int
1431 RefineAndCoarsenList<dim, spacedim>::coarsen_callback(
1432 typename internal::p4est::types<dim>::forest *forest,
1433 typename internal::p4est::types<dim>::topidx coarse_cell_index,
1434 typename internal::p4est::types<dim>::quadrant *children[])
1435 {
1436 RefineAndCoarsenList<dim, spacedim> *this_object =
1437 reinterpret_cast<RefineAndCoarsenList<dim, spacedim> *>(
1438 forest->user_pointer);
1439
1440 // if there are no more cells in our list the current cell can't be
1441 // flagged for coarsening
1442 if (this_object->current_coarsen_pointer == this_object->coarsen_list.end())
1443 return 0;
1444
1445 Assert(coarse_cell_index <=
1446 this_object->current_coarsen_pointer->p.which_tree,
1448
1449 // if p4est hasn't yet reached the tree of the next flagged cell the
1450 // current cell can't be flagged for coarsening
1451 if (coarse_cell_index < this_object->current_coarsen_pointer->p.which_tree)
1452 return 0;
1453
1454 // now we're in the right tree in the forest
1455 Assert(coarse_cell_index <=
1456 this_object->current_coarsen_pointer->p.which_tree,
1458
1459 // make sure that the p4est loop over cells hasn't gotten ahead of our own
1460 // pointer
1462 children[0], &*this_object->current_coarsen_pointer) <= 0,
1464
1465 // now, if the p4est cell is one in the list, it is supposed to be
1466 // coarsened
1468 children[0], &*this_object->current_coarsen_pointer))
1469 {
1470 // move current pointer one up
1471 ++this_object->current_coarsen_pointer;
1472
1473 // note that the next 3 cells in our list need to correspond to the
1474 // other siblings of the cell we have just found
1475 for (unsigned int c = 1; c < GeometryInfo<dim>::max_children_per_cell;
1476 ++c)
1477 {
1479 children[c], &*this_object->current_coarsen_pointer),
1481 ++this_object->current_coarsen_pointer;
1482 }
1483
1484 return 1;
1485 }
1486
1487 // p4est cell is not in list
1488 return 0;
1489 }
1490
1491
1492
1499 template <int dim, int spacedim>
1500 class PartitionWeights
1501 {
1502 public:
1508 explicit PartitionWeights(const std::vector<unsigned int> &cell_weights);
1509
1517 static int
1518 cell_weight(typename internal::p4est::types<dim>::forest *forest,
1519 typename internal::p4est::types<dim>::topidx coarse_cell_index,
1520 typename internal::p4est::types<dim>::quadrant *quadrant);
1521
1522 private:
1523 std::vector<unsigned int> cell_weights_list;
1524 std::vector<unsigned int>::const_iterator current_pointer;
1525 };
1526
1527
1528 template <int dim, int spacedim>
1529 PartitionWeights<dim, spacedim>::PartitionWeights(
1530 const std::vector<unsigned int> &cell_weights)
1531 : cell_weights_list(cell_weights)
1532 {
1533 // set the current pointer to the first element of the list, given that
1534 // we will walk through it sequentially
1535 current_pointer = cell_weights_list.begin();
1536 }
1537
1538
1539 template <int dim, int spacedim>
1540 int
1541 PartitionWeights<dim, spacedim>::cell_weight(
1542 typename internal::p4est::types<dim>::forest *forest,
1545 {
1546 // the function gets two additional arguments, but we don't need them
1547 // since we know in which order p4est will walk through the cells
1548 // and have already built our weight lists in this order
1549
1550 PartitionWeights<dim, spacedim> *this_object =
1551 reinterpret_cast<PartitionWeights<dim, spacedim> *>(forest->user_pointer);
1552
1553 Assert(this_object->current_pointer >=
1554 this_object->cell_weights_list.begin(),
1556 Assert(this_object->current_pointer < this_object->cell_weights_list.end(),
1558
1559 // Get the weight, increment the pointer, and return the weight. Also
1560 // make sure that we don't exceed the 'int' data type that p4est uses
1561 // to represent weights
1562 const unsigned int weight = *this_object->current_pointer;
1563 ++this_object->current_pointer;
1564
1565 Assert(weight < static_cast<unsigned int>(std::numeric_limits<int>::max()),
1566 ExcMessage("p4est uses 'signed int' to represent the partition "
1567 "weights for cells. The weight provided here exceeds "
1568 "the maximum value represented as a 'signed int'."));
1569 return static_cast<int>(weight);
1570 }
1571
1572 template <int dim, int spacedim>
1573 using cell_relation_t = typename std::pair<
1574 typename ::Triangulation<dim, spacedim>::cell_iterator,
1575 CellStatus>;
1576
1586 template <int dim, int spacedim>
1587 inline void
1588 add_single_cell_relation(
1589 std::vector<cell_relation_t<dim, spacedim>> &cell_rel,
1590 const typename ::internal::p4est::types<dim>::tree &tree,
1591 const unsigned int idx,
1592 const typename Triangulation<dim, spacedim>::cell_iterator &dealii_cell,
1593 const CellStatus status)
1594 {
1595 const unsigned int local_quadrant_index = tree.quadrants_offset + idx;
1596
1597 // check if we will be writing into valid memory
1598 Assert(local_quadrant_index < cell_rel.size(), ExcInternalError());
1599
1600 // store relation
1601 cell_rel[local_quadrant_index] = std::make_pair(dealii_cell, status);
1602 }
1603
1604
1605
1615 template <int dim, int spacedim>
1616 void
1617 update_cell_relations_recursively(
1618 std::vector<cell_relation_t<dim, spacedim>> &cell_rel,
1619 const typename ::internal::p4est::types<dim>::tree &tree,
1620 const typename Triangulation<dim, spacedim>::cell_iterator &dealii_cell,
1621 const typename ::internal::p4est::types<dim>::quadrant &p4est_cell)
1622 {
1623 // find index of p4est_cell in the quadrants array of the corresponding tree
1624 const int idx = sc_array_bsearch(
1625 const_cast<sc_array_t *>(&tree.quadrants),
1626 &p4est_cell,
1628 if (idx == -1 &&
1630 const_cast<typename ::internal::p4est::types<dim>::tree *>(
1631 &tree),
1632 &p4est_cell) == false))
1633 // this quadrant and none of its children belong to us.
1634 return;
1635
1636 // recurse further if both p4est and dealii still have children
1637 const bool p4est_has_children = (idx == -1);
1638 if (p4est_has_children && dealii_cell->has_children())
1639 {
1640 // recurse further
1641 typename ::internal::p4est::types<dim>::quadrant
1643
1644 for (unsigned int c = 0; c < GeometryInfo<dim>::max_children_per_cell;
1645 ++c)
1646 switch (dim)
1647 {
1648 case 2:
1649 P4EST_QUADRANT_INIT(&p4est_child[c]);
1650 break;
1651 case 3:
1652 P8EST_QUADRANT_INIT(&p4est_child[c]);
1653 break;
1654 default:
1656 }
1657
1659 &p4est_cell, p4est_child);
1660
1661 for (unsigned int c = 0; c < GeometryInfo<dim>::max_children_per_cell;
1662 ++c)
1663 {
1664 update_cell_relations_recursively<dim, spacedim>(
1665 cell_rel, tree, dealii_cell->child(c), p4est_child[c]);
1666 }
1667 }
1668 else if (!p4est_has_children && !dealii_cell->has_children())
1669 {
1670 // this active cell didn't change
1671 // save pair into corresponding position
1672 add_single_cell_relation<dim, spacedim>(
1673 cell_rel, tree, idx, dealii_cell, CellStatus::cell_will_persist);
1674 }
1675 else if (p4est_has_children) // based on the conditions above, we know that
1676 // dealii_cell has no children
1677 {
1678 // this cell got refined in p4est, but the dealii_cell has not yet been
1679 // refined
1680
1681 // this quadrant is not active
1682 // generate its children, and store information in those
1683 typename ::internal::p4est::types<dim>::quadrant
1685 for (unsigned int c = 0; c < GeometryInfo<dim>::max_children_per_cell;
1686 ++c)
1687 switch (dim)
1688 {
1689 case 2:
1690 P4EST_QUADRANT_INIT(&p4est_child[c]);
1691 break;
1692 case 3:
1693 P8EST_QUADRANT_INIT(&p4est_child[c]);
1694 break;
1695 default:
1697 }
1698
1700 &p4est_cell, p4est_child);
1701
1702 // mark first child with CellStatus::cell_will_be_refined and the
1703 // remaining children with CellStatus::cell_invalid, but associate them
1704 // all with the parent cell unpack algorithm will be called only on
1705 // CellStatus::cell_will_be_refined flagged quadrant
1706 int child_idx;
1707 CellStatus cell_status;
1708 for (unsigned int i = 0; i < GeometryInfo<dim>::max_children_per_cell;
1709 ++i)
1710 {
1711 child_idx = sc_array_bsearch(
1712 const_cast<sc_array_t *>(&tree.quadrants),
1713 &p4est_child[i],
1715
1716 cell_status = (i == 0) ? CellStatus::cell_will_be_refined :
1718
1719 add_single_cell_relation<dim, spacedim>(
1720 cell_rel, tree, child_idx, dealii_cell, cell_status);
1721 }
1722 }
1723 else // based on the conditions above, we know that p4est_cell has no
1724 // children, and the dealii_cell does
1725 {
1726 // its children got coarsened into this cell in p4est,
1727 // but the dealii_cell still has its children
1728 add_single_cell_relation<dim, spacedim>(
1729 cell_rel,
1730 tree,
1731 idx,
1732 dealii_cell,
1734 }
1735 }
1736} // namespace
1737
1738
1739
1740namespace parallel
1741{
1742 namespace distributed
1743 {
1744 /*----------------- class Triangulation<dim,spacedim> ---------------*/
1745 template <int dim, int spacedim>
1748 const MPI_Comm mpi_communicator,
1749 const typename ::Triangulation<dim, spacedim>::MeshSmoothing
1750 smooth_grid,
1751 const Settings settings)
1752 : // Do not check for distorted cells.
1753 // For multigrid, we need limit_level_difference_at_vertices
1754 // to make sure the transfer operators only need to consider two levels.
1755 ::parallel::DistributedTriangulationBase<dim, spacedim>(
1756 mpi_communicator,
1757 (settings & construct_multigrid_hierarchy) ?
1758 static_cast<
1759 typename ::Triangulation<dim, spacedim>::MeshSmoothing>(
1760 smooth_grid |
1761 Triangulation<dim, spacedim>::limit_level_difference_at_vertices) :
1762 smooth_grid,
1763 false)
1764 , settings(settings)
1765 , triangulation_has_content(false)
1766 , connectivity(nullptr)
1767 , parallel_forest(nullptr)
1768 {
1769 parallel_ghost = nullptr;
1770 }
1771
1772
1773
1774 template <int dim, int spacedim>
1776 Triangulation<dim, spacedim>::~Triangulation()
1777 {
1778 // virtual functions called in constructors and destructors never use the
1779 // override in a derived class
1780 // for clarity be explicit on which function is called
1781 try
1782 {
1784 }
1785 catch (...)
1786 {}
1787
1788 AssertNothrow(triangulation_has_content == false, ExcInternalError());
1789 AssertNothrow(connectivity == nullptr, ExcInternalError());
1790 AssertNothrow(parallel_forest == nullptr, ExcInternalError());
1791 }
1792
1793
1794
1795 template <int dim, int spacedim>
1797 void Triangulation<dim, spacedim>::create_triangulation(
1798 const std::vector<Point<spacedim>> &vertices,
1799 const std::vector<CellData<dim>> &cells,
1800 const SubCellData &subcelldata)
1801 {
1802 try
1803 {
1805 vertices, cells, subcelldata);
1806 }
1807 catch (
1808 const typename ::Triangulation<dim, spacedim>::DistortedCellList
1809 &)
1810 {
1811 // the underlying triangulation should not be checking for distorted
1812 // cells
1814 }
1815
1816 Assert(
1817 this->all_reference_cells_are_hyper_cube(),
1818 ExcMessage(
1819 "The class parallel::distributed::Triangulation only supports meshes "
1820 "consisting only of hypercube-like cells."));
1821
1822 // note that now we have some content in the p4est objects and call the
1823 // functions that do the actual work (which are dimension dependent, so
1824 // separate)
1825 triangulation_has_content = true;
1826
1827 setup_coarse_cell_to_p4est_tree_permutation();
1828
1829 copy_new_triangulation_to_p4est(std::integral_constant<int, dim>());
1830
1831 try
1832 {
1833 copy_local_forest_to_triangulation();
1834 }
1835 catch (const typename Triangulation<dim>::DistortedCellList &)
1836 {
1837 // the underlying triangulation should not be checking for distorted
1838 // cells
1840 }
1841
1842 this->update_periodic_face_map();
1843 this->update_number_cache();
1844 }
1845
1846
1847
1848 template <int dim, int spacedim>
1850 void Triangulation<dim, spacedim>::create_triangulation(
1851 const TriangulationDescription::Description<dim, spacedim>
1852 & /*construction_data*/)
1853 {
1855 }
1856
1857
1858
1859 template <int dim, int spacedim>
1861 void Triangulation<dim, spacedim>::clear()
1862 {
1863 triangulation_has_content = false;
1864
1865 if (parallel_ghost != nullptr)
1866 {
1868 parallel_ghost);
1869 parallel_ghost = nullptr;
1870 }
1871
1872 if (parallel_forest != nullptr)
1873 {
1875 parallel_forest = nullptr;
1876 }
1877
1878 if (connectivity != nullptr)
1879 {
1881 connectivity);
1882 connectivity = nullptr;
1883 }
1884
1885 coarse_cell_to_p4est_tree_permutation.resize(0);
1886 p4est_tree_to_coarse_cell_permutation.resize(0);
1887
1889
1890 this->update_number_cache();
1891 }
1892
1893
1894
1895 template <int dim, int spacedim>
1897 bool Triangulation<dim, spacedim>::is_multilevel_hierarchy_constructed()
1898 const
1899 {
1900 return settings &
1902 }
1903
1904
1905
1906 template <int dim, int spacedim>
1908 bool Triangulation<dim, spacedim>::are_vertices_communicated_to_p4est()
1909 const
1910 {
1911 return settings &
1913 }
1914
1915
1916
1917 template <int dim, int spacedim>
1919 void Triangulation<dim, spacedim>::execute_transfer(
1920 const typename ::internal::p4est::types<dim>::forest
1921 *parallel_forest,
1922 const typename ::internal::p4est::types<dim>::gloidx
1923 *previous_global_first_quadrant)
1924 {
1925 Assert(this->data_serializer.sizes_fixed_cumulative.size() > 0,
1926 ExcMessage("No data has been packed!"));
1927
1928 // Resize memory according to the data that we will receive.
1929 this->data_serializer.dest_data_fixed.resize(
1930 parallel_forest->local_num_quadrants *
1931 this->data_serializer.sizes_fixed_cumulative.back());
1932
1933 // Execute non-blocking fixed size transfer.
1934 typename ::internal::p4est::types<dim>::transfer_context
1935 *tf_context;
1936 tf_context =
1938 parallel_forest->global_first_quadrant,
1939 previous_global_first_quadrant,
1940 parallel_forest->mpicomm,
1941 0,
1942 this->data_serializer.dest_data_fixed.data(),
1943 this->data_serializer.src_data_fixed.data(),
1944 this->data_serializer.sizes_fixed_cumulative.back());
1945
1946 if (this->data_serializer.variable_size_data_stored)
1947 {
1948 // Resize memory according to the data that we will receive.
1949 this->data_serializer.dest_sizes_variable.resize(
1950 parallel_forest->local_num_quadrants);
1951
1952 // Execute fixed size transfer of data sizes for variable size
1953 // transfer.
1955 parallel_forest->global_first_quadrant,
1956 previous_global_first_quadrant,
1957 parallel_forest->mpicomm,
1958 1,
1959 this->data_serializer.dest_sizes_variable.data(),
1960 this->data_serializer.src_sizes_variable.data(),
1961 sizeof(unsigned int));
1962 }
1963
1965
1966 // Release memory of previously packed data.
1967 this->data_serializer.src_data_fixed.clear();
1968 this->data_serializer.src_data_fixed.shrink_to_fit();
1969
1970 if (this->data_serializer.variable_size_data_stored)
1971 {
1972 // Resize memory according to the data that we will receive.
1973 this->data_serializer.dest_data_variable.resize(
1974 std::accumulate(this->data_serializer.dest_sizes_variable.begin(),
1975 this->data_serializer.dest_sizes_variable.end(),
1976 std::vector<int>::size_type(0)));
1977
1978 // Execute variable size transfer.
1980 parallel_forest->global_first_quadrant,
1981 previous_global_first_quadrant,
1982 parallel_forest->mpicomm,
1983 1,
1984 this->data_serializer.dest_data_variable.data(),
1985 this->data_serializer.dest_sizes_variable.data(),
1986 this->data_serializer.src_data_variable.data(),
1987 this->data_serializer.src_sizes_variable.data());
1988
1989 // Release memory of previously packed data.
1990 this->data_serializer.src_sizes_variable.clear();
1991 this->data_serializer.src_sizes_variable.shrink_to_fit();
1992 this->data_serializer.src_data_variable.clear();
1993 this->data_serializer.src_data_variable.shrink_to_fit();
1994 }
1995 }
1996
1997
1998
1999 template <int dim, int spacedim>
2001 void Triangulation<dim,
2002 spacedim>::setup_coarse_cell_to_p4est_tree_permutation()
2003 {
2004 DynamicSparsityPattern cell_connectivity;
2006 cell_connectivity);
2007 coarse_cell_to_p4est_tree_permutation.resize(this->n_cells(0));
2009 cell_connectivity, coarse_cell_to_p4est_tree_permutation);
2010
2011 p4est_tree_to_coarse_cell_permutation =
2012 Utilities::invert_permutation(coarse_cell_to_p4est_tree_permutation);
2013 }
2014
2015
2016
2017 template <int dim, int spacedim>
2019 void Triangulation<dim, spacedim>::write_mesh_vtk(
2020 const std::string &file_basename) const
2021 {
2022 Assert(parallel_forest != nullptr,
2023 ExcMessage("Can't produce output when no forest is created yet."));
2024
2025 AssertThrow(are_vertices_communicated_to_p4est(),
2026 ExcMessage(
2027 "To use this function the triangulation's flag "
2028 "Settings::communicate_vertices_to_p4est must be set."));
2029
2031 parallel_forest, nullptr, file_basename.c_str());
2032 }
2033
2034
2035
2036 template <int dim, int spacedim>
2038 void Triangulation<dim, spacedim>::save(
2039 const std::string &file_basename) const
2040 {
2041 Assert(
2042 this->cell_attached_data.n_attached_deserialize == 0,
2043 ExcMessage(
2044 "Not all SolutionTransfer objects have been deserialized after the last call to load()."));
2045 Assert(this->n_cells() > 0,
2046 ExcMessage("Can not save() an empty Triangulation."));
2047
2048 const int myrank =
2049 Utilities::MPI::this_mpi_process(this->mpi_communicator);
2050
2051 // signal that serialization is going to happen
2052 this->signals.pre_distributed_save();
2053
2054 if (this->my_subdomain == 0)
2055 {
2056 std::string fname = file_basename + ".info";
2057 std::ofstream f(fname);
2058 f << "version nproc n_attached_fixed_size_objs n_attached_variable_size_objs n_coarse_cells"
2059 << std::endl
2060 << 5 << " "
2061 << Utilities::MPI::n_mpi_processes(this->mpi_communicator) << " "
2062 << this->cell_attached_data.pack_callbacks_fixed.size() << " "
2063 << this->cell_attached_data.pack_callbacks_variable.size() << " "
2064 << this->n_cells(0) << std::endl;
2065 }
2066
2067 // each cell should have been flagged `CellStatus::cell_will_persist`
2068 for ([[maybe_unused]] const auto &cell_rel : this->local_cell_relations)
2069 {
2070 Assert((cell_rel.second == // cell_status
2073 }
2074
2075 // Save cell attached data.
2076 this->save_attached_data(parallel_forest->global_first_quadrant[myrank],
2077 parallel_forest->global_num_quadrants,
2078 file_basename);
2079
2080 ::internal::p4est::functions<dim>::save(file_basename.c_str(),
2081 parallel_forest,
2082 false);
2083
2084 // signal that serialization has finished
2085 this->signals.post_distributed_save();
2086 }
2087
2088
2089
2090 template <int dim, int spacedim>
2092 void Triangulation<dim, spacedim>::load(const std::string &file_basename)
2093 {
2094 Assert(
2095 this->n_cells() > 0,
2096 ExcMessage(
2097 "load() only works if the Triangulation already contains a coarse mesh!"));
2098 Assert(
2099 this->n_levels() == 1,
2100 ExcMessage(
2101 "Triangulation may only contain coarse cells when calling load()."));
2102
2103 const int myrank =
2104 Utilities::MPI::this_mpi_process(this->mpi_communicator);
2105
2106 // signal that de-serialization is going to happen
2107 this->signals.pre_distributed_load();
2108
2109 if (parallel_ghost != nullptr)
2110 {
2112 parallel_ghost);
2113 parallel_ghost = nullptr;
2114 }
2116 parallel_forest = nullptr;
2118 connectivity);
2119 connectivity = nullptr;
2120
2121 unsigned int version, numcpus, attached_count_fixed,
2122 attached_count_variable, n_coarse_cells;
2123 {
2124 std::string fname = std::string(file_basename) + ".info";
2125 std::ifstream f(fname);
2126 AssertThrow(f.fail() == false, ExcIO());
2127 std::string firstline;
2128 getline(f, firstline); // skip first line
2129 f >> version >> numcpus >> attached_count_fixed >>
2130 attached_count_variable >> n_coarse_cells;
2131 }
2132
2133 AssertThrow(version == 5,
2134 ExcMessage("Incompatible version found in .info file."));
2135 Assert(this->n_cells(0) == n_coarse_cells,
2136 ExcMessage("Number of coarse cells differ!"));
2137
2138 // clear all of the callback data, as explained in the documentation of
2139 // register_data_attach()
2140 this->cell_attached_data.n_attached_data_sets = 0;
2141 this->cell_attached_data.n_attached_deserialize =
2142 attached_count_fixed + attached_count_variable;
2143
2145 file_basename.c_str(),
2146 this->mpi_communicator,
2147 0,
2148 0,
2149 1,
2150 0,
2151 this,
2152 &connectivity);
2153
2154 // We partition the p4est mesh that it conforms to the requirements of the
2155 // deal.II mesh, i.e., partition for coarsening.
2156 // This function call is optional.
2158 parallel_forest,
2159 /* prepare coarsening */ 1,
2160 /* weight_callback */ nullptr);
2161
2162 try
2163 {
2164 copy_local_forest_to_triangulation();
2165 }
2166 catch (const typename Triangulation<dim>::DistortedCellList &)
2167 {
2168 // the underlying triangulation should not be checking for distorted
2169 // cells
2171 }
2172
2173 // Load attached cell data, if any was stored.
2174 this->load_attached_data(parallel_forest->global_first_quadrant[myrank],
2175 parallel_forest->global_num_quadrants,
2176 parallel_forest->local_num_quadrants,
2177 file_basename,
2178 attached_count_fixed,
2179 attached_count_variable);
2180
2181 // signal that de-serialization is finished
2182 this->signals.post_distributed_load();
2183
2184 this->update_periodic_face_map();
2185 this->update_number_cache();
2186 }
2187
2188
2189
2190 template <int dim, int spacedim>
2192 void Triangulation<dim, spacedim>::load(
2193 const typename ::internal::p4est::types<dim>::forest *forest)
2194 {
2195 Assert(this->n_cells() > 0,
2196 ExcMessage(
2197 "load() only works if the Triangulation already contains "
2198 "a coarse mesh!"));
2199 Assert(this->n_cells() == forest->trees->elem_count,
2200 ExcMessage(
2201 "Coarse mesh of the Triangulation does not match the one "
2202 "of the provided forest!"));
2203
2204 // clear the old forest
2205 if (parallel_ghost != nullptr)
2206 {
2208 parallel_ghost);
2209 parallel_ghost = nullptr;
2210 }
2212 parallel_forest = nullptr;
2213
2214 // note: we can keep the connectivity, since the coarse grid does not
2215 // change
2216
2217 // create deep copy of the new forest
2218 typename ::internal::p4est::types<dim>::forest *temp =
2219 const_cast<typename ::internal::p4est::types<dim>::forest *>(
2220 forest);
2221 parallel_forest =
2223 parallel_forest->connectivity = connectivity;
2224 parallel_forest->user_pointer = this;
2225
2226 try
2227 {
2228 copy_local_forest_to_triangulation();
2229 }
2230 catch (const typename Triangulation<dim>::DistortedCellList &)
2231 {
2232 // the underlying triangulation should not be checking for distorted
2233 // cells
2235 }
2236
2237 this->update_periodic_face_map();
2238 this->update_number_cache();
2239 }
2240
2241
2242
2243 template <int dim, int spacedim>
2245 unsigned int Triangulation<dim, spacedim>::get_checksum() const
2246 {
2247 Assert(parallel_forest != nullptr,
2248 ExcMessage(
2249 "Can't produce a check sum when no forest is created yet."));
2250
2251 auto checksum =
2253
2254# if !DEAL_II_P4EST_VERSION_GTE(2, 8, 6, 0)
2255 /*
2256 * p4est prior to 2.8.6 returns the proper checksum only on rank 0
2257 * and simply "0" on all other ranks. This is not really what we
2258 * want, thus broadcast the correct value to all other ranks:
2259 */
2260 checksum = Utilities::MPI::broadcast(this->mpi_communicator,
2261 checksum,
2262 /*root_process*/ 0);
2263# endif
2264
2265 return checksum;
2266 }
2267
2268
2269
2270 template <int dim, int spacedim>
2272 const typename ::internal::p4est::types<dim>::forest
2274 {
2275 Assert(parallel_forest != nullptr,
2276 ExcMessage("The forest has not been allocated yet."));
2277 return parallel_forest;
2278 }
2279
2280
2281
2282 template <int dim, int spacedim>
2284 typename ::internal::p4est::types<dim>::tree
2286 const int dealii_coarse_cell_index) const
2287 {
2288 const unsigned int tree_index =
2289 coarse_cell_to_p4est_tree_permutation[dealii_coarse_cell_index];
2290 typename ::internal::p4est::types<dim>::tree *tree =
2291 static_cast<typename ::internal::p4est::types<dim>::tree *>(
2292 sc_array_index(parallel_forest->trees, tree_index));
2293
2294 return tree;
2295 }
2296
2297
2298
2299 // Note: this has been added here to prevent that these functions
2300 // appear in the Doxygen documentation of ::Triangulation
2301# ifndef DOXYGEN
2302
2303 template <>
2304 void
2306 std::integral_constant<int, 2>)
2307 {
2308 const unsigned int dim = 2, spacedim = 2;
2309 Assert(this->n_cells(0) > 0, ExcInternalError());
2310 Assert(this->n_levels() == 1, ExcInternalError());
2311
2312 // data structures that counts how many cells touch each vertex
2313 // (vertex_touch_count), and which cells touch a given vertex (together
2314 // with the local numbering of that vertex within the cells that touch
2315 // it)
2316 std::vector<unsigned int> vertex_touch_count;
2317 std::vector<
2318 std::list<std::pair<Triangulation<dim, spacedim>::active_cell_iterator,
2319 unsigned int>>>
2320 vertex_to_cell;
2321 get_vertex_to_cell_mappings(*this, vertex_touch_count, vertex_to_cell);
2322 const ::internal::p4est::types<2>::locidx num_vtt =
2323 std::accumulate(vertex_touch_count.begin(),
2324 vertex_touch_count.end(),
2325 0u);
2326
2327 // now create a connectivity object with the right sizes for all
2328 // arrays. set vertex information only in debug mode (saves a few bytes
2329 // in optimized mode)
2330 const bool set_vertex_info = this->are_vertices_communicated_to_p4est();
2331
2333 (set_vertex_info == true ? this->n_vertices() : 0),
2334 this->n_cells(0),
2335 this->n_vertices(),
2336 num_vtt);
2337
2338 set_vertex_and_cell_info(*this,
2339 vertex_touch_count,
2340 vertex_to_cell,
2341 coarse_cell_to_p4est_tree_permutation,
2342 set_vertex_info,
2343 connectivity);
2344
2345 Assert(p4est_connectivity_is_valid(connectivity) == 1,
2347
2348 // now create a forest out of the connectivity data structure
2350 this->mpi_communicator,
2351 connectivity,
2352 /* minimum initial number of quadrants per tree */ 0,
2353 /* minimum level of upfront refinement */ 0,
2354 /* use uniform upfront refinement */ 1,
2355 /* user_data_size = */ 0,
2356 /* user_data_constructor = */ nullptr,
2357 /* user_pointer */ this);
2358 }
2359
2360
2361
2362 // TODO: This is a verbatim copy of the 2,2 case. However, we can't just
2363 // specialize the dim template argument, but let spacedim open
2364 template <>
2365 void
2367 std::integral_constant<int, 2>)
2368 {
2369 const unsigned int dim = 2, spacedim = 3;
2370 Assert(this->n_cells(0) > 0, ExcInternalError());
2371 Assert(this->n_levels() == 1, ExcInternalError());
2372
2373 // data structures that counts how many cells touch each vertex
2374 // (vertex_touch_count), and which cells touch a given vertex (together
2375 // with the local numbering of that vertex within the cells that touch
2376 // it)
2377 std::vector<unsigned int> vertex_touch_count;
2378 std::vector<
2379 std::list<std::pair<Triangulation<dim, spacedim>::active_cell_iterator,
2380 unsigned int>>>
2381 vertex_to_cell;
2382 get_vertex_to_cell_mappings(*this, vertex_touch_count, vertex_to_cell);
2383 const ::internal::p4est::types<2>::locidx num_vtt =
2384 std::accumulate(vertex_touch_count.begin(),
2385 vertex_touch_count.end(),
2386 0u);
2387
2388 // now create a connectivity object with the right sizes for all
2389 // arrays. set vertex information only in debug mode (saves a few bytes
2390 // in optimized mode)
2391 const bool set_vertex_info = this->are_vertices_communicated_to_p4est();
2392
2394 (set_vertex_info == true ? this->n_vertices() : 0),
2395 this->n_cells(0),
2396 this->n_vertices(),
2397 num_vtt);
2398
2399 set_vertex_and_cell_info(*this,
2400 vertex_touch_count,
2401 vertex_to_cell,
2402 coarse_cell_to_p4est_tree_permutation,
2403 set_vertex_info,
2404 connectivity);
2405
2406 Assert(p4est_connectivity_is_valid(connectivity) == 1,
2408
2409 // now create a forest out of the connectivity data structure
2411 this->mpi_communicator,
2412 connectivity,
2413 /* minimum initial number of quadrants per tree */ 0,
2414 /* minimum level of upfront refinement */ 0,
2415 /* use uniform upfront refinement */ 1,
2416 /* user_data_size = */ 0,
2417 /* user_data_constructor = */ nullptr,
2418 /* user_pointer */ this);
2419 }
2420
2421
2422
2423 template <>
2424 void
2426 std::integral_constant<int, 3>)
2427 {
2428 const int dim = 3, spacedim = 3;
2429 Assert(this->n_cells(0) > 0, ExcInternalError());
2430 Assert(this->n_levels() == 1, ExcInternalError());
2431
2432 // data structures that counts how many cells touch each vertex
2433 // (vertex_touch_count), and which cells touch a given vertex (together
2434 // with the local numbering of that vertex within the cells that touch
2435 // it)
2436 std::vector<unsigned int> vertex_touch_count;
2437 std::vector<std::list<
2438 std::pair<Triangulation<3>::active_cell_iterator, unsigned int>>>
2439 vertex_to_cell;
2440 get_vertex_to_cell_mappings(*this, vertex_touch_count, vertex_to_cell);
2441 const ::internal::p4est::types<2>::locidx num_vtt =
2442 std::accumulate(vertex_touch_count.begin(),
2443 vertex_touch_count.end(),
2444 0u);
2445
2446 std::vector<unsigned int> edge_touch_count;
2447 std::vector<std::list<
2448 std::pair<Triangulation<3>::active_cell_iterator, unsigned int>>>
2449 edge_to_cell;
2450 get_edge_to_cell_mappings(*this, edge_touch_count, edge_to_cell);
2451 const ::internal::p4est::types<2>::locidx num_ett =
2452 std::accumulate(edge_touch_count.begin(), edge_touch_count.end(), 0u);
2453
2454 // now create a connectivity object with the right sizes for all arrays
2455 const bool set_vertex_info = this->are_vertices_communicated_to_p4est();
2456
2458 (set_vertex_info == true ? this->n_vertices() : 0),
2459 this->n_cells(0),
2460 this->n_active_lines(),
2461 num_ett,
2462 this->n_vertices(),
2463 num_vtt);
2464
2465 set_vertex_and_cell_info(*this,
2466 vertex_touch_count,
2467 vertex_to_cell,
2468 coarse_cell_to_p4est_tree_permutation,
2469 set_vertex_info,
2470 connectivity);
2471
2472 // next to tree-to-edge
2473 // data. note that in p4est lines
2474 // are ordered as follows
2475 // *---3---* *---3---*
2476 // /| | / /|
2477 // 6 | 11 6 7 11
2478 // / 10 | / / |
2479 // * | | *---2---* |
2480 // | *---1---* | | *
2481 // | / / | 9 /
2482 // 8 4 5 8 | 5
2483 // |/ / | |/
2484 // *---0---* *---0---*
2485 // whereas in deal.II they are like this:
2486 // *---7---* *---7---*
2487 // /| | / /|
2488 // 4 | 11 4 5 11
2489 // / 10 | / / |
2490 // * | | *---6---* |
2491 // | *---3---* | | *
2492 // | / / | 9 /
2493 // 8 0 1 8 | 1
2494 // |/ / | |/
2495 // *---2---* *---2---*
2496
2497 const unsigned int deal_to_p4est_line_index[12] = {
2498 4, 5, 0, 1, 6, 7, 2, 3, 8, 9, 10, 11};
2499
2501 this->begin_active();
2502 cell != this->end();
2503 ++cell)
2504 {
2505 const unsigned int index =
2506 coarse_cell_to_p4est_tree_permutation[cell->index()];
2507 for (unsigned int e = 0; e < GeometryInfo<3>::lines_per_cell; ++e)
2508 connectivity->tree_to_edge[index * GeometryInfo<3>::lines_per_cell +
2509 deal_to_p4est_line_index[e]] =
2510 cell->line(e)->index();
2511 }
2512
2513 // now also set edge-to-tree
2514 // information
2515 connectivity->ett_offset[0] = 0;
2516 std::partial_sum(edge_touch_count.begin(),
2517 edge_touch_count.end(),
2518 &connectivity->ett_offset[1]);
2519
2520 Assert(connectivity->ett_offset[this->n_active_lines()] == num_ett,
2522
2523 for (unsigned int v = 0; v < this->n_active_lines(); ++v)
2524 {
2525 Assert(edge_to_cell[v].size() == edge_touch_count[v],
2527
2528 std::list<
2529 std::pair<Triangulation<dim, spacedim>::active_cell_iterator,
2530 unsigned int>>::const_iterator p =
2531 edge_to_cell[v].begin();
2532 for (unsigned int c = 0; c < edge_touch_count[v]; ++c, ++p)
2533 {
2534 connectivity->edge_to_tree[connectivity->ett_offset[v] + c] =
2535 coarse_cell_to_p4est_tree_permutation[p->first->index()];
2536 connectivity->edge_to_edge[connectivity->ett_offset[v] + c] =
2537 deal_to_p4est_line_index[p->second];
2538 }
2539 }
2540
2541 Assert(p8est_connectivity_is_valid(connectivity) == 1,
2543
2544 // now create a forest out of the connectivity data structure
2546 this->mpi_communicator,
2547 connectivity,
2548 /* minimum initial number of quadrants per tree */ 0,
2549 /* minimum level of upfront refinement */ 0,
2550 /* use uniform upfront refinement */ 1,
2551 /* user_data_size = */ 0,
2552 /* user_data_constructor = */ nullptr,
2553 /* user_pointer */ this);
2554 }
2555# endif
2556
2557
2558
2559 namespace
2560 {
2561 // ensures the 2:1 mesh balance for periodic boundary conditions in the
2562 // artificial cell layer (the active cells are taken care of by p4est)
2563 template <int dim, int spacedim>
2564 bool
2565 enforce_mesh_balance_over_periodic_boundaries(
2567 {
2568 if (tria.get_periodic_face_map().empty())
2569 return false;
2570
2571 std::vector<bool> flags_before[2];
2572 tria.save_coarsen_flags(flags_before[0]);
2573 tria.save_refine_flags(flags_before[1]);
2574
2575 std::vector<unsigned int> topological_vertex_numbering(
2576 tria.n_vertices());
2577 for (unsigned int i = 0; i < topological_vertex_numbering.size(); ++i)
2578 topological_vertex_numbering[i] = i;
2579 // combine vertices that have different locations (and thus, different
2580 // vertex_index) but represent the same topological entity over
2581 // periodic boundaries. The vector topological_vertex_numbering
2582 // contains a linear map from 0 to n_vertices at input and at output
2583 // relates periodic vertices with only one vertex index. The output is
2584 // used to always identify the same vertex according to the
2585 // periodicity, e.g. when finding the maximum cell level around a
2586 // vertex.
2587 //
2588 // Example: On a 3d cell with vertices numbered from 0 to 7 and
2589 // periodic boundary conditions in x direction, the vector
2590 // topological_vertex_numbering will contain the numbers
2591 // {0,0,2,2,4,4,6,6} (because the vertex pairs {0,1}, {2,3}, {4,5},
2592 // {6,7} belong together, respectively). If periodicity is set in x
2593 // and z direction, the output is {0,0,2,2,0,0,2,2}, and if
2594 // periodicity is in all directions, the output is simply
2595 // {0,0,0,0,0,0,0,0}.
2596 using cell_iterator =
2598 for (const auto &it : tria.get_periodic_face_map())
2599 {
2600 const cell_iterator &cell_1 = it.first.first;
2601 const unsigned int face_no_1 = it.first.second;
2602 const cell_iterator &cell_2 = it.second.first.first;
2603 const unsigned int face_no_2 = it.second.first.second;
2604 const unsigned char combined_orientation = it.second.second;
2605
2606 if (cell_1->level() == cell_2->level())
2607 {
2608 for (const unsigned int v :
2609 cell_1->face(face_no_1)->vertex_indices())
2610 {
2611 // take possible non-standard orientation of face on
2612 // cell[0] into account
2613 const unsigned int vface1 =
2614 cell_1->reference_cell().standard_to_real_face_vertex(
2615 v, face_no_1, combined_orientation);
2616 const unsigned int vi1 =
2617 topological_vertex_numbering[cell_1->face(face_no_1)
2618 ->vertex_index(vface1)];
2619 const unsigned int vi2 =
2620 topological_vertex_numbering[cell_2->face(face_no_2)
2621 ->vertex_index(v)];
2622 const unsigned int min_index = std::min(vi1, vi2);
2623 topological_vertex_numbering[cell_1->face(face_no_1)
2624 ->vertex_index(vface1)] =
2625 topological_vertex_numbering[cell_2->face(face_no_2)
2626 ->vertex_index(v)] =
2627 min_index;
2628 }
2629 }
2630 }
2631
2632# ifdef DEBUG
2633 // There must not be any chains!
2634 for (unsigned int i = 0; i < topological_vertex_numbering.size(); ++i)
2635 {
2636 const unsigned int j = topological_vertex_numbering[i];
2637 Assert(j == i || topological_vertex_numbering[j] == j,
2638 ExcMessage("Got inconclusive constraints with chain: " +
2639 std::to_string(i) + " vs " + std::to_string(j) +
2640 " which should be equal to " +
2641 std::to_string(topological_vertex_numbering[j])));
2642 }
2643# endif
2644
2645
2646 // this code is replicated from grid/tria.cc but using an indirection
2647 // for periodic boundary conditions
2648 bool continue_iterating = true;
2649 std::vector<int> vertex_level(tria.n_vertices());
2650 while (continue_iterating)
2651 {
2652 // store highest level one of the cells adjacent to a vertex
2653 // belongs to
2654 std::fill(vertex_level.begin(), vertex_level.end(), 0);
2656 cell = tria.begin_active(),
2657 endc = tria.end();
2658 for (; cell != endc; ++cell)
2659 {
2660 if (cell->refine_flag_set())
2661 for (const unsigned int vertex :
2663 vertex_level[topological_vertex_numbering
2664 [cell->vertex_index(vertex)]] =
2665 std::max(vertex_level[topological_vertex_numbering
2666 [cell->vertex_index(vertex)]],
2667 cell->level() + 1);
2668 else if (!cell->coarsen_flag_set())
2669 for (const unsigned int vertex :
2671 vertex_level[topological_vertex_numbering
2672 [cell->vertex_index(vertex)]] =
2673 std::max(vertex_level[topological_vertex_numbering
2674 [cell->vertex_index(vertex)]],
2675 cell->level());
2676 else
2677 {
2678 // if coarsen flag is set then tentatively assume
2679 // that the cell will be coarsened. this isn't
2680 // always true (the coarsen flag could be removed
2681 // again) and so we may make an error here. we try
2682 // to correct this by iterating over the entire
2683 // process until we are converged
2684 Assert(cell->coarsen_flag_set(), ExcInternalError());
2685 for (const unsigned int vertex :
2687 vertex_level[topological_vertex_numbering
2688 [cell->vertex_index(vertex)]] =
2689 std::max(vertex_level[topological_vertex_numbering
2690 [cell->vertex_index(vertex)]],
2691 cell->level() - 1);
2692 }
2693 }
2694
2695 continue_iterating = false;
2696
2697 // loop over all cells in reverse order. do so because we
2698 // can then update the vertex levels on the adjacent
2699 // vertices and maybe already flag additional cells in this
2700 // loop
2701 //
2702 // note that not only may we have to add additional
2703 // refinement flags, but we will also have to remove
2704 // coarsening flags on cells adjacent to vertices that will
2705 // see refinement
2706 for (cell = tria.last_active(); cell != endc; --cell)
2707 if (cell->refine_flag_set() == false)
2708 {
2709 for (const unsigned int vertex :
2711 if (vertex_level[topological_vertex_numbering
2712 [cell->vertex_index(vertex)]] >=
2713 cell->level() + 1)
2714 {
2715 // remove coarsen flag...
2716 cell->clear_coarsen_flag();
2717
2718 // ...and if necessary also refine the current
2719 // cell, at the same time updating the level
2720 // information about vertices
2721 if (vertex_level[topological_vertex_numbering
2722 [cell->vertex_index(vertex)]] >
2723 cell->level() + 1)
2724 {
2725 cell->set_refine_flag();
2726 continue_iterating = true;
2727
2728 for (const unsigned int v :
2730 vertex_level[topological_vertex_numbering
2731 [cell->vertex_index(v)]] =
2732 std::max(
2733 vertex_level[topological_vertex_numbering
2734 [cell->vertex_index(v)]],
2735 cell->level() + 1);
2736 }
2737
2738 // continue and see whether we may, for example,
2739 // go into the inner 'if' above based on a
2740 // different vertex
2741 }
2742 }
2743
2744 // clear coarsen flag if not all children were marked
2745 for (const auto &cell : tria.cell_iterators())
2746 {
2747 // nothing to do if we are already on the finest level
2748 if (cell->is_active())
2749 continue;
2750
2751 const unsigned int n_children = cell->n_children();
2752 unsigned int flagged_children = 0;
2753 for (unsigned int child = 0; child < n_children; ++child)
2754 if (cell->child(child)->is_active() &&
2755 cell->child(child)->coarsen_flag_set())
2756 ++flagged_children;
2757
2758 // if not all children were flagged for coarsening, remove
2759 // coarsen flags
2760 if (flagged_children < n_children)
2761 for (unsigned int child = 0; child < n_children; ++child)
2762 if (cell->child(child)->is_active())
2763 cell->child(child)->clear_coarsen_flag();
2764 }
2765 }
2766 std::vector<bool> flags_after[2];
2767 tria.save_coarsen_flags(flags_after[0]);
2768 tria.save_refine_flags(flags_after[1]);
2769 return ((flags_before[0] != flags_after[0]) ||
2770 (flags_before[1] != flags_after[1]));
2771 }
2772 } // namespace
2773
2774
2775
2776 template <int dim, int spacedim>
2778 bool Triangulation<dim, spacedim>::prepare_coarsening_and_refinement()
2779 {
2780 // First exchange coarsen/refinement flags on ghost cells. After this
2781 // collective communication call all flags on ghost cells match the
2782 // flags set by the user on the owning rank.
2785
2786 // Now we can call the sequential version to apply mesh smoothing and
2787 // other modifications:
2788 const bool any_changes = this->::Triangulation<dim, spacedim>::
2790 return any_changes;
2791 }
2792
2793
2794
2795 template <int dim, int spacedim>
2797 void Triangulation<dim, spacedim>::copy_local_forest_to_triangulation()
2798 {
2799 // Disable mesh smoothing for recreating the deal.II triangulation,
2800 // otherwise we might not be able to reproduce the p4est mesh
2801 // exactly. We restore the original smoothing at the end of this
2802 // function. Note that the smoothing flag is used in the normal
2803 // refinement process.
2804 typename Triangulation<dim, spacedim>::MeshSmoothing save_smooth =
2805 this->smooth_grid;
2806
2807 // We will refine manually to match the p4est further down, which
2808 // obeys a level difference of 2 at each vertex (see the balance call
2809 // to p4est). We can disable this here so we store fewer artificial
2810 // cells (in some cases).
2811 // For geometric multigrid it turns out that
2812 // we will miss level cells at shared vertices if we ignore this.
2813 // See tests/mpi/mg_06. In particular, the flag is still necessary
2814 // even though we force it for the original smooth_grid in the
2815 // constructor.
2816 if (settings & construct_multigrid_hierarchy)
2817 this->smooth_grid =
2818 ::Triangulation<dim,
2819 spacedim>::limit_level_difference_at_vertices;
2820 else
2821 this->smooth_grid = ::Triangulation<dim, spacedim>::none;
2822
2823 bool mesh_changed = false;
2824
2825 // Remove all deal.II refinements. Note that we could skip this and
2826 // start from our current state, because the algorithm later coarsens as
2827 // necessary. This has the advantage of being faster when large parts
2828 // of the local partition changes (likely) and gives a deterministic
2829 // ordering of the cells (useful for snapshot/resume).
2830 // TODO: is there a more efficient way to do this?
2831 if (settings & mesh_reconstruction_after_repartitioning)
2832 while (this->n_levels() > 1)
2833 {
2834 // Instead of marking all active cells, we slice off the finest
2835 // level, one level at a time. This takes the same number of
2836 // iterations but solves an issue where not all cells on a
2837 // periodic boundary are indeed coarsened and we run into an
2838 // irrelevant Assert() in update_periodic_face_map().
2839 for (const auto &cell :
2840 this->active_cell_iterators_on_level(this->n_levels() - 1))
2841 {
2842 cell->set_coarsen_flag();
2843 }
2844 try
2845 {
2848 }
2849 catch (
2851 {
2852 // the underlying triangulation should not be checking for
2853 // distorted cells
2855 }
2856 }
2857
2858
2859 // query p4est for the ghost cells
2860 if (parallel_ghost != nullptr)
2861 {
2863 parallel_ghost);
2864 parallel_ghost = nullptr;
2865 }
2867 parallel_forest,
2868 (dim == 2 ? typename ::internal::p4est::types<dim>::balance_type(
2869 P4EST_CONNECT_CORNER) :
2870 typename ::internal::p4est::types<dim>::balance_type(
2871 P8EST_CONNECT_CORNER)));
2872
2873 Assert(parallel_ghost, ExcInternalError());
2874
2875
2876 // set all cells to artificial. we will later set it to the correct
2877 // subdomain in match_tree_recursively
2878 for (const auto &cell : this->cell_iterators_on_level(0))
2879 cell->recursively_set_subdomain_id(numbers::artificial_subdomain_id);
2880
2881 do
2882 {
2883 for (const auto &cell : this->cell_iterators_on_level(0))
2884 {
2885 // if this processor stores no part of the forest that comes out
2886 // of this coarse grid cell, then we need to delete all children
2887 // of this cell (the coarse grid cell remains)
2888 if (tree_exists_locally<dim, spacedim>(
2889 parallel_forest,
2890 coarse_cell_to_p4est_tree_permutation[cell->index()]) ==
2891 false)
2892 {
2893 delete_all_children<dim, spacedim>(cell);
2894 if (cell->is_active())
2895 cell->set_subdomain_id(numbers::artificial_subdomain_id);
2896 }
2897
2898 else
2899 {
2900 // this processor stores at least a part of the tree that
2901 // comes out of this cell.
2902
2903 typename ::internal::p4est::types<dim>::quadrant
2904 p4est_coarse_cell;
2905 typename ::internal::p4est::types<dim>::tree *tree =
2906 init_tree(cell->index());
2907
2908 ::internal::p4est::init_coarse_quadrant<dim>(
2909 p4est_coarse_cell);
2910
2911 match_tree_recursively<dim, spacedim>(*tree,
2912 cell,
2913 p4est_coarse_cell,
2914 *parallel_forest,
2915 this->my_subdomain);
2916 }
2917 }
2918
2919 // check mesh for ghost cells, refine as necessary. iterate over
2920 // every ghostquadrant, find corresponding deal coarsecell and
2921 // recurse.
2922 typename ::internal::p4est::types<dim>::quadrant *quadr;
2923 types::subdomain_id ghost_owner = 0;
2924 typename ::internal::p4est::types<dim>::topidx ghost_tree = 0;
2925
2926 for (unsigned int g_idx = 0;
2927 g_idx < parallel_ghost->ghosts.elem_count;
2928 ++g_idx)
2929 {
2930 while (g_idx >= static_cast<unsigned int>(
2931 parallel_ghost->proc_offsets[ghost_owner + 1]))
2932 ++ghost_owner;
2933 while (g_idx >= static_cast<unsigned int>(
2934 parallel_ghost->tree_offsets[ghost_tree + 1]))
2935 ++ghost_tree;
2936
2937 quadr = static_cast<
2938 typename ::internal::p4est::types<dim>::quadrant *>(
2939 sc_array_index(&parallel_ghost->ghosts, g_idx));
2940
2941 unsigned int coarse_cell_index =
2942 p4est_tree_to_coarse_cell_permutation[ghost_tree];
2943
2944 match_quadrant<dim, spacedim>(this,
2945 coarse_cell_index,
2946 *quadr,
2947 ghost_owner);
2948 }
2949
2950 // Fix all the flags to make sure we have a consistent local
2951 // mesh. For some reason periodic boundaries involving artificial
2952 // cells are not obeying the 2:1 ratio that we require (and that is
2953 // enforced by p4est between active cells). So, here we will loop
2954 // refining across periodic boundaries until 2:1 is satisfied. Note
2955 // that we are using the base class (sequential) prepare and execute
2956 // calls here, not involving communication, because we are only
2957 // trying to recreate a local triangulation from the p4est data.
2958 {
2959 bool mesh_changed = true;
2960 unsigned int loop_counter = 0;
2961
2962 do
2963 {
2966
2967 this->update_periodic_face_map();
2968
2969 mesh_changed =
2970 enforce_mesh_balance_over_periodic_boundaries(*this);
2971
2972 // We can't be sure that we won't run into a situation where we
2973 // can not reconcile mesh smoothing and balancing of periodic
2974 // faces. As we don't know what else to do, at least abort with
2975 // an error message.
2976 ++loop_counter;
2977
2979 loop_counter < 32,
2980 ExcMessage(
2981 "Infinite loop in "
2982 "parallel::distributed::Triangulation::copy_local_forest_to_triangulation() "
2983 "for periodic boundaries detected. Aborting."));
2984 }
2985 while (mesh_changed);
2986 }
2987
2988 // see if any flags are still set
2989 mesh_changed =
2990 std::any_of(this->begin_active(),
2991 active_cell_iterator{this->end()},
2992 [](const CellAccessor<dim, spacedim> &cell) {
2993 return cell.refine_flag_set() ||
2994 cell.coarsen_flag_set();
2995 });
2996
2997 // actually do the refinement to change the local mesh by
2998 // calling the base class refinement function directly
2999 try
3000 {
3003 }
3004 catch (
3006 {
3007 // the underlying triangulation should not be checking for
3008 // distorted cells
3010 }
3011 }
3012 while (mesh_changed);
3013
3014# ifdef DEBUG
3015 // check if correct number of ghosts is created
3016 unsigned int num_ghosts = 0;
3017
3018 for (const auto &cell : this->active_cell_iterators())
3019 {
3020 if (cell->subdomain_id() != this->my_subdomain &&
3021 cell->subdomain_id() != numbers::artificial_subdomain_id)
3022 ++num_ghosts;
3023 }
3024
3025 Assert(num_ghosts == parallel_ghost->ghosts.elem_count,
3027# endif
3028
3029
3030
3031 // fill level_subdomain_ids for geometric multigrid
3032 // the level ownership of a cell is defined as the owner if the cell is
3033 // active or as the owner of child(0) we need this information for all
3034 // our ancestors and the same-level neighbors of our own cells (=level
3035 // ghosts)
3036 if (settings & construct_multigrid_hierarchy)
3037 {
3038 // step 1: We set our own ids all the way down and all the others to
3039 // -1. Note that we do not fill other cells we could figure out the
3040 // same way, because we might accidentally set an id for a cell that
3041 // is not a ghost cell.
3042 for (unsigned int lvl = this->n_levels(); lvl > 0;)
3043 {
3044 --lvl;
3045 for (const auto &cell : this->cell_iterators_on_level(lvl))
3046 {
3047 if ((cell->is_active() &&
3048 cell->subdomain_id() ==
3049 this->locally_owned_subdomain()) ||
3050 (cell->has_children() &&
3051 cell->child(0)->level_subdomain_id() ==
3052 this->locally_owned_subdomain()))
3053 cell->set_level_subdomain_id(
3054 this->locally_owned_subdomain());
3055 else
3056 {
3057 // not our cell
3058 cell->set_level_subdomain_id(
3060 }
3061 }
3062 }
3063
3064 // step 2: make sure all the neighbors to our level_cells exist.
3065 // Need to look up in p4est...
3066 std::vector<std::vector<bool>> marked_vertices(this->n_levels());
3067 for (unsigned int lvl = 0; lvl < this->n_levels(); ++lvl)
3068 marked_vertices[lvl] = mark_locally_active_vertices_on_level(lvl);
3069
3070 for (const auto &cell : this->cell_iterators_on_level(0))
3071 {
3072 typename ::internal::p4est::types<dim>::quadrant
3073 p4est_coarse_cell;
3074 const unsigned int tree_index =
3075 coarse_cell_to_p4est_tree_permutation[cell->index()];
3076 typename ::internal::p4est::types<dim>::tree *tree =
3077 init_tree(cell->index());
3078
3079 ::internal::p4est::init_coarse_quadrant<dim>(
3080 p4est_coarse_cell);
3081
3082 determine_level_subdomain_id_recursively<dim, spacedim>(
3083 *tree,
3084 tree_index,
3085 cell,
3086 p4est_coarse_cell,
3087 *parallel_forest,
3088 this->my_subdomain,
3089 marked_vertices);
3090 }
3091
3092 // step 3: make sure we have the parent of our level cells
3093 for (unsigned int lvl = this->n_levels(); lvl > 0;)
3094 {
3095 --lvl;
3096 for (const auto &cell : this->cell_iterators_on_level(lvl))
3097 {
3098 if (cell->has_children())
3099 for (unsigned int c = 0;
3100 c < GeometryInfo<dim>::max_children_per_cell;
3101 ++c)
3102 {
3103 if (cell->child(c)->level_subdomain_id() ==
3104 this->locally_owned_subdomain())
3105 {
3106 // at least one of the children belongs to us, so
3107 // make sure we set the level subdomain id
3108 const types::subdomain_id mark =
3109 cell->child(0)->level_subdomain_id();
3111 ExcInternalError()); // we should know the
3112 // child(0)
3113 cell->set_level_subdomain_id(mark);
3114 break;
3115 }
3116 }
3117 }
3118 }
3119 }
3120
3121
3122
3123# ifdef DEBUG
3124 // check that our local copy has exactly as many cells as the p4est
3125 // original (at least if we are on only one processor); for parallel
3126 // computations, we want to check that we have at least as many as p4est
3127 // stores locally (in the future we should check that we have exactly as
3128 // many non-artificial cells as parallel_forest->local_num_quadrants)
3129 {
3130 const unsigned int total_local_cells = this->n_active_cells();
3131
3132
3133 if (Utilities::MPI::n_mpi_processes(this->mpi_communicator) == 1)
3134 {
3135 Assert(static_cast<unsigned int>(
3136 parallel_forest->local_num_quadrants) == total_local_cells,
3138 }
3139 else
3140 {
3141 Assert(static_cast<unsigned int>(
3142 parallel_forest->local_num_quadrants) <= total_local_cells,
3144 }
3145
3146 // count the number of owned, active cells and compare with p4est.
3147 unsigned int n_owned = 0;
3148 for (const auto &cell : this->active_cell_iterators())
3149 {
3150 if (cell->subdomain_id() == this->my_subdomain)
3151 ++n_owned;
3152 }
3153
3154 Assert(static_cast<unsigned int>(
3155 parallel_forest->local_num_quadrants) == n_owned,
3157 }
3158# endif
3159
3160 this->smooth_grid = save_smooth;
3161
3162 // finally, after syncing the parallel_forest with the triangulation,
3163 // also update the cell_relations, which will be used for
3164 // repartitioning, further refinement/coarsening, and unpacking
3165 // of stored or transferred data.
3166 update_cell_relations();
3167 }
3168
3169
3170
3171 template <int dim, int spacedim>
3175 {
3176 // Call the other function
3177 std::vector<Point<dim>> point{p};
3178 std::vector<types::subdomain_id> owner = find_point_owner_rank(point);
3179
3180 return owner[0];
3181 }
3182
3183
3184
3185 template <int dim, int spacedim>
3187 std::vector<types::subdomain_id> Triangulation<dim, spacedim>::
3188 find_point_owner_rank(const std::vector<Point<dim>> &points)
3189 {
3190 // We can only use this function if vertices are communicated to p4est
3191 AssertThrow(this->are_vertices_communicated_to_p4est(),
3192 ExcMessage(
3193 "Vertices need to be communicated to p4est to use this "
3194 "function. This must explicitly be turned on in the "
3195 "settings of the triangulation's constructor."));
3196
3197 // We can only use this function if all manifolds are flat
3198 for (const auto &manifold_id : this->get_manifold_ids())
3199 {
3201 manifold_id == numbers::flat_manifold_id,
3202 ExcMessage(
3203 "This function can only be used if the triangulation "
3204 "has no other manifold than a Cartesian (flat) manifold attached."));
3205 }
3206
3207 // Create object for callback
3208 PartitionSearch<dim> partition_search;
3209
3210 // Pointer should be this triangulation before we set it to something else
3211 Assert(parallel_forest->user_pointer == this, ExcInternalError());
3212
3213 // re-assign p4est's user pointer
3214 parallel_forest->user_pointer = &partition_search;
3215
3216 //
3217 // Copy points into p4est internal array data struct
3218 //
3219 // pointer to an array of points.
3220 sc_array_t *point_sc_array;
3221 // allocate memory for a number of dim-dimensional points including their
3222 // MPI rank, i.e., dim + 1 fields
3223 point_sc_array =
3224 sc_array_new_count(sizeof(double[dim + 1]), points.size());
3225
3226 // now assign the actual value
3227 for (size_t i = 0; i < points.size(); ++i)
3228 {
3229 // alias
3230 const Point<dim> &p = points[i];
3231 // get a non-const view of the array
3232 double *this_sc_point =
3233 static_cast<double *>(sc_array_index_ssize_t(point_sc_array, i));
3234 // fill this with the point data
3235 for (unsigned int d = 0; d < dim; ++d)
3236 {
3237 this_sc_point[d] = p(d);
3238 }
3239 this_sc_point[dim] = -1.0; // owner rank
3240 }
3241
3243 parallel_forest,
3244 /* execute quadrant function when leaving quadrant */
3245 static_cast<int>(false),
3246 &PartitionSearch<dim>::local_quadrant_fn,
3247 &PartitionSearch<dim>::local_point_fn,
3248 point_sc_array);
3249
3250 // copy the points found to an std::array
3251 std::vector<types::subdomain_id> owner_rank(
3252 points.size(), numbers::invalid_subdomain_id);
3253
3254 // fill the array
3255 for (size_t i = 0; i < points.size(); ++i)
3256 {
3257 // get a non-const view of the array
3258 double *this_sc_point =
3259 static_cast<double *>(sc_array_index_ssize_t(point_sc_array, i));
3260 Assert(this_sc_point[dim] >= 0. || this_sc_point[dim] == -1.,
3262 if (this_sc_point[dim] < 0.)
3263 owner_rank[i] = numbers::invalid_subdomain_id;
3264 else
3265 owner_rank[i] =
3266 static_cast<types::subdomain_id>(this_sc_point[dim]);
3267 }
3268
3269 // reset the internal pointer to this triangulation
3270 parallel_forest->user_pointer = this;
3271
3272 // release the memory (otherwise p4est will complain)
3273 sc_array_destroy_null(&point_sc_array);
3274
3275 return owner_rank;
3276 }
3277
3278
3279
3280 template <int dim, int spacedim>
3282 void Triangulation<dim, spacedim>::execute_coarsening_and_refinement()
3283 {
3284 // do not allow anisotropic refinement
3285# ifdef DEBUG
3286 for (const auto &cell : this->active_cell_iterators())
3287 if (cell->is_locally_owned() && cell->refine_flag_set())
3288 Assert(cell->refine_flag_set() ==
3290 ExcMessage(
3291 "This class does not support anisotropic refinement"));
3292# endif
3293
3294
3295 // safety check: p4est has an upper limit on the level of a cell
3296 if (this->n_levels() ==
3298 {
3300 cell = this->begin_active(
3302 cell !=
3304 1);
3305 ++cell)
3306 {
3308 !(cell->refine_flag_set()),
3309 ExcMessage(
3310 "Fatal Error: maximum refinement level of p4est reached."));
3311 }
3312 }
3313
3314 this->prepare_coarsening_and_refinement();
3315
3316 // signal that refinement is going to happen
3317 this->signals.pre_distributed_refinement();
3318
3319 // now do the work we're supposed to do when we are in charge
3320 // make sure all flags are cleared on cells we don't own, since nothing
3321 // good can come of that if they are still around
3322 for (const auto &cell : this->active_cell_iterators())
3323 if (cell->is_ghost() || cell->is_artificial())
3324 {
3325 cell->clear_refine_flag();
3326 cell->clear_coarsen_flag();
3327 }
3328
3329
3330 // count how many cells will be refined and coarsened, and allocate that
3331 // much memory
3332 RefineAndCoarsenList<dim, spacedim> refine_and_coarsen_list(
3333 *this, p4est_tree_to_coarse_cell_permutation, this->my_subdomain);
3334
3335 // copy refine and coarsen flags into p4est and execute the refinement
3336 // and coarsening. this uses the refine_and_coarsen_list just built,
3337 // which is communicated to the callback functions through
3338 // p4est's user_pointer object
3339 Assert(parallel_forest->user_pointer == this, ExcInternalError());
3340 parallel_forest->user_pointer = &refine_and_coarsen_list;
3341
3342 if (parallel_ghost != nullptr)
3343 {
3345 parallel_ghost);
3346 parallel_ghost = nullptr;
3347 }
3349 parallel_forest,
3350 /* refine_recursive */ false,
3351 &RefineAndCoarsenList<dim, spacedim>::refine_callback,
3352 /*init_callback=*/nullptr);
3354 parallel_forest,
3355 /* coarsen_recursive */ false,
3356 &RefineAndCoarsenList<dim, spacedim>::coarsen_callback,
3357 /*init_callback=*/nullptr);
3358
3359 // make sure all cells in the lists have been consumed
3360 Assert(refine_and_coarsen_list.pointers_are_at_end(), ExcInternalError());
3361
3362 // reset the pointer
3363 parallel_forest->user_pointer = this;
3364
3365 // enforce 2:1 hanging node condition
3367 parallel_forest,
3368 /* face and corner balance */
3369 (dim == 2 ? typename ::internal::p4est::types<dim>::balance_type(
3370 P4EST_CONNECT_FULL) :
3371 typename ::internal::p4est::types<dim>::balance_type(
3372 P8EST_CONNECT_FULL)),
3373 /*init_callback=*/nullptr);
3374
3375 // since refinement and/or coarsening on the parallel forest
3376 // has happened, we need to update the quadrant cell relations
3377 update_cell_relations();
3378
3379 // signals that parallel_forest has been refined and cell relations have
3380 // been updated
3381 this->signals.post_p4est_refinement();
3382
3383 // before repartitioning the mesh, save a copy of the current positions
3384 // of quadrants only if data needs to be transferred later
3385 std::vector<typename ::internal::p4est::types<dim>::gloidx>
3386 previous_global_first_quadrant;
3387
3388 if (this->cell_attached_data.n_attached_data_sets > 0)
3389 {
3390 previous_global_first_quadrant.resize(parallel_forest->mpisize + 1);
3391 std::memcpy(previous_global_first_quadrant.data(),
3392 parallel_forest->global_first_quadrant,
3393 sizeof(
3394 typename ::internal::p4est::types<dim>::gloidx) *
3395 (parallel_forest->mpisize + 1));
3396 }
3397
3398 if (!(settings & no_automatic_repartitioning))
3399 {
3400 // partition the new mesh between all processors. If cell weights
3401 // have not been given balance the number of cells.
3402 if (this->signals.weight.empty())
3404 parallel_forest,
3405 /* prepare coarsening */ 1,
3406 /* weight_callback */ nullptr);
3407 else
3408 {
3409 // get cell weights for a weighted repartitioning.
3410 const std::vector<unsigned int> cell_weights = get_cell_weights();
3411
3412 // verify that the global sum of weights is larger than 0
3413 Assert(Utilities::MPI::sum(std::accumulate(cell_weights.begin(),
3414 cell_weights.end(),
3415 std::uint64_t(0)),
3416 this->mpi_communicator) > 0,
3417 ExcMessage(
3418 "The global sum of weights over all active cells "
3419 "is zero. Please verify how you generate weights."));
3420
3421 PartitionWeights<dim, spacedim> partition_weights(cell_weights);
3422
3423 // attach (temporarily) a pointer to the cell weights through
3424 // p4est's user_pointer object
3425 Assert(parallel_forest->user_pointer == this, ExcInternalError());
3426 parallel_forest->user_pointer = &partition_weights;
3427
3429 parallel_forest,
3430 /* prepare coarsening */ 1,
3431 /* weight_callback */
3432 &PartitionWeights<dim, spacedim>::cell_weight);
3433
3434 // release data
3436 parallel_forest, 0, nullptr, nullptr);
3437 // reset the user pointer to its previous state
3438 parallel_forest->user_pointer = this;
3439 }
3440 }
3441
3442 // pack data before triangulation gets updated
3443 if (this->cell_attached_data.n_attached_data_sets > 0)
3444 {
3445 this->data_serializer.pack_data(
3446 this->local_cell_relations,
3447 this->cell_attached_data.pack_callbacks_fixed,
3448 this->cell_attached_data.pack_callbacks_variable,
3449 this->get_mpi_communicator());
3450 }
3451
3452 // finally copy back from local part of tree to deal.II
3453 // triangulation. before doing so, make sure there are no refine or
3454 // coarsen flags pending
3455 for (const auto &cell : this->active_cell_iterators())
3456 {
3457 cell->clear_refine_flag();
3458 cell->clear_coarsen_flag();
3459 }
3460
3461 try
3462 {
3463 copy_local_forest_to_triangulation();
3464 }
3465 catch (const typename Triangulation<dim>::DistortedCellList &)
3466 {
3467 // the underlying triangulation should not be checking for distorted
3468 // cells
3470 }
3471
3472 // transfer data after triangulation got updated
3473 if (this->cell_attached_data.n_attached_data_sets > 0)
3474 {
3475 this->execute_transfer(parallel_forest,
3476 previous_global_first_quadrant.data());
3477
3478 // also update the CellStatus information on the new mesh
3479 this->data_serializer.unpack_cell_status(this->local_cell_relations);
3480 }
3481
3482# ifdef DEBUG
3483 // Check that we know the level subdomain ids of all our neighbors. This
3484 // also involves coarser cells that share a vertex if they are active.
3485 //
3486 // Example (M= my, O=other):
3487 // *------*
3488 // | |
3489 // | O |
3490 // | |
3491 // *---*---*------*
3492 // | M | M |
3493 // *---*---*
3494 // | | M |
3495 // *---*---*
3496 // ^- the parent can be owned by somebody else, so O is not a neighbor
3497 // one level coarser
3498 if (settings & construct_multigrid_hierarchy)
3499 {
3500 for (unsigned int lvl = 0; lvl < this->n_global_levels(); ++lvl)
3501 {
3502 std::vector<bool> active_verts =
3503 this->mark_locally_active_vertices_on_level(lvl);
3504
3505 const unsigned int maybe_coarser_lvl =
3506 (lvl > 0) ? (lvl - 1) : lvl;
3508 cell = this->begin(maybe_coarser_lvl),
3509 endc = this->end(lvl);
3510 for (; cell != endc; ++cell)
3511 if (cell->level() == static_cast<int>(lvl) || cell->is_active())
3512 {
3513 const bool is_level_artificial =
3514 (cell->level_subdomain_id() ==
3516 bool need_to_know = false;
3517 for (const unsigned int vertex :
3519 if (active_verts[cell->vertex_index(vertex)])
3520 {
3521 need_to_know = true;
3522 break;
3523 }
3524
3525 Assert(
3526 !need_to_know || !is_level_artificial,
3527 ExcMessage(
3528 "Internal error: the owner of cell" +
3529 cell->id().to_string() +
3530 " is unknown even though it is needed for geometric multigrid."));
3531 }
3532 }
3533 }
3534# endif
3535
3536 this->update_periodic_face_map();
3537 this->update_number_cache();
3538
3539 // signal that refinement is finished
3540 this->signals.post_distributed_refinement();
3541 }
3542
3543
3544
3545 template <int dim, int spacedim>
3547 void Triangulation<dim, spacedim>::repartition()
3548 {
3549# ifdef DEBUG
3550 for (const auto &cell : this->active_cell_iterators())
3551 if (cell->is_locally_owned())
3552 Assert(
3553 !cell->refine_flag_set() && !cell->coarsen_flag_set(),
3554 ExcMessage(
3555 "Error: There shouldn't be any cells flagged for coarsening/refinement when calling repartition()."));
3556# endif
3557
3558 // signal that repartitioning is going to happen
3559 this->signals.pre_distributed_repartition();
3560
3561 // before repartitioning the mesh, save a copy of the current positions
3562 // of quadrants only if data needs to be transferred later
3563 std::vector<typename ::internal::p4est::types<dim>::gloidx>
3564 previous_global_first_quadrant;
3565
3566 if (this->cell_attached_data.n_attached_data_sets > 0)
3567 {
3568 previous_global_first_quadrant.resize(parallel_forest->mpisize + 1);
3569 std::memcpy(previous_global_first_quadrant.data(),
3570 parallel_forest->global_first_quadrant,
3571 sizeof(
3572 typename ::internal::p4est::types<dim>::gloidx) *
3573 (parallel_forest->mpisize + 1));
3574 }
3575
3576 if (this->signals.weight.empty())
3577 {
3578 // no cell weights given -- call p4est's 'partition' without a
3579 // callback for cell weights
3581 parallel_forest,
3582 /* prepare coarsening */ 1,
3583 /* weight_callback */ nullptr);
3584 }
3585 else
3586 {
3587 // get cell weights for a weighted repartitioning.
3588 const std::vector<unsigned int> cell_weights = get_cell_weights();
3589
3590 // verify that the global sum of weights is larger than 0
3591 Assert(Utilities::MPI::sum(std::accumulate(cell_weights.begin(),
3592 cell_weights.end(),
3593 std::uint64_t(0)),
3594 this->mpi_communicator) > 0,
3595 ExcMessage(
3596 "The global sum of weights over all active cells "
3597 "is zero. Please verify how you generate weights."));
3598
3599 PartitionWeights<dim, spacedim> partition_weights(cell_weights);
3600
3601 // attach (temporarily) a pointer to the cell weights through
3602 // p4est's user_pointer object
3603 Assert(parallel_forest->user_pointer == this, ExcInternalError());
3604 parallel_forest->user_pointer = &partition_weights;
3605
3607 parallel_forest,
3608 /* prepare coarsening */ 1,
3609 /* weight_callback */
3610 &PartitionWeights<dim, spacedim>::cell_weight);
3611
3612 // reset the user pointer to its previous state
3613 parallel_forest->user_pointer = this;
3614 }
3615
3616 // pack data before triangulation gets updated
3617 if (this->cell_attached_data.n_attached_data_sets > 0)
3618 {
3619 this->data_serializer.pack_data(
3620 this->local_cell_relations,
3621 this->cell_attached_data.pack_callbacks_fixed,
3622 this->cell_attached_data.pack_callbacks_variable,
3623 this->get_mpi_communicator());
3624 }
3625
3626 try
3627 {
3628 copy_local_forest_to_triangulation();
3629 }
3630 catch (const typename Triangulation<dim>::DistortedCellList &)
3631 {
3632 // the underlying triangulation should not be checking for distorted
3633 // cells
3635 }
3636
3637 // transfer data after triangulation got updated
3638 if (this->cell_attached_data.n_attached_data_sets > 0)
3639 {
3640 this->execute_transfer(parallel_forest,
3641 previous_global_first_quadrant.data());
3642 }
3643
3644 this->update_periodic_face_map();
3645
3646 // update how many cells, edges, etc, we store locally
3647 this->update_number_cache();
3648
3649 // signal that repartitioning is finished
3650 this->signals.post_distributed_repartition();
3651 }
3652
3653
3654
3655 template <int dim, int spacedim>
3657 const std::vector<types::global_dof_index>
3659 const
3660 {
3661 return p4est_tree_to_coarse_cell_permutation;
3662 }
3663
3664
3665
3666 template <int dim, int spacedim>
3668 const std::vector<types::global_dof_index>
3670 const
3671 {
3672 return coarse_cell_to_p4est_tree_permutation;
3673 }
3674
3675
3676
3677 template <int dim, int spacedim>
3679 std::vector<bool> Triangulation<dim, spacedim>::
3680 mark_locally_active_vertices_on_level(const int level) const
3681 {
3682 Assert(dim > 1, ExcNotImplemented());
3683
3684 std::vector<bool> marked_vertices(this->n_vertices(), false);
3685 for (const auto &cell : this->cell_iterators_on_level(level))
3686 if (cell->level_subdomain_id() == this->locally_owned_subdomain())
3687 for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
3688 marked_vertices[cell->vertex_index(v)] = true;
3689
3695 // When a connectivity in the code below is detected, the assignment
3696 // 'marked_vertices[v1] = marked_vertices[v2] = true' makes sure that
3697 // the information about the periodicity propagates back to vertices on
3698 // cells that are not owned locally. However, in the worst case we want
3699 // to connect to a vertex that is 'dim' hops away from the locally owned
3700 // cell. Depending on the order of the periodic face map, we might
3701 // connect to that point by chance or miss it. However, after looping
3702 // through all the periodic directions (which are at most as many as
3703 // the number of space dimensions) we can be sure that all connections
3704 // to vertices have been created.
3705 for (unsigned int repetition = 0; repetition < dim; ++repetition)
3706 for (const auto &it : this->get_periodic_face_map())
3707 {
3708 const cell_iterator &cell_1 = it.first.first;
3709 const unsigned int face_no_1 = it.first.second;
3710 const cell_iterator &cell_2 = it.second.first.first;
3711 const unsigned int face_no_2 = it.second.first.second;
3712 const unsigned char combined_orientation = it.second.second;
3713 const auto [orientation, rotation, flip] =
3714 ::internal::split_face_orientation(combined_orientation);
3715
3716 if (cell_1->level() == level && cell_2->level() == level)
3717 {
3718 for (unsigned int v = 0;
3719 v < GeometryInfo<dim - 1>::vertices_per_cell;
3720 ++v)
3721 {
3722 // take possible non-standard orientation of faces into
3723 // account
3724 const unsigned int vface0 =
3726 v, orientation, flip, rotation);
3727 if (marked_vertices[cell_1->face(face_no_1)->vertex_index(
3728 vface0)] ||
3729 marked_vertices[cell_2->face(face_no_2)->vertex_index(
3730 v)])
3731 marked_vertices[cell_1->face(face_no_1)->vertex_index(
3732 vface0)] =
3733 marked_vertices[cell_2->face(face_no_2)->vertex_index(
3734 v)] = true;
3735 }
3736 }
3737 }
3738
3739 return marked_vertices;
3740 }
3741
3742
3743
3744 template <int dim, int spacedim>
3746 unsigned int Triangulation<dim, spacedim>::
3747 coarse_cell_id_to_coarse_cell_index(
3748 const types::coarse_cell_id coarse_cell_id) const
3749 {
3750 return p4est_tree_to_coarse_cell_permutation[coarse_cell_id];
3751 }
3752
3753
3754
3755 template <int dim, int spacedim>
3759 const unsigned int coarse_cell_index) const
3760 {
3761 return coarse_cell_to_p4est_tree_permutation[coarse_cell_index];
3762 }
3763
3764
3765
3766 template <int dim, int spacedim>
3768 void Triangulation<dim, spacedim>::add_periodicity(
3769 const std::vector<::GridTools::PeriodicFacePair<cell_iterator>>
3770 &periodicity_vector)
3771 {
3772 Assert(triangulation_has_content == true,
3773 ExcMessage("The triangulation is empty!"));
3774 Assert(this->n_levels() == 1,
3775 ExcMessage("The triangulation is refined!"));
3776
3777 // call the base class for storing the periodicity information; we must
3778 // do this before going to p4est and rebuilding the triangulation to get
3779 // the level subdomain ids correct in the multigrid case
3781
3782 const auto reference_cell = ReferenceCells::get_hypercube<dim>();
3783 const auto face_reference_cell = ReferenceCells::get_hypercube<dim - 1>();
3784 for (const auto &face_pair : periodicity_vector)
3785 {
3786 const cell_iterator first_cell = face_pair.cell[0];
3787 const cell_iterator second_cell = face_pair.cell[1];
3788 const unsigned int face_left = face_pair.face_idx[0];
3789 const unsigned int face_right = face_pair.face_idx[1];
3790
3791 // respective cells of the matching faces in p4est
3792 const unsigned int tree_left =
3793 coarse_cell_to_p4est_tree_permutation[first_cell->index()];
3794 const unsigned int tree_right =
3795 coarse_cell_to_p4est_tree_permutation[second_cell->index()];
3796
3797 // p4est wants to know which corner the first corner on the face with
3798 // the lower id is mapped to on the face with with the higher id. For
3799 // d==2 there are only two possibilities: i.e., face_pair.orientation
3800 // must be 0 or 1. For d==3 we have to use a lookup table. The result
3801 // is given below.
3802
3803 unsigned int p4est_orientation = 0;
3804 if (dim == 2)
3805 {
3806 AssertIndexRange(face_pair.orientation, 2);
3807 p4est_orientation =
3808 face_pair.orientation ==
3810 0u :
3811 1u;
3812 }
3813 else
3814 {
3815 const unsigned int face_idx_list[] = {face_left, face_right};
3816 const cell_iterator cell_list[] = {first_cell, second_cell};
3817 unsigned int lower_idx, higher_idx;
3818 unsigned char orientation;
3819 if (face_left <= face_right)
3820 {
3821 higher_idx = 1;
3822 lower_idx = 0;
3823 orientation =
3824 face_reference_cell.get_inverse_combined_orientation(
3825 face_pair.orientation);
3826 }
3827 else
3828 {
3829 higher_idx = 0;
3830 lower_idx = 1;
3831 orientation = face_pair.orientation;
3832 }
3833
3834 // get the cell index of the first index on the face with the
3835 // lower id
3836 unsigned int first_p4est_idx_on_cell =
3837 p8est_face_corners[face_idx_list[lower_idx]][0];
3838 unsigned int first_dealii_idx_on_face =
3840 for (unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_face;
3841 ++i)
3842 {
3843 const unsigned int first_dealii_idx_on_cell =
3845 face_idx_list[lower_idx],
3846 i,
3847 cell_list[lower_idx]->face_orientation(
3848 face_idx_list[lower_idx]),
3849 cell_list[lower_idx]->face_flip(face_idx_list[lower_idx]),
3850 cell_list[lower_idx]->face_rotation(
3851 face_idx_list[lower_idx]));
3852 if (first_p4est_idx_on_cell == first_dealii_idx_on_cell)
3853 {
3854 first_dealii_idx_on_face = i;
3855 break;
3856 }
3857 }
3858 Assert(first_dealii_idx_on_face != numbers::invalid_unsigned_int,
3860
3861 // Now map dealii_idx_on_face according to the orientation.
3862 const unsigned int second_dealii_idx_on_face =
3863 reference_cell.standard_to_real_face_vertex(
3864 first_dealii_idx_on_face,
3865 face_idx_list[lower_idx],
3866 orientation);
3867 const unsigned int second_dealii_idx_on_cell =
3868 reference_cell.face_to_cell_vertices(
3869 face_idx_list[higher_idx],
3870 second_dealii_idx_on_face,
3871 cell_list[higher_idx]->combined_face_orientation(
3872 face_idx_list[higher_idx]));
3873 // map back to p4est
3874 const unsigned int second_p4est_idx_on_face =
3875 p8est_corner_face_corners[second_dealii_idx_on_cell]
3876 [face_idx_list[higher_idx]];
3877 p4est_orientation = second_p4est_idx_on_face;
3878 }
3879
3881 connectivity,
3882 tree_left,
3883 tree_right,
3884 face_left,
3885 face_right,
3886 p4est_orientation);
3887 }
3888
3889
3891 connectivity) == 1,
3893
3894 // now create a forest out of the connectivity data structure
3897 this->mpi_communicator,
3898 connectivity,
3899 /* minimum initial number of quadrants per tree */ 0,
3900 /* minimum level of upfront refinement */ 0,
3901 /* use uniform upfront refinement */ 1,
3902 /* user_data_size = */ 0,
3903 /* user_data_constructor = */ nullptr,
3904 /* user_pointer */ this);
3905
3906 try
3907 {
3908 copy_local_forest_to_triangulation();
3909 }
3910 catch (const typename Triangulation<dim>::DistortedCellList &)
3911 {
3912 // the underlying triangulation should not be checking for distorted
3913 // cells
3915 }
3916
3917 // The range of ghost_owners might have changed so update that
3918 // information
3919 this->update_number_cache();
3920 }
3921
3922
3923
3924 template <int dim, int spacedim>
3926 std::size_t Triangulation<dim, spacedim>::memory_consumption() const
3927 {
3928 std::size_t mem =
3931 MemoryConsumption::memory_consumption(triangulation_has_content) +
3933 MemoryConsumption::memory_consumption(parallel_forest) +
3935 this->cell_attached_data.n_attached_data_sets) +
3936 // MemoryConsumption::memory_consumption(cell_attached_data.pack_callbacks_fixed)
3937 // +
3938 // MemoryConsumption::memory_consumption(cell_attached_data.pack_callbacks_variable)
3939 // +
3940 // TODO[TH]: how?
3942 coarse_cell_to_p4est_tree_permutation) +
3944 p4est_tree_to_coarse_cell_permutation) +
3945 memory_consumption_p4est();
3946
3947 return mem;
3948 }
3949
3950
3951
3952 template <int dim, int spacedim>
3954 std::size_t Triangulation<dim, spacedim>::memory_consumption_p4est() const
3955 {
3956 return ::internal::p4est::functions<dim>::forest_memory_used(
3957 parallel_forest) +
3959 connectivity);
3960 }
3961
3962
3963
3964 template <int dim, int spacedim>
3966 void Triangulation<dim, spacedim>::copy_triangulation(
3967 const ::Triangulation<dim, spacedim> &other_tria)
3968 {
3969 Assert(
3970 (dynamic_cast<
3971 const ::parallel::distributed::Triangulation<dim, spacedim> *>(
3972 &other_tria)) ||
3973 (other_tria.n_global_levels() == 1),
3975
3977
3978 try
3979 {
3981 copy_triangulation(other_tria);
3982 }
3983 catch (
3984 const typename ::Triangulation<dim, spacedim>::DistortedCellList
3985 &)
3986 {
3987 // the underlying triangulation should not be checking for distorted
3988 // cells
3990 }
3991
3992 if (const ::parallel::distributed::Triangulation<dim, spacedim>
3993 *other_distributed =
3994 dynamic_cast<const ::parallel::distributed::
3995 Triangulation<dim, spacedim> *>(&other_tria))
3996 {
3997 // copy parallel distributed specifics
3998 settings = other_distributed->settings;
3999 triangulation_has_content =
4000 other_distributed->triangulation_has_content;
4001 coarse_cell_to_p4est_tree_permutation =
4002 other_distributed->coarse_cell_to_p4est_tree_permutation;
4003 p4est_tree_to_coarse_cell_permutation =
4004 other_distributed->p4est_tree_to_coarse_cell_permutation;
4005
4006 // create deep copy of connectivity graph
4007 typename ::internal::p4est::types<dim>::connectivity
4008 *temp_connectivity = const_cast<
4009 typename ::internal::p4est::types<dim>::connectivity *>(
4010 other_distributed->connectivity);
4011 connectivity =
4012 ::internal::p4est::copy_connectivity<dim>(temp_connectivity);
4013
4014 // create deep copy of parallel forest
4015 typename ::internal::p4est::types<dim>::forest *temp_forest =
4016 const_cast<typename ::internal::p4est::types<dim>::forest *>(
4017 other_distributed->parallel_forest);
4018 parallel_forest =
4020 false);
4021 parallel_forest->connectivity = connectivity;
4022 parallel_forest->user_pointer = this;
4023 }
4024 else
4025 {
4026 triangulation_has_content = true;
4027 setup_coarse_cell_to_p4est_tree_permutation();
4028 copy_new_triangulation_to_p4est(std::integral_constant<int, dim>());
4029 }
4030
4031 try
4032 {
4033 copy_local_forest_to_triangulation();
4034 }
4035 catch (const typename Triangulation<dim>::DistortedCellList &)
4036 {
4037 // the underlying triangulation should not be checking for distorted
4038 // cells
4040 }
4041
4042 this->update_periodic_face_map();
4043 this->update_number_cache();
4044 }
4045
4046
4047
4048 template <int dim, int spacedim>
4050 void Triangulation<dim, spacedim>::update_cell_relations()
4051 {
4052 // reorganize memory for local_cell_relations
4053 this->local_cell_relations.resize(parallel_forest->local_num_quadrants);
4054 this->local_cell_relations.shrink_to_fit();
4055
4056 // recurse over p4est
4057 for (const auto &cell : this->cell_iterators_on_level(0))
4058 {
4059 // skip coarse cells that are not ours
4060 if (tree_exists_locally<dim, spacedim>(
4061 parallel_forest,
4062 coarse_cell_to_p4est_tree_permutation[cell->index()]) == false)
4063 continue;
4064
4065 // initialize auxiliary top level p4est quadrant
4066 typename ::internal::p4est::types<dim>::quadrant
4067 p4est_coarse_cell;
4068 ::internal::p4est::init_coarse_quadrant<dim>(p4est_coarse_cell);
4069
4070 // determine tree to start recursion on
4071 typename ::internal::p4est::types<dim>::tree *tree =
4072 init_tree(cell->index());
4073
4074 update_cell_relations_recursively<dim, spacedim>(
4075 this->local_cell_relations, *tree, cell, p4est_coarse_cell);
4076 }
4077 }
4078
4079
4080
4081 template <int dim, int spacedim>
4083 std::vector<unsigned int> Triangulation<dim, spacedim>::get_cell_weights()
4084 const
4085 {
4086 // check if local_cell_relations have been previously gathered
4087 // correctly
4088 Assert(this->local_cell_relations.size() ==
4089 static_cast<unsigned int>(parallel_forest->local_num_quadrants),
4091
4092 // Allocate the space for the weights. We reserve an integer for each
4093 // locally owned quadrant on the already refined p4est object.
4094 std::vector<unsigned int> weights;
4095 weights.reserve(this->local_cell_relations.size());
4096
4097 // Iterate over p4est and Triangulation relations
4098 // to find refined/coarsened/kept
4099 // cells. Then append weight.
4100 // Note that we need to follow the p4est ordering
4101 // instead of the deal.II ordering to get the weights
4102 // in the same order p4est will encounter them during repartitioning.
4103 for (const auto &cell_rel : this->local_cell_relations)
4104 {
4105 const auto &cell_it = cell_rel.first;
4106 const auto &cell_status = cell_rel.second;
4107
4108 weights.push_back(this->signals.weight(cell_it, cell_status));
4109 }
4110
4111 return weights;
4112 }
4113
4114
4115
4116 template <int spacedim>
4119 const MPI_Comm mpi_communicator,
4120 const typename ::Triangulation<1, spacedim>::MeshSmoothing
4121 smooth_grid,
4122 const Settings /*settings*/)
4123 : ::parallel::DistributedTriangulationBase<1, spacedim>(
4124 mpi_communicator,
4125 smooth_grid,
4126 false)
4127 {
4129 }
4130
4131
4132 template <int spacedim>
4135 {
4137 }
4138
4139
4140
4141 template <int spacedim>
4143 const std::vector<types::global_dof_index>
4145 const
4146 {
4147 static std::vector<types::global_dof_index> a;
4148 return a;
4149 }
4150
4151
4152
4153 template <int spacedim>
4155 std::map<unsigned int,
4156 std::set<::types::subdomain_id>> Triangulation<1, spacedim>::
4158 const unsigned int /*level*/) const
4159 {
4161
4162 return std::map<unsigned int, std::set<::types::subdomain_id>>();
4163 }
4164
4165
4166
4167 template <int spacedim>
4169 std::vector<bool> Triangulation<1, spacedim>::
4170 mark_locally_active_vertices_on_level(const unsigned int) const
4171 {
4173 return std::vector<bool>();
4174 }
4175
4176
4177
4178 template <int spacedim>
4180 unsigned int Triangulation<1, spacedim>::
4181 coarse_cell_id_to_coarse_cell_index(const types::coarse_cell_id) const
4182 {
4184 return 0;
4185 }
4186
4187
4188
4189 template <int spacedim>
4193 const unsigned int) const
4194 {
4196 return 0;
4197 }
4198
4199
4200
4201 template <int spacedim>
4203 void Triangulation<1, spacedim>::load(const std::string &)
4204 {
4206 }
4207
4208
4209
4210 template <int spacedim>
4212 void Triangulation<1, spacedim>::save(const std::string &) const
4213 {
4215 }
4216
4217
4218
4219 template <int spacedim>
4221 bool Triangulation<1, spacedim>::is_multilevel_hierarchy_constructed() const
4222 {
4224 return false;
4225 }
4226
4227
4228
4229 template <int spacedim>
4231 bool Triangulation<1, spacedim>::are_vertices_communicated_to_p4est() const
4232 {
4234 return false;
4235 }
4236
4237
4238
4239 template <int spacedim>
4241 void Triangulation<1, spacedim>::update_cell_relations()
4242 {
4244 }
4245
4246 } // namespace distributed
4247} // namespace parallel
4248
4249
4250#endif // DEAL_II_WITH_P4EST
4251
4252
4253
4254namespace parallel
4255{
4256 namespace distributed
4257 {
4258 template <int dim, int spacedim>
4261 : distributed_tria(
4262 dynamic_cast<
4263 ::parallel::distributed::Triangulation<dim, spacedim> *>(
4264 &tria))
4265 {
4266#ifdef DEAL_II_WITH_P4EST
4267 if (distributed_tria != nullptr)
4268 {
4269 // Save the current set of refinement flags, and adjust the
4270 // refinement flags to be consistent with the p4est oracle.
4271 distributed_tria->save_coarsen_flags(saved_coarsen_flags);
4272 distributed_tria->save_refine_flags(saved_refine_flags);
4273
4274 for (const auto &pair : distributed_tria->local_cell_relations)
4275 {
4276 const auto &cell = pair.first;
4277 const auto &status = pair.second;
4278
4279 switch (status)
4280 {
4282 // cell remains unchanged
4283 cell->clear_refine_flag();
4284 cell->clear_coarsen_flag();
4285 break;
4286
4288 // cell will be refined
4289 cell->clear_coarsen_flag();
4290 cell->set_refine_flag();
4291 break;
4292
4294 // children of this cell will be coarsened
4295 for (const auto &child : cell->child_iterators())
4296 {
4297 child->clear_refine_flag();
4298 child->set_coarsen_flag();
4299 }
4300 break;
4301
4303 // do nothing as cell does not exist yet
4304 break;
4305
4306 default:
4308 break;
4309 }
4310 }
4311 }
4312#endif
4313 }
4314
4315
4316
4317 template <int dim, int spacedim>
4319 {
4320#ifdef DEAL_II_WITH_P4EST
4321 if (distributed_tria)
4322 {
4323 // Undo the refinement flags modification.
4324 distributed_tria->load_coarsen_flags(saved_coarsen_flags);
4325 distributed_tria->load_refine_flags(saved_refine_flags);
4326 }
4327#else
4328 // pretend that this destructor does something to silence clang-tidy
4329 (void)distributed_tria;
4330#endif
4331 }
4332 } // namespace distributed
4333} // namespace parallel
4334
4335
4336
4337/*-------------- Explicit Instantiations -------------------------------*/
4338#include "tria.inst"
4339
4340
CellStatus
Definition cell_status.h:31
@ cell_will_be_refined
@ children_will_be_coarsened
Definition point.h:111
static constexpr unsigned char default_combined_face_orientation()
virtual void add_periodicity(const std::vector< GridTools::PeriodicFacePair< cell_iterator > > &)
virtual void clear()
const std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, unsigned char > > & get_periodic_face_map() const
active_cell_iterator last_active() const
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
cell_iterator end() const
virtual types::coarse_cell_id coarse_cell_index_to_coarse_cell_id(const unsigned int coarse_cell_index) const
virtual void execute_coarsening_and_refinement()
virtual bool prepare_coarsening_and_refinement()
void save_refine_flags(std::ostream &out) const
unsigned int n_vertices() const
void save_coarsen_flags(std::ostream &out) const
active_cell_iterator begin_active(const unsigned int level=0) const
virtual std::size_t memory_consumption() const override
Definition tria_base.cc:92
virtual void clear() override
Definition tria_base.cc:687
virtual void copy_triangulation(const ::Triangulation< dim, spacedim > &old_tria) override
Definition tria_base.cc:67
const ObserverPointer< ::parallel::distributed::Triangulation< dim, spacedim > > distributed_tria
Definition tria.h:1180
virtual void clear() override
Definition tria.cc:1861
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:175
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
unsigned int level
Definition grid_out.cc:4626
unsigned int vertex_indices[2]
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertNothrow(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
typename ::Triangulation< dim, spacedim >::cell_iterator cell_iterator
Definition tria.h:324
typename ::Triangulation< dim, spacedim >::active_cell_iterator active_cell_iterator
Definition tria.h:345
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
std::size_t size
Definition mpi.cc:734
void get_vertex_connectivity_of_cells(const Triangulation< dim, spacedim > &triangulation, DynamicSparsityPattern &connectivity)
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:191
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
constexpr const ReferenceCell & get_hypercube()
void reorder_hierarchical(const DynamicSparsityPattern &sparsity, std::vector< DynamicSparsityPattern::size_type > &new_indices)
VectorType::value_type * end(VectorType &V)
T sum(const T &t, const MPI_Comm mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
Definition mpi.cc:92
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
Definition mpi.cc:107
T broadcast(const MPI_Comm comm, const T &object_to_send, const unsigned int root_process=0)
std::vector< Integer > invert_permutation(const std::vector< Integer > &permutation)
Definition utilities.h:1699
bool tree_exists_locally(const typename types< dim >::forest *parallel_forest, const typename types< dim >::topidx coarse_grid_cell)
void exchange_refinement_flags(::parallel::distributed::Triangulation< dim, spacedim > &tria)
Definition tria.cc:56
std::tuple< bool, bool, bool > split_face_orientation(const unsigned char combined_face_orientation)
const types::subdomain_id artificial_subdomain_id
Definition types.h:362
const types::subdomain_id invalid_subdomain_id
Definition types.h:341
static const unsigned int invalid_unsigned_int
Definition types.h:220
const types::manifold_id flat_manifold_id
Definition types.h:325
STL namespace.
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
Definition types.h:32
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
static unsigned int standard_to_real_face_vertex(const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices()
static bool is_inside_unit_cell(const Point< dim > &p)
static Point< dim > unit_cell_vertex(const unsigned int vertex)