Reference documentation for deal.II version 9.3.3
\(\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\}}\)
p4est_wrappers.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2008 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16
19
21
22#ifdef DEAL_II_WITH_P4EST
23
24namespace internal
25{
26 namespace p4est
27 {
28 namespace
29 {
30 template <int dim, int spacedim>
31 typename ::Triangulation<dim, spacedim>::cell_iterator
32 cell_from_quad(
33 const ::parallel::distributed::Triangulation<dim, spacedim>
35 const typename ::internal::p4est::types<dim>::topidx treeidx,
36 const typename ::internal::p4est::types<dim>::quadrant &quad)
37 {
38 int i, l = quad.level;
39 ::types::global_dof_index dealii_index =
40 triangulation->get_p4est_tree_to_coarse_cell_permutation()[treeidx];
41
42 for (i = 0; i < l; i++)
43 {
44 typename ::Triangulation<dim, spacedim>::cell_iterator cell(
45 triangulation, i, dealii_index);
46 const int child_id =
48 &quad, i + 1);
49 Assert(cell->has_children(),
50 ExcMessage("p4est quadrant does not correspond to a cell!"));
51 dealii_index = cell->child_index(child_id);
52 }
53
54 typename ::Triangulation<dim, spacedim>::cell_iterator out_cell(
55 triangulation, l, dealii_index);
56
57 return out_cell;
58 }
59
64 template <int dim, int spacedim>
65 struct FindGhosts
66 {
68 spacedim>
70 sc_array_t *subids;
71 std::map<unsigned int, std::set<::types::subdomain_id>>
73 };
74
75
81 template <int dim, int spacedim>
82 void
83 find_ghosts_corner(
84 typename ::internal::p4est::iter<dim>::corner_info *info,
85 void * user_data)
86 {
87 int i, j;
88 int nsides = info->sides.elem_count;
89 auto *sides = reinterpret_cast<
90 typename ::internal::p4est::iter<dim>::corner_side *>(
91 info->sides.array);
92 FindGhosts<dim, spacedim> *fg =
93 static_cast<FindGhosts<dim, spacedim> *>(user_data);
94 sc_array_t *subids = fg->subids;
95 const ::parallel::distributed::Triangulation<dim, spacedim>
96 * triangulation = fg->triangulation;
97 int nsubs;
98 ::types::subdomain_id *subdomain_ids;
99 std::map<unsigned int, std::set<::types::subdomain_id>>
100 *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
101
102 subids->elem_count = 0;
103 for (i = 0; i < nsides; i++)
104 {
105 if (sides[i].is_ghost)
106 {
107 typename ::parallel::distributed::
108 Triangulation<dim, spacedim>::cell_iterator cell =
109 cell_from_quad(triangulation,
110 sides[i].treeid,
111 *(sides[i].quad));
112 Assert(cell->is_ghost(),
113 ExcMessage("ghost quad did not find ghost cell"));
114 ::types::subdomain_id *subid =
115 static_cast<::types::subdomain_id *>(
116 sc_array_push(subids));
117 *subid = cell->subdomain_id();
118 }
119 }
120
121 if (!subids->elem_count)
122 {
123 return;
124 }
125
126 nsubs = static_cast<int>(subids->elem_count);
127 subdomain_ids =
128 reinterpret_cast<::types::subdomain_id *>(subids->array);
129
130 for (i = 0; i < nsides; i++)
131 {
132 if (!sides[i].is_ghost)
133 {
134 typename ::parallel::distributed::
135 Triangulation<dim, spacedim>::cell_iterator cell =
136 cell_from_quad(triangulation,
137 sides[i].treeid,
138 *(sides[i].quad));
139
140 Assert(!cell->is_ghost(),
141 ExcMessage("local quad found ghost cell"));
142
143 for (j = 0; j < nsubs; j++)
144 {
145 (*vertices_with_ghost_neighbors)[cell->vertex_index(
146 sides[i].corner)]
147 .insert(subdomain_ids[j]);
148 }
149 }
150 }
151
152 subids->elem_count = 0;
153 }
154
158 template <int dim, int spacedim>
159 void
160 find_ghosts_edge(
161 typename ::internal::p4est::iter<dim>::edge_info *info,
162 void * user_data)
163 {
164 int i, j, k;
165 int nsides = info->sides.elem_count;
166 auto *sides = reinterpret_cast<
167 typename ::internal::p4est::iter<dim>::edge_side *>(
168 info->sides.array);
169 auto * fg = static_cast<FindGhosts<dim, spacedim> *>(user_data);
170 sc_array_t *subids = fg->subids;
171 const ::parallel::distributed::Triangulation<dim, spacedim>
172 * triangulation = fg->triangulation;
173 int nsubs;
174 ::types::subdomain_id *subdomain_ids;
175 std::map<unsigned int, std::set<::types::subdomain_id>>
176 *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
177
178 subids->elem_count = 0;
179 for (i = 0; i < nsides; i++)
180 {
181 if (sides[i].is_hanging)
182 {
183 for (j = 0; j < 2; j++)
184 {
185 if (sides[i].is.hanging.is_ghost[j])
186 {
187 typename ::parallel::distributed::
188 Triangulation<dim, spacedim>::cell_iterator cell =
189 cell_from_quad(triangulation,
190 sides[i].treeid,
191 *(sides[i].is.hanging.quad[j]));
192 ::types::subdomain_id *subid =
193 static_cast<::types::subdomain_id *>(
194 sc_array_push(subids));
195 *subid = cell->subdomain_id();
196 }
197 }
198 }
199 }
200
201 if (!subids->elem_count)
202 {
203 return;
204 }
205
206 nsubs = static_cast<int>(subids->elem_count);
207 subdomain_ids =
208 reinterpret_cast<::types::subdomain_id *>(subids->array);
209
210 for (i = 0; i < nsides; i++)
211 {
212 if (sides[i].is_hanging)
213 {
214 for (j = 0; j < 2; j++)
215 {
216 if (!sides[i].is.hanging.is_ghost[j])
217 {
218 typename ::parallel::distributed::
219 Triangulation<dim, spacedim>::cell_iterator cell =
220 cell_from_quad(triangulation,
221 sides[i].treeid,
222 *(sides[i].is.hanging.quad[j]));
223
224 for (k = 0; k < nsubs; k++)
225 {
226 (*vertices_with_ghost_neighbors)
227 [cell->vertex_index(
228 p8est_edge_corners[sides[i].edge][1 ^ j])]
229 .insert(subdomain_ids[k]);
230 }
231 }
232 }
233 }
234 }
235
236 subids->elem_count = 0;
237 }
238
242 template <int dim, int spacedim>
243 void
244 find_ghosts_face(
245 typename ::internal::p4est::iter<dim>::face_info *info,
246 void * user_data)
247 {
248 int i, j, k;
249 int nsides = info->sides.elem_count;
250 auto *sides = reinterpret_cast<
251 typename ::internal::p4est::iter<dim>::face_side *>(
252 info->sides.array);
253 FindGhosts<dim, spacedim> *fg =
254 static_cast<FindGhosts<dim, spacedim> *>(user_data);
255 sc_array_t *subids = fg->subids;
256 const ::parallel::distributed::Triangulation<dim, spacedim>
257 * triangulation = fg->triangulation;
258 int nsubs;
259 ::types::subdomain_id *subdomain_ids;
260 std::map<unsigned int, std::set<::types::subdomain_id>>
261 * vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
262 int limit = (dim == 2) ? 2 : 4;
263
264 subids->elem_count = 0;
265 for (i = 0; i < nsides; i++)
266 {
267 if (sides[i].is_hanging)
268 {
269 for (j = 0; j < limit; j++)
270 {
271 if (sides[i].is.hanging.is_ghost[j])
272 {
273 typename ::parallel::distributed::
274 Triangulation<dim, spacedim>::cell_iterator cell =
275 cell_from_quad(triangulation,
276 sides[i].treeid,
277 *(sides[i].is.hanging.quad[j]));
278 ::types::subdomain_id *subid =
279 static_cast<::types::subdomain_id *>(
280 sc_array_push(subids));
281 *subid = cell->subdomain_id();
282 }
283 }
284 }
285 }
286
287 if (!subids->elem_count)
288 {
289 return;
290 }
291
292 nsubs = static_cast<int>(subids->elem_count);
293 subdomain_ids =
294 reinterpret_cast<::types::subdomain_id *>(subids->array);
295
296 for (i = 0; i < nsides; i++)
297 {
298 if (sides[i].is_hanging)
299 {
300 for (j = 0; j < limit; j++)
301 {
302 if (!sides[i].is.hanging.is_ghost[j])
303 {
304 typename ::parallel::distributed::
305 Triangulation<dim, spacedim>::cell_iterator cell =
306 cell_from_quad(triangulation,
307 sides[i].treeid,
308 *(sides[i].is.hanging.quad[j]));
309
310 for (k = 0; k < nsubs; k++)
311 {
312 if (dim == 2)
313 {
314 (*vertices_with_ghost_neighbors)
315 [cell->vertex_index(
316 p4est_face_corners[sides[i].face]
317 [(limit - 1) ^ j])]
318 .insert(subdomain_ids[k]);
319 }
320 else
321 {
322 (*vertices_with_ghost_neighbors)
323 [cell->vertex_index(
324 p8est_face_corners[sides[i].face]
325 [(limit - 1) ^ j])]
326 .insert(subdomain_ids[k]);
327 }
328 }
329 }
330 }
331 }
332 }
333
334 subids->elem_count = 0;
335 }
336 } // namespace
337
338
339 int (&functions<2>::quadrant_compare)(const void *v1, const void *v2) =
340 p4est_quadrant_compare;
341
343 types<2>::quadrant c[]) =
344 p4est_quadrant_childrenv;
345
347 const types<2>::quadrant *q) =
348 p4est_quadrant_overlaps_tree;
349
351 int level,
352 uint64_t id) =
353 p4est_quadrant_set_morton;
354
356 const types<2>::quadrant *q2) =
357 p4est_quadrant_is_equal;
358
360 const types<2>::quadrant *q2) =
361 p4est_quadrant_is_sibling;
362
364 const types<2>::quadrant *q2) =
365 p4est_quadrant_is_ancestor;
366
368 int level) =
369 p4est_quadrant_ancestor_id;
370
372 const types<2>::locidx which_tree,
373 const types<2>::quadrant *q,
374 const int guess) =
375 p4est_comm_find_owner;
376
378 types<2>::topidx num_vertices,
379 types<2>::topidx num_trees,
380 types<2>::topidx num_corners,
381 types<2>::topidx num_vtt) = p4est_connectivity_new;
382
384 types<2>::topidx num_vertices,
385 types<2>::topidx num_trees,
386 types<2>::topidx num_corners,
387 const double * vertices,
388 const types<2>::topidx *ttv,
389 const types<2>::topidx *ttt,
390 const int8_t * ttf,
391 const types<2>::topidx *ttc,
392 const types<2>::topidx *coff,
393 const types<2>::topidx *ctt,
394 const int8_t * ctc) = p4est_connectivity_new_copy;
395
397 types<2>::topidx tree_left,
398 types<2>::topidx tree_right,
399 int face_left,
400 int face_right,
401 int orientation) =
402 p4est_connectivity_join_faces;
403
405 p4est_connectivity_t *connectivity) = p4est_connectivity_destroy;
406
408 MPI_Comm mpicomm,
409 types<2>::connectivity *connectivity,
410 types<2>::locidx min_quadrants,
411 int min_level,
412 int fill_uniform,
413 std::size_t data_size,
414 p4est_init_t init_fn,
415 void * user_pointer) = p4est_new_ext;
416
418 int copy_data) = p4est_copy;
419
420 void (&functions<2>::destroy)(types<2>::forest *p4est) = p4est_destroy;
421
423 int refine_recursive,
424 p4est_refine_t refine_fn,
425 p4est_init_t init_fn) = p4est_refine;
426
428 int coarsen_recursive,
429 p4est_coarsen_t coarsen_fn,
430 p4est_init_t init_fn) = p4est_coarsen;
431
434 p4est_init_t init_fn) = p4est_balance;
435
437 int partition_for_coarsening,
438 p4est_weight_t weight_fn) =
439 p4est_partition_ext;
440
441 void (&functions<2>::save)(const char * filename,
442 types<2>::forest *p4est,
443 int save_data) = p4est_save;
444
446 const char * filename,
447 MPI_Comm mpicomm,
448 std::size_t data_size,
449 int load_data,
450 int autopartition,
451 int broadcasthead,
452 void * user_pointer,
453 types<2>::connectivity **p4est) = p4est_load_ext;
454
456 const char * filename,
457 types<2>::connectivity *connectivity) = p4est_connectivity_save;
458
460 types<2>::connectivity *connectivity) = p4est_connectivity_is_valid;
461
463 const char * filename,
464 std::size_t *length) = p4est_connectivity_load;
465
466 unsigned int (&functions<2>::checksum)(types<2>::forest *p4est) =
467 p4est_checksum;
468
470 p4est_geometry_t *,
471 const char *baseName) =
472 p4est_vtk_write_file;
473
476 p4est_ghost_new;
477
479 p4est_ghost_destroy;
480
482 std::size_t data_size,
483 p4est_init_t init_fn,
484 void *user_pointer) = p4est_reset_data;
485
487 p4est_memory_used;
488
490 types<2>::connectivity *p4est) = p4est_connectivity_memory_used;
491
492 constexpr unsigned int functions<2>::max_level;
493
494 void (&functions<2>::transfer_fixed)(const types<2>::gloidx *dest_gfq,
495 const types<2>::gloidx *src_gfq,
496 MPI_Comm mpicomm,
497 int tag,
498 void * dest_data,
499 const void * src_data,
500 std::size_t data_size) =
501 p4est_transfer_fixed;
502
504 const types<2>::gloidx *dest_gfq,
505 const types<2>::gloidx *src_gfq,
506 MPI_Comm mpicomm,
507 int tag,
508 void * dest_data,
509 const void * src_data,
510 std::size_t data_size) = p4est_transfer_fixed_begin;
511
513 p4est_transfer_fixed_end;
514
515 void (&functions<2>::transfer_custom)(const types<2>::gloidx *dest_gfq,
516 const types<2>::gloidx *src_gfq,
517 MPI_Comm mpicomm,
518 int tag,
519 void * dest_data,
520 const int * dest_sizes,
521 const void * src_data,
522 const int * src_sizes) =
523 p4est_transfer_custom;
524
526 const types<2>::gloidx *dest_gfq,
527 const types<2>::gloidx *src_gfq,
528 MPI_Comm mpicomm,
529 int tag,
530 void * dest_data,
531 const int * dest_sizes,
532 const void * src_data,
533 const int * src_sizes) = p4est_transfer_custom_begin;
534
536 p4est_transfer_custom_end;
537
538
539
540 int (&functions<3>::quadrant_compare)(const void *v1, const void *v2) =
541 p8est_quadrant_compare;
542
544 types<3>::quadrant c[]) =
545 p8est_quadrant_childrenv;
546
548 const types<3>::quadrant *q) =
549 p8est_quadrant_overlaps_tree;
550
552 int level,
553 uint64_t id) =
554 p8est_quadrant_set_morton;
555
557 const types<3>::quadrant *q2) =
558 p8est_quadrant_is_equal;
559
561 const types<3>::quadrant *q2) =
562 p8est_quadrant_is_sibling;
563
565 const types<3>::quadrant *q2) =
566 p8est_quadrant_is_ancestor;
567
569 int level) =
570 p8est_quadrant_ancestor_id;
571
573 const types<3>::locidx which_tree,
574 const types<3>::quadrant *q,
575 const int guess) =
576 p8est_comm_find_owner;
577
579 types<3>::topidx num_vertices,
580 types<3>::topidx num_trees,
581 types<3>::topidx num_edges,
582 types<3>::topidx num_ett,
583 types<3>::topidx num_corners,
584 types<3>::topidx num_ctt) = p8est_connectivity_new;
585
587 types<3>::topidx num_vertices,
588 types<3>::topidx num_trees,
589 types<3>::topidx num_edges,
590 types<3>::topidx num_corners,
591 const double * vertices,
592 const types<3>::topidx *ttv,
593 const types<3>::topidx *ttt,
594 const int8_t * ttf,
595 const types<3>::topidx *tte,
596 const types<3>::topidx *eoff,
597 const types<3>::topidx *ett,
598 const int8_t * ete,
599 const types<3>::topidx *ttc,
600 const types<3>::topidx *coff,
601 const types<3>::topidx *ctt,
602 const int8_t * ctc) = p8est_connectivity_new_copy;
603
605 p8est_connectivity_t *connectivity) = p8est_connectivity_destroy;
606
608 types<3>::topidx tree_left,
609 types<3>::topidx tree_right,
610 int face_left,
611 int face_right,
612 int orientation) =
613 p8est_connectivity_join_faces;
614
616 MPI_Comm mpicomm,
617 types<3>::connectivity *connectivity,
618 types<3>::locidx min_quadrants,
619 int min_level,
620 int fill_uniform,
621 std::size_t data_size,
622 p8est_init_t init_fn,
623 void * user_pointer) = p8est_new_ext;
624
626 int copy_data) = p8est_copy;
627
628 void (&functions<3>::destroy)(types<3>::forest *p8est) = p8est_destroy;
629
631 int refine_recursive,
632 p8est_refine_t refine_fn,
633 p8est_init_t init_fn) = p8est_refine;
634
636 int coarsen_recursive,
637 p8est_coarsen_t coarsen_fn,
638 p8est_init_t init_fn) = p8est_coarsen;
639
642 p8est_init_t init_fn) = p8est_balance;
643
645 int partition_for_coarsening,
646 p8est_weight_t weight_fn) =
647 p8est_partition_ext;
648
649 void (&functions<3>::save)(const char * filename,
650 types<3>::forest *p4est,
651 int save_data) = p8est_save;
652
654 const char * filename,
655 MPI_Comm mpicomm,
656 std::size_t data_size,
657 int load_data,
658 int autopartition,
659 int broadcasthead,
660 void * user_pointer,
661 types<3>::connectivity **p4est) = p8est_load_ext;
662
664 const char * filename,
665 types<3>::connectivity *connectivity) = p8est_connectivity_save;
666
668 types<3>::connectivity *connectivity) = p8est_connectivity_is_valid;
669
671 const char * filename,
672 std::size_t *length) = p8est_connectivity_load;
673
674 unsigned int (&functions<3>::checksum)(types<3>::forest *p8est) =
675 p8est_checksum;
676
678 p8est_geometry_t *,
679 const char *baseName) =
680 p8est_vtk_write_file;
681
684 p8est_ghost_new;
685
687 p8est_ghost_destroy;
688
690 std::size_t data_size,
691 p8est_init_t init_fn,
692 void *user_pointer) = p8est_reset_data;
693
695 p8est_memory_used;
696
698 types<3>::connectivity *p4est) = p8est_connectivity_memory_used;
699
700 constexpr unsigned int functions<3>::max_level;
701
702 void (&functions<3>::transfer_fixed)(const types<3>::gloidx *dest_gfq,
703 const types<3>::gloidx *src_gfq,
704 MPI_Comm mpicomm,
705 int tag,
706 void * dest_data,
707 const void * src_data,
708 std::size_t data_size) =
709 p8est_transfer_fixed;
710
712 const types<3>::gloidx *dest_gfq,
713 const types<3>::gloidx *src_gfq,
714 MPI_Comm mpicomm,
715 int tag,
716 void * dest_data,
717 const void * src_data,
718 std::size_t data_size) = p8est_transfer_fixed_begin;
719
721 p8est_transfer_fixed_end;
722
723 void (&functions<3>::transfer_custom)(const types<3>::gloidx *dest_gfq,
724 const types<3>::gloidx *src_gfq,
725 MPI_Comm mpicomm,
726 int tag,
727 void * dest_data,
728 const int * dest_sizes,
729 const void * src_data,
730 const int * src_sizes) =
731 p8est_transfer_custom;
732
734 const types<3>::gloidx *dest_gfq,
735 const types<3>::gloidx *src_gfq,
736 MPI_Comm mpicomm,
737 int tag,
738 void * dest_data,
739 const int * dest_sizes,
740 const void * src_data,
741 const int * src_sizes) = p8est_transfer_custom_begin;
742
744 p8est_transfer_custom_end;
745
746
747
748 template <int dim>
749 void
751 const typename types<dim>::quadrant &p4est_cell,
752 typename types<dim>::quadrant (
754 {
755 for (unsigned int c = 0;
756 c < ::GeometryInfo<dim>::max_children_per_cell;
757 ++c)
758 switch (dim)
759 {
760 case 2:
761 P4EST_QUADRANT_INIT(&p4est_children[c]);
762 break;
763 case 3:
764 P8EST_QUADRANT_INIT(&p4est_children[c]);
765 break;
766 default:
767 Assert(false, ExcNotImplemented());
768 }
769
770
771 functions<dim>::quadrant_childrenv(&p4est_cell, p4est_children);
772 }
773
774 template <int dim>
775 void
777 {
778 switch (dim)
779 {
780 case 2:
781 P4EST_QUADRANT_INIT(&quad);
782 break;
783 case 3:
784 P8EST_QUADRANT_INIT(&quad);
785 break;
786 default:
787 Assert(false, ExcNotImplemented());
788 }
790 /*level=*/0,
791 /*index=*/0);
792 }
793
794 template <int dim>
795 bool
797 const typename types<dim>::quadrant &q2)
798 {
799 return functions<dim>::quadrant_is_equal(&q1, &q2);
800 }
801
802
803
804 template <int dim>
805 bool
807 const typename types<dim>::quadrant &q2)
808 {
810 }
811
812 template <int dim>
813 bool
814 tree_exists_locally(const typename types<dim>::forest *parallel_forest,
815 const typename types<dim>::topidx coarse_grid_cell)
816 {
817 Assert(coarse_grid_cell < parallel_forest->connectivity->num_trees,
819 return ((coarse_grid_cell >= parallel_forest->first_local_tree) &&
820 (coarse_grid_cell <= parallel_forest->last_local_tree));
821 }
822
823
824
825 // template specializations
826
827 template <>
829 copy_connectivity<2>(const typename types<2>::connectivity *connectivity)
830 {
832 connectivity->num_vertices,
833 connectivity->num_trees,
834 connectivity->num_corners,
835 connectivity->vertices,
836 connectivity->tree_to_vertex,
837 connectivity->tree_to_tree,
838 connectivity->tree_to_face,
839 connectivity->tree_to_corner,
840 connectivity->ctt_offset,
841 connectivity->corner_to_tree,
842 connectivity->corner_to_corner);
843 }
844
845 template <>
847 copy_connectivity<3>(const typename types<3>::connectivity *connectivity)
848 {
850 connectivity->num_vertices,
851 connectivity->num_trees,
852 connectivity->num_edges,
853 connectivity->num_corners,
854 connectivity->vertices,
855 connectivity->tree_to_vertex,
856 connectivity->tree_to_tree,
857 connectivity->tree_to_face,
858 connectivity->tree_to_edge,
859 connectivity->ett_offset,
860 connectivity->edge_to_tree,
861 connectivity->edge_to_edge,
862 connectivity->tree_to_corner,
863 connectivity->ctt_offset,
864 connectivity->corner_to_tree,
865 connectivity->corner_to_corner);
866 }
867
868
869
870 template <>
871 bool
872 quadrant_is_equal<1>(const typename types<1>::quadrant &q1,
873 const typename types<1>::quadrant &q2)
874 {
875 return q1 == q2;
876 }
877
878
879
880 template <>
882 types<1>::quadrant const &q2)
883 {
884 // determine level of quadrants
885 const int level_1 = (q1 << types<1>::max_n_child_indices_bits) >>
887 const int level_2 = (q2 << types<1>::max_n_child_indices_bits) >>
889
890 // q1 can be an ancestor of q2 if q1's level is smaller
891 if (level_1 >= level_2)
892 return false;
893
894 // extract path of quadrants up to level of possible ancestor q1
895 const int truncated_id_1 = (q1 >> (types<1>::n_bits - 1 - level_1))
896 << (types<1>::n_bits - 1 - level_1);
897 const int truncated_id_2 = (q2 >> (types<1>::n_bits - 1 - level_1))
898 << (types<1>::n_bits - 1 - level_1);
899
900 // compare paths
901 return truncated_id_1 == truncated_id_2;
902 }
903
904
905
906 template <>
907 void
909 const typename types<1>::quadrant &q,
910 typename types<1>::quadrant (
912 {
913 // determine the current level of quadrant
914 const int level_parent = (q << types<1>::max_n_child_indices_bits) >>
916 const int level_child = level_parent + 1;
917
918 // left child: only n_child_indices has to be incremented
919 p4est_children[0] = (q + 1);
920
921 // right child: increment and set a bit to 1 indicating that it is a right
922 // child
923 p4est_children[1] = (q + 1) | (1 << (types<1>::n_bits - 1 - level_child));
924 }
925
926
927
928 template <>
930 {
931 quad = 0;
932 }
933
934 } // namespace p4est
935} // namespace internal
936
937#endif // DEAL_II_WITH_P4EST
938
939/*-------------- Explicit Instantiations -------------------------------*/
940#include "p4est_wrappers.inst"
941
942
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
Point< 3 > vertices[4]
unsigned int level
Definition: grid_out.cc:4590
const unsigned int v1
Definition: grid_tools.cc:963
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
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)
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 quadrant_is_ancestor< 1 >(types< 1 >::quadrant const &q1, types< 1 >::quadrant const &q2)
bool tree_exists_locally(const typename types< dim >::forest *parallel_forest, const typename types< dim >::topidx coarse_grid_cell)
TriangulationBase< dim, spacedim > Triangulation
Definition: tria_base.h:396
sc_array_t * subids
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
std::map< unsigned int, std::set<::types::subdomain_id > > * vertices_with_ghost_neighbors