deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
p4est_wrappers.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2016 - 2025 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
19#ifdef DEAL_II_WITH_P4EST
20# include <p4est.h>
21# include <p8est.h>
22# include <sc_containers.h>
23
24// Below, we will use the P4EST_QUADRANT_INIT and P8EST_QUADRANT_INIT
25// function-like macros. If we are building the library based on
26// header files, we get these from the <p4est.h> and <p8est.h> header
27// inclusions. But if we build a C++20 module, we only import
28// declarations, not preprocessor macros. As a consequence, let us
29// duplicate these macros here, hoping that at some point, the p4est
30// library folks add regular functions that can do the job.
31# ifndef P4EST_QUADRANT_INIT
32# define P4EST_QUADRANT_INIT(q) \
33 ((void)std::memset((q), -1, sizeof(p4est_quadrant_t)))
34# endif
35
36# ifndef P8EST_QUADRANT_INIT
37# define P8EST_QUADRANT_INIT(q) \
38 ((void)std::memset((q), -1, sizeof(p8est_quadrant_t)))
39# endif
40
41#endif
42
43
45
46#ifdef DEAL_II_WITH_P4EST
47
48namespace internal
49{
50 namespace p4est
51 {
52 namespace
53 {
54 template <int dim, int spacedim>
55 typename ::Triangulation<dim, spacedim>::cell_iterator
56 cell_from_quad(
57 const ::parallel::distributed::Triangulation<dim, spacedim>
58 *triangulation,
59 const typename ::internal::p4est::types<dim>::topidx treeidx,
60 const typename ::internal::p4est::types<dim>::quadrant &quad)
61 {
62 int i, l = quad.level;
63 ::types::global_dof_index dealii_index =
64 triangulation->get_p4est_tree_to_coarse_cell_permutation()[treeidx];
65
66 for (i = 0; i < l; ++i)
67 {
68 typename ::Triangulation<dim, spacedim>::cell_iterator cell(
69 triangulation, i, dealii_index);
70 const int child_id =
72 &quad, i + 1);
73 Assert(cell->has_children(),
74 ExcMessage("p4est quadrant does not correspond to a cell!"));
75 dealii_index = cell->child_index(child_id);
76 }
77
78 typename ::Triangulation<dim, spacedim>::cell_iterator out_cell(
79 triangulation, l, dealii_index);
80
81 return out_cell;
82 }
83
88 template <int dim, int spacedim>
89 struct FindGhosts
90 {
91 const typename ::parallel::distributed::Triangulation<dim,
92 spacedim>
93 *triangulation;
94 sc_array_t *subids;
95 std::map<unsigned int, std::set<::types::subdomain_id>>
96 *vertices_with_ghost_neighbors;
97 };
98
99
105 template <int dim, int spacedim>
106 void
107 find_ghosts_corner(
108 typename ::internal::p4est::iter<dim>::corner_info *info,
109 void *user_data)
110 {
111 int i, j;
112 int nsides = info->sides.elem_count;
113 auto *sides = reinterpret_cast<
114 typename ::internal::p4est::iter<dim>::corner_side *>(
115 info->sides.array);
116 FindGhosts<dim, spacedim> *fg =
117 static_cast<FindGhosts<dim, spacedim> *>(user_data);
118 sc_array_t *subids = fg->subids;
119 const ::parallel::distributed::Triangulation<dim, spacedim>
120 *triangulation = fg->triangulation;
121 int nsubs;
122 ::types::subdomain_id *subdomain_ids;
123 std::map<unsigned int, std::set<::types::subdomain_id>>
124 *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
125
126 subids->elem_count = 0;
127 for (i = 0; i < nsides; ++i)
128 {
129 if (sides[i].is_ghost)
130 {
131 typename ::parallel::distributed::
132 Triangulation<dim, spacedim>::cell_iterator cell =
133 cell_from_quad(triangulation,
134 sides[i].treeid,
135 *(sides[i].quad));
136 Assert(cell->is_ghost(),
137 ExcMessage("ghost quad did not find ghost cell"));
138 ::types::subdomain_id *subid =
139 static_cast<::types::subdomain_id *>(
140 sc_array_push(subids));
141 *subid = cell->subdomain_id();
142 }
143 }
144
145 if (!subids->elem_count)
146 {
147 return;
148 }
149
150 nsubs = static_cast<int>(subids->elem_count);
151 subdomain_ids =
152 reinterpret_cast<::types::subdomain_id *>(subids->array);
153
154 for (i = 0; i < nsides; ++i)
155 {
156 if (!sides[i].is_ghost)
157 {
158 typename ::parallel::distributed::
159 Triangulation<dim, spacedim>::cell_iterator cell =
160 cell_from_quad(triangulation,
161 sides[i].treeid,
162 *(sides[i].quad));
163
164 Assert(!cell->is_ghost(),
165 ExcMessage("local quad found ghost cell"));
166
167 for (j = 0; j < nsubs; ++j)
168 {
169 (*vertices_with_ghost_neighbors)[cell->vertex_index(
170 sides[i].corner)]
171 .insert(subdomain_ids[j]);
172 }
173 }
174 }
175
176 subids->elem_count = 0;
177 }
178
182 template <int dim, int spacedim>
183 void
184 find_ghosts_edge(
185 typename ::internal::p4est::iter<dim>::edge_info *info,
186 void *user_data)
187 {
188 int i, j, k;
189 int nsides = info->sides.elem_count;
190 auto *sides = reinterpret_cast<
191 typename ::internal::p4est::iter<dim>::edge_side *>(
192 info->sides.array);
193 auto *fg = static_cast<FindGhosts<dim, spacedim> *>(user_data);
194 sc_array_t *subids = fg->subids;
195 const ::parallel::distributed::Triangulation<dim, spacedim>
196 *triangulation = fg->triangulation;
197 int nsubs;
198 ::types::subdomain_id *subdomain_ids;
199 std::map<unsigned int, std::set<::types::subdomain_id>>
200 *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
201
202 subids->elem_count = 0;
203 for (i = 0; i < nsides; ++i)
204 {
205 if (sides[i].is_hanging)
206 {
207 for (j = 0; j < 2; ++j)
208 {
209 if (sides[i].is.hanging.is_ghost[j])
210 {
211 typename ::parallel::distributed::
212 Triangulation<dim, spacedim>::cell_iterator cell =
213 cell_from_quad(triangulation,
214 sides[i].treeid,
215 *(sides[i].is.hanging.quad[j]));
216 ::types::subdomain_id *subid =
217 static_cast<::types::subdomain_id *>(
218 sc_array_push(subids));
219 *subid = cell->subdomain_id();
220 }
221 }
222 }
223 }
224
225 if (!subids->elem_count)
226 {
227 return;
228 }
229
230 nsubs = static_cast<int>(subids->elem_count);
231 subdomain_ids =
232 reinterpret_cast<::types::subdomain_id *>(subids->array);
233
234 for (i = 0; i < nsides; ++i)
235 {
236 if (sides[i].is_hanging)
237 {
238 for (j = 0; j < 2; ++j)
239 {
240 if (!sides[i].is.hanging.is_ghost[j])
241 {
242 typename ::parallel::distributed::
243 Triangulation<dim, spacedim>::cell_iterator cell =
244 cell_from_quad(triangulation,
245 sides[i].treeid,
246 *(sides[i].is.hanging.quad[j]));
247
248 for (k = 0; k < nsubs; ++k)
249 {
250 (*vertices_with_ghost_neighbors)
251 [cell->vertex_index(
252 p8est_edge_corners[sides[i].edge][1 ^ j])]
253 .insert(subdomain_ids[k]);
254 }
255 }
256 }
257 }
258 }
259
260 subids->elem_count = 0;
261 }
262
266 template <int dim, int spacedim>
267 void
268 find_ghosts_face(
269 typename ::internal::p4est::iter<dim>::face_info *info,
270 void *user_data)
271 {
272 int i, j, k;
273 int nsides = info->sides.elem_count;
274 auto *sides = reinterpret_cast<
275 typename ::internal::p4est::iter<dim>::face_side *>(
276 info->sides.array);
277 FindGhosts<dim, spacedim> *fg =
278 static_cast<FindGhosts<dim, spacedim> *>(user_data);
279 sc_array_t *subids = fg->subids;
280 const ::parallel::distributed::Triangulation<dim, spacedim>
281 *triangulation = fg->triangulation;
282 int nsubs;
283 ::types::subdomain_id *subdomain_ids;
284 std::map<unsigned int, std::set<::types::subdomain_id>>
285 *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
286 int limit = (dim == 2) ? 2 : 4;
287
288 subids->elem_count = 0;
289 for (i = 0; i < nsides; ++i)
290 {
291 if (sides[i].is_hanging)
292 {
293 for (j = 0; j < limit; ++j)
294 {
295 if (sides[i].is.hanging.is_ghost[j])
296 {
297 typename ::parallel::distributed::
298 Triangulation<dim, spacedim>::cell_iterator cell =
299 cell_from_quad(triangulation,
300 sides[i].treeid,
301 *(sides[i].is.hanging.quad[j]));
302 ::types::subdomain_id *subid =
303 static_cast<::types::subdomain_id *>(
304 sc_array_push(subids));
305 *subid = cell->subdomain_id();
306 }
307 }
308 }
309 }
310
311 if (!subids->elem_count)
312 {
313 return;
314 }
315
316 nsubs = static_cast<int>(subids->elem_count);
317 subdomain_ids =
318 reinterpret_cast<::types::subdomain_id *>(subids->array);
319
320 for (i = 0; i < nsides; ++i)
321 {
322 if (sides[i].is_hanging)
323 {
324 for (j = 0; j < limit; ++j)
325 {
326 if (!sides[i].is.hanging.is_ghost[j])
327 {
328 typename ::parallel::distributed::
329 Triangulation<dim, spacedim>::cell_iterator cell =
330 cell_from_quad(triangulation,
331 sides[i].treeid,
332 *(sides[i].is.hanging.quad[j]));
333
334 for (k = 0; k < nsubs; ++k)
335 {
336 if (dim == 2)
337 {
338 (*vertices_with_ghost_neighbors)
339 [cell->vertex_index(
340 p4est_face_corners[sides[i].face]
341 [(limit - 1) ^ j])]
342 .insert(subdomain_ids[k]);
343 }
344 else
345 {
346 (*vertices_with_ghost_neighbors)
347 [cell->vertex_index(
348 p8est_face_corners[sides[i].face]
349 [(limit - 1) ^ j])]
350 .insert(subdomain_ids[k]);
351 }
352 }
353 }
354 }
355 }
356 }
357
358 subids->elem_count = 0;
359 }
360 } // namespace
361
362
363 int (&functions<2>::quadrant_compare)(const void *v1, const void *v2) =
364 p4est_quadrant_compare;
365
367 types<2>::quadrant c[]) =
368 p4est_quadrant_childrenv;
369
371 const types<2>::quadrant *q) =
372 p4est_quadrant_overlaps_tree;
373
375 int level,
376 std::uint64_t id) =
377 p4est_quadrant_set_morton;
378
379 void
384
386 const types<2>::quadrant *q2) =
387 p4est_quadrant_is_equal;
388
390 const types<2>::quadrant *q2) =
391 p4est_quadrant_is_sibling;
392
394 const types<2>::quadrant *q2) =
395 p4est_quadrant_is_ancestor;
396
398 int level) =
399 p4est_quadrant_ancestor_id;
400
402 const types<2>::locidx which_tree,
403 const types<2>::quadrant *q,
404 const int guess) =
405 p4est_comm_find_owner;
406
408 types<2>::topidx num_vertices,
409 types<2>::topidx num_trees,
410 types<2>::topidx num_corners,
411 types<2>::topidx num_vtt) = p4est_connectivity_new;
412
414 types<2>::topidx num_vertices,
415 types<2>::topidx num_trees,
416 types<2>::topidx num_corners,
417 const double *vertices,
418 const types<2>::topidx *ttv,
419 const types<2>::topidx *ttt,
420 const int8_t *ttf,
421 const types<2>::topidx *ttc,
422 const types<2>::topidx *coff,
423 const types<2>::topidx *ctt,
424 const int8_t *ctc) = p4est_connectivity_new_copy;
425
427 types<2>::topidx tree_left,
428 types<2>::topidx tree_right,
429 int face_left,
430 int face_right,
431 int orientation) =
432 p4est_connectivity_join_faces;
433
435 p4est_connectivity_t *connectivity) = p4est_connectivity_destroy;
436
438 MPI_Comm mpicomm,
439 types<2>::connectivity *connectivity,
440 types<2>::locidx min_quadrants,
441 int min_level,
442 int fill_uniform,
443 std::size_t data_size,
444 p4est_init_t init_fn,
445 void *user_pointer) = p4est_new_ext;
446
448 int copy_data) = p4est_copy;
449
450 void (&functions<2>::destroy)(types<2>::forest *p4est) = p4est_destroy;
451
453 int refine_recursive,
454 p4est_refine_t refine_fn,
455 p4est_init_t init_fn) = p4est_refine;
456
458 int coarsen_recursive,
459 p4est_coarsen_t coarsen_fn,
460 p4est_init_t init_fn) = p4est_coarsen;
461
464 p4est_init_t init_fn) = p4est_balance;
465
467 int partition_for_coarsening,
468 p4est_weight_t weight_fn) =
469 p4est_partition_ext;
470
471 void (&functions<2>::save)(const char *filename,
473 int save_data) = p4est_save;
474
476 const char *filename,
477 MPI_Comm mpicomm,
478 std::size_t data_size,
479 int load_data,
480 int autopartition,
481 int broadcasthead,
482 void *user_pointer,
483 types<2>::connectivity **p4est) = p4est_load_ext;
484
486 const char *filename,
487 types<2>::connectivity *connectivity) = p4est_connectivity_save;
488
490 types<2>::connectivity *connectivity) = p4est_connectivity_is_valid;
491
493 const char *filename,
494 std::size_t *length) = p4est_connectivity_load;
495
497 p4est_checksum;
498
500 p4est_geometry_t *,
501 const char *baseName) =
502 p4est_vtk_write_file;
503
506 p4est_ghost_new;
507
509 p4est_ghost_destroy;
510
512 std::size_t data_size,
513 p4est_init_t init_fn,
514 void *user_pointer) = p4est_reset_data;
515
517 p4est_memory_used;
518
520 types<2>::connectivity *p4est) = p4est_connectivity_memory_used;
521
522 constexpr unsigned int functions<2>::max_level;
523
524 void (&functions<2>::transfer_fixed)(const types<2>::gloidx *dest_gfq,
525 const types<2>::gloidx *src_gfq,
526 MPI_Comm mpicomm,
527 int tag,
528 void *dest_data,
529 const void *src_data,
530 std::size_t data_size) =
531 p4est_transfer_fixed;
532
534 const types<2>::gloidx *dest_gfq,
535 const types<2>::gloidx *src_gfq,
536 MPI_Comm mpicomm,
537 int tag,
538 void *dest_data,
539 const void *src_data,
540 std::size_t data_size) = p4est_transfer_fixed_begin;
541
543 p4est_transfer_fixed_end;
544
545 void (&functions<2>::transfer_custom)(const types<2>::gloidx *dest_gfq,
546 const types<2>::gloidx *src_gfq,
547 MPI_Comm mpicomm,
548 int tag,
549 void *dest_data,
550 const int *dest_sizes,
551 const void *src_data,
552 const int *src_sizes) =
553 p4est_transfer_custom;
554
556 const types<2>::gloidx *dest_gfq,
557 const types<2>::gloidx *src_gfq,
558 MPI_Comm mpicomm,
559 int tag,
560 void *dest_data,
561 const int *dest_sizes,
562 const void *src_data,
563 const int *src_sizes) = p4est_transfer_custom_begin;
564
566 p4est_transfer_custom_end;
567
570 int call_post,
573 sc_array_t *points) = p4est_search_partition;
574
576 types<2>::connectivity *connectivity,
577 types<2>::topidx treeid,
580 double vxyz[3]) = p4est_qcoord_to_vertex;
581
582 int (&functions<3>::quadrant_compare)(const void *v1, const void *v2) =
583 p8est_quadrant_compare;
584
586 types<3>::quadrant c[]) =
587 p8est_quadrant_childrenv;
588
590 const types<3>::quadrant *q) =
591 p8est_quadrant_overlaps_tree;
592
594 int level,
595 std::uint64_t id) =
596 p8est_quadrant_set_morton;
597
598 void
603
605 const types<3>::quadrant *q2) =
606 p8est_quadrant_is_equal;
607
609 const types<3>::quadrant *q2) =
610 p8est_quadrant_is_sibling;
611
613 const types<3>::quadrant *q2) =
614 p8est_quadrant_is_ancestor;
615
617 int level) =
618 p8est_quadrant_ancestor_id;
619
621 const types<3>::locidx which_tree,
622 const types<3>::quadrant *q,
623 const int guess) =
624 p8est_comm_find_owner;
625
627 types<3>::topidx num_vertices,
628 types<3>::topidx num_trees,
629 types<3>::topidx num_edges,
630 types<3>::topidx num_ett,
631 types<3>::topidx num_corners,
632 types<3>::topidx num_ctt) = p8est_connectivity_new;
633
635 types<3>::topidx num_vertices,
636 types<3>::topidx num_trees,
637 types<3>::topidx num_edges,
638 types<3>::topidx num_corners,
639 const double *vertices,
640 const types<3>::topidx *ttv,
641 const types<3>::topidx *ttt,
642 const int8_t *ttf,
643 const types<3>::topidx *tte,
644 const types<3>::topidx *eoff,
645 const types<3>::topidx *ett,
646 const int8_t *ete,
647 const types<3>::topidx *ttc,
648 const types<3>::topidx *coff,
649 const types<3>::topidx *ctt,
650 const int8_t *ctc) = p8est_connectivity_new_copy;
651
653 p8est_connectivity_t *connectivity) = p8est_connectivity_destroy;
654
656 types<3>::topidx tree_left,
657 types<3>::topidx tree_right,
658 int face_left,
659 int face_right,
660 int orientation) =
661 p8est_connectivity_join_faces;
662
664 MPI_Comm mpicomm,
665 types<3>::connectivity *connectivity,
666 types<3>::locidx min_quadrants,
667 int min_level,
668 int fill_uniform,
669 std::size_t data_size,
670 p8est_init_t init_fn,
671 void *user_pointer) = p8est_new_ext;
672
674 int copy_data) = p8est_copy;
675
676 void (&functions<3>::destroy)(types<3>::forest *p8est) = p8est_destroy;
677
679 int refine_recursive,
680 p8est_refine_t refine_fn,
681 p8est_init_t init_fn) = p8est_refine;
682
684 int coarsen_recursive,
685 p8est_coarsen_t coarsen_fn,
686 p8est_init_t init_fn) = p8est_coarsen;
687
690 p8est_init_t init_fn) = p8est_balance;
691
693 int partition_for_coarsening,
694 p8est_weight_t weight_fn) =
695 p8est_partition_ext;
696
697 void (&functions<3>::save)(const char *filename,
699 int save_data) = p8est_save;
700
702 const char *filename,
703 MPI_Comm mpicomm,
704 std::size_t data_size,
705 int load_data,
706 int autopartition,
707 int broadcasthead,
708 void *user_pointer,
709 types<3>::connectivity **p4est) = p8est_load_ext;
710
712 const char *filename,
713 types<3>::connectivity *connectivity) = p8est_connectivity_save;
714
716 types<3>::connectivity *connectivity) = p8est_connectivity_is_valid;
717
719 const char *filename,
720 std::size_t *length) = p8est_connectivity_load;
721
722 unsigned int (&functions<3>::checksum)(types<3>::forest *p8est) =
723 p8est_checksum;
724
726 p8est_geometry_t *,
727 const char *baseName) =
728 p8est_vtk_write_file;
729
732 p8est_ghost_new;
733
735 p8est_ghost_destroy;
736
738 std::size_t data_size,
739 p8est_init_t init_fn,
740 void *user_pointer) = p8est_reset_data;
741
743 p8est_memory_used;
744
746 types<3>::connectivity *p4est) = p8est_connectivity_memory_used;
747
748 constexpr unsigned int functions<3>::max_level;
749
750 void (&functions<3>::transfer_fixed)(const types<3>::gloidx *dest_gfq,
751 const types<3>::gloidx *src_gfq,
752 MPI_Comm mpicomm,
753 int tag,
754 void *dest_data,
755 const void *src_data,
756 std::size_t data_size) =
757 p8est_transfer_fixed;
758
760 const types<3>::gloidx *dest_gfq,
761 const types<3>::gloidx *src_gfq,
762 MPI_Comm mpicomm,
763 int tag,
764 void *dest_data,
765 const void *src_data,
766 std::size_t data_size) = p8est_transfer_fixed_begin;
767
769 p8est_transfer_fixed_end;
770
771 void (&functions<3>::transfer_custom)(const types<3>::gloidx *dest_gfq,
772 const types<3>::gloidx *src_gfq,
773 MPI_Comm mpicomm,
774 int tag,
775 void *dest_data,
776 const int *dest_sizes,
777 const void *src_data,
778 const int *src_sizes) =
779 p8est_transfer_custom;
780
782 const types<3>::gloidx *dest_gfq,
783 const types<3>::gloidx *src_gfq,
784 MPI_Comm mpicomm,
785 int tag,
786 void *dest_data,
787 const int *dest_sizes,
788 const void *src_data,
789 const int *src_sizes) = p8est_transfer_custom_begin;
790
792 p8est_transfer_custom_end;
793
796 int call_post,
799 sc_array_t *points) = p8est_search_partition;
800
802 types<3>::connectivity *connectivity,
803 types<3>::topidx treeid,
807 double vxyz[3]) = p8est_qcoord_to_vertex;
808
809 template <int dim>
810 void
812 const typename types<dim>::quadrant &p4est_cell,
813 typename types<dim>::quadrant (
815 {
816 for (unsigned int c = 0;
817 c < ::GeometryInfo<dim>::max_children_per_cell;
818 ++c)
819 functions<dim>::quadrant_init(p4est_children[c]);
820
821 functions<dim>::quadrant_childrenv(&p4est_cell, p4est_children);
822 }
823
824 template <int dim>
825 void
827 {
830 /*level=*/0,
831 /*index=*/0);
832 }
833
834 template <int dim>
835 bool
837 const typename types<dim>::quadrant &q2)
838 {
839 return functions<dim>::quadrant_is_equal(&q1, &q2);
840 }
841
842
843
844 template <int dim>
845 bool
847 const typename types<dim>::quadrant &q2)
848 {
850 }
851
852 template <int dim>
853 bool
854 tree_exists_locally(const typename types<dim>::forest *parallel_forest,
855 const typename types<dim>::topidx coarse_grid_cell)
856 {
857 Assert(coarse_grid_cell < parallel_forest->connectivity->num_trees,
859 return ((coarse_grid_cell >= parallel_forest->first_local_tree) &&
860 (coarse_grid_cell <= parallel_forest->last_local_tree));
861 }
862
863
864
865 // template specializations
866
867 template <>
869 copy_connectivity<2>(const typename types<2>::connectivity *connectivity)
870 {
872 connectivity->num_vertices,
873 connectivity->num_trees,
874 connectivity->num_corners,
875 connectivity->vertices,
876 connectivity->tree_to_vertex,
877 connectivity->tree_to_tree,
878 connectivity->tree_to_face,
879 connectivity->tree_to_corner,
880 connectivity->ctt_offset,
881 connectivity->corner_to_tree,
882 connectivity->corner_to_corner);
883 }
884
885 template <>
887 copy_connectivity<3>(const typename types<3>::connectivity *connectivity)
888 {
890 connectivity->num_vertices,
891 connectivity->num_trees,
892 connectivity->num_edges,
893 connectivity->num_corners,
894 connectivity->vertices,
895 connectivity->tree_to_vertex,
896 connectivity->tree_to_tree,
897 connectivity->tree_to_face,
898 connectivity->tree_to_edge,
899 connectivity->ett_offset,
900 connectivity->edge_to_tree,
901 connectivity->edge_to_edge,
902 connectivity->tree_to_corner,
903 connectivity->ctt_offset,
904 connectivity->corner_to_tree,
905 connectivity->corner_to_corner);
906 }
907
908
909
910 template <>
911 bool
912 quadrant_is_equal<1>(const typename types<1>::quadrant &q1,
913 const typename types<1>::quadrant &q2)
914 {
915 return q1 == q2;
916 }
917
918
919
920 template <>
921 bool
923 const types<1>::quadrant &q2)
924 {
925 // determine level of quadrants
926 const int level_1 = (q1 << types<1>::max_n_child_indices_bits) >>
928 const int level_2 = (q2 << types<1>::max_n_child_indices_bits) >>
930
931 // q1 can be an ancestor of q2 if q1's level is smaller
932 if (level_1 >= level_2)
933 return false;
934
935 // extract path of quadrants up to level of possible ancestor q1
936 const int truncated_id_1 = (q1 >> (types<1>::n_bits - 1 - level_1))
937 << (types<1>::n_bits - 1 - level_1);
938 const int truncated_id_2 = (q2 >> (types<1>::n_bits - 1 - level_1))
939 << (types<1>::n_bits - 1 - level_1);
940
941 // compare paths
942 return truncated_id_1 == truncated_id_2;
943 }
944
945
946
947 template <>
948 void
950 const typename types<1>::quadrant &q,
951 typename types<1>::quadrant (
953 {
954 // determine the current level of quadrant
955 const int level_parent = (q << types<1>::max_n_child_indices_bits) >>
957 const int level_child = level_parent + 1;
958
959 // left child: only n_child_indices has to be incremented
960 p4est_children[0] = (q + 1);
961
962 // right child: increment and set a bit to 1 indicating that it is a right
963 // child
964 p4est_children[1] = (q + 1) | (1 << (types<1>::n_bits - 1 - level_child));
965 }
966
967
968
969 template <>
970 void
972 {
973 quad = 0;
974 }
975
976 } // namespace p4est
977} // namespace internal
978
979#endif // DEAL_II_WITH_P4EST
980
981/*-------------- Explicit Instantiations -------------------------------*/
982#include "distributed/p4est_wrappers.inst"
983
984
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
void init_quadrant_children< 1 >(const typename types< 1 >::quadrant &q, typename types< 1 >::quadrant(&p4est_children)[::GeometryInfo< 1 >::max_children_per_cell])
void init_coarse_quadrant(typename types< dim >::quadrant &quad)
bool quadrant_is_equal< 1 >(const typename types< 1 >::quadrant &q1, const typename types< 1 >::quadrant &q2)
types< 2 >::connectivity * copy_connectivity< 2 >(const typename types< 2 >::connectivity *connectivity)
bool quadrant_is_ancestor< 1 >(const types< 1 >::quadrant &q1, const types< 1 >::quadrant &q2)
types< 3 >::connectivity * copy_connectivity< 3 >(const typename types< 3 >::connectivity *connectivity)
bool quadrant_is_equal(const typename types< dim >::quadrant &q1, const typename types< dim >::quadrant &q2)
void init_quadrant_children(const typename types< dim >::quadrant &p4est_cell, typename types< dim >::quadrant(&p4est_children)[::GeometryInfo< dim >::max_children_per_cell])
bool quadrant_is_ancestor(const typename types< dim >::quadrant &q1, const typename types< dim >::quadrant &q2)
void init_coarse_quadrant< 1 >(typename types< 1 >::quadrant &quad)
bool tree_exists_locally(const typename types< dim >::forest *parallel_forest, const typename types< dim >::topidx coarse_grid_cell)
unsigned int subdomain_id
Definition types.h:52
unsigned int global_dof_index
Definition types.h:94
#define P8EST_QUADRANT_INIT(q)
#define P4EST_QUADRANT_INIT(q)
static constexpr unsigned int max_children_per_cell