22 #ifdef DEAL_II_WITH_P4EST
30 template <
int dim,
int spacedim>
31 typename ::Triangulation<dim, spacedim>::cell_iterator
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)
38 int i,
l = quad.level;
40 triangulation->get_p4est_tree_to_coarse_cell_permutation()[treeidx];
42 for (i = 0; i <
l; i++)
44 typename ::Triangulation<dim, spacedim>::cell_iterator cell(
49 Assert(cell->has_children(),
50 ExcMessage(
"p4est quadrant does not correspond to a cell!"));
51 dealii_index = cell->child_index(child_id);
54 typename ::Triangulation<dim, spacedim>::cell_iterator out_cell(
64 template <
int dim,
int spacedim>
71 std::map<unsigned int, std::set<::types::subdomain_id>>
81 template <
int dim,
int spacedim>
84 typename ::internal::p4est::iter<dim>::corner_info *info,
88 int nsides = info->sides.elem_count;
89 auto *sides =
reinterpret_cast<
90 typename ::internal::p4est::iter<dim>::corner_side *
>(
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>
99 std::map<unsigned int, std::set<::types::subdomain_id>>
103 for (i = 0; i < nsides; i++)
105 if (sides[i].is_ghost)
107 typename ::parallel::distributed::
108 Triangulation<dim, spacedim>::cell_iterator cell =
113 ExcMessage(
"ghost quad did not find ghost cell"));
117 *subid = cell->subdomain_id();
126 nsubs =
static_cast<int>(
subids->elem_count);
130 for (i = 0; i < nsides; i++)
132 if (!sides[i].is_ghost)
134 typename ::parallel::distributed::
135 Triangulation<dim, spacedim>::cell_iterator cell =
143 for (j = 0; j < nsubs; j++)
145 (*vertices_with_ghost_neighbors)[cell->vertex_index(
147 .insert(subdomain_ids[j]);
158 template <
int dim,
int spacedim>
161 typename ::internal::p4est::iter<dim>::edge_info *info,
165 int nsides = info->sides.elem_count;
166 auto *sides =
reinterpret_cast<
167 typename ::internal::p4est::iter<dim>::edge_side *
>(
169 auto * fg =
static_cast<FindGhosts<dim, spacedim> *
>(user_data);
170 sc_array_t *
subids = fg->subids;
171 const ::parallel::distributed::Triangulation<dim, spacedim>
175 std::map<unsigned int, std::set<::types::subdomain_id>>
179 for (i = 0; i < nsides; i++)
181 if (sides[i].is_hanging)
183 for (j = 0; j < 2; j++)
185 if (sides[i].is.hanging.is_ghost[j])
187 typename ::parallel::distributed::
188 Triangulation<dim, spacedim>::cell_iterator cell =
191 *(sides[i].is.hanging.quad[j]));
195 *subid = cell->subdomain_id();
206 nsubs =
static_cast<int>(
subids->elem_count);
210 for (i = 0; i < nsides; i++)
212 if (sides[i].is_hanging)
214 for (j = 0; j < 2; j++)
216 if (!sides[i].is.hanging.is_ghost[j])
218 typename ::parallel::distributed::
219 Triangulation<dim, spacedim>::cell_iterator cell =
222 *(sides[i].is.hanging.quad[j]));
224 for (k = 0; k < nsubs; k++)
226 (*vertices_with_ghost_neighbors)
228 p8est_edge_corners[sides[i].edge][1 ^ j])]
229 .insert(subdomain_ids[k]);
242 template <
int dim,
int spacedim>
245 typename ::internal::p4est::iter<dim>::face_info *info,
249 int nsides = info->sides.elem_count;
250 auto *sides =
reinterpret_cast<
251 typename ::internal::p4est::iter<dim>::face_side *
>(
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>
260 std::map<unsigned int, std::set<::types::subdomain_id>>
262 int limit = (dim == 2) ? 2 : 4;
265 for (i = 0; i < nsides; i++)
267 if (sides[i].is_hanging)
269 for (j = 0; j < limit; j++)
271 if (sides[i].is.hanging.is_ghost[j])
273 typename ::parallel::distributed::
274 Triangulation<dim, spacedim>::cell_iterator cell =
277 *(sides[i].is.hanging.quad[j]));
281 *subid = cell->subdomain_id();
292 nsubs =
static_cast<int>(
subids->elem_count);
296 for (i = 0; i < nsides; i++)
298 if (sides[i].is_hanging)
300 for (j = 0; j < limit; j++)
302 if (!sides[i].is.hanging.is_ghost[j])
304 typename ::parallel::distributed::
305 Triangulation<dim, spacedim>::cell_iterator cell =
308 *(sides[i].is.hanging.quad[j]));
310 for (k = 0; k < nsubs; k++)
314 (*vertices_with_ghost_neighbors)
316 p4est_face_corners[sides[i].face]
318 .insert(subdomain_ids[k]);
322 (*vertices_with_ghost_neighbors)
324 p8est_face_corners[sides[i].face]
326 .insert(subdomain_ids[k]);
340 p4est_quadrant_compare;
344 p4est_quadrant_childrenv;
348 p4est_quadrant_overlaps_tree;
353 p4est_quadrant_set_morton;
357 p4est_quadrant_is_equal;
361 p4est_quadrant_is_sibling;
365 p4est_quadrant_is_ancestor;
369 p4est_quadrant_ancestor_id;
375 p4est_comm_find_owner;
389 p4est_connectivity_join_faces;
392 p4est_connectivity_t *connectivity) = p4est_connectivity_destroy;
400 std::size_t data_size,
401 p4est_init_t init_fn,
402 void * user_pointer) = p4est_new_ext;
407 int refine_recursive,
408 p4est_refine_t refine_fn,
409 p4est_init_t init_fn) = p4est_refine;
412 int coarsen_recursive,
413 p4est_coarsen_t coarsen_fn,
414 p4est_init_t init_fn) = p4est_coarsen;
418 p4est_init_t init_fn) = p4est_balance;
421 int partition_for_coarsening,
422 p4est_weight_t weight_fn) =
427 int save_data) = p4est_save;
430 const char * filename,
432 std::size_t data_size,
440 const char * filename,
447 const char * filename,
448 std::size_t *length) = p4est_connectivity_load;
455 const char *baseName) =
456 p4est_vtk_write_file;
466 std::size_t data_size,
467 p4est_init_t init_fn,
468 void *user_pointer) = p4est_reset_data;
483 const void * src_data,
484 std::size_t data_size) =
485 p4est_transfer_fixed;
493 const void * src_data,
494 std::size_t data_size) = p4est_transfer_fixed_begin;
497 p4est_transfer_fixed_end;
504 const int * dest_sizes,
505 const void * src_data,
506 const int * src_sizes) =
507 p4est_transfer_custom;
515 const int * dest_sizes,
516 const void * src_data,
517 const int * src_sizes) = p4est_transfer_custom_begin;
520 p4est_transfer_custom_end;
525 p8est_quadrant_compare;
529 p8est_quadrant_childrenv;
533 p8est_quadrant_overlaps_tree;
538 p8est_quadrant_set_morton;
542 p8est_quadrant_is_equal;
546 p8est_quadrant_is_sibling;
550 p8est_quadrant_is_ancestor;
554 p8est_quadrant_ancestor_id;
560 p8est_comm_find_owner;
571 p8est_connectivity_t *connectivity) = p8est_connectivity_destroy;
579 p8est_connectivity_join_faces;
587 std::size_t data_size,
588 p8est_init_t init_fn,
589 void * user_pointer) = p8est_new_ext;
594 int refine_recursive,
595 p8est_refine_t refine_fn,
596 p8est_init_t init_fn) = p8est_refine;
599 int coarsen_recursive,
600 p8est_coarsen_t coarsen_fn,
601 p8est_init_t init_fn) = p8est_coarsen;
605 p8est_init_t init_fn) = p8est_balance;
608 int partition_for_coarsening,
609 p8est_weight_t weight_fn) =
614 int save_data) = p8est_save;
617 const char * filename,
619 std::size_t data_size,
627 const char * filename,
634 const char * filename,
635 std::size_t *length) = p8est_connectivity_load;
642 const char *baseName) =
643 p8est_vtk_write_file;
653 std::size_t data_size,
654 p8est_init_t init_fn,
655 void *user_pointer) = p8est_reset_data;
670 const void * src_data,
671 std::size_t data_size) =
672 p8est_transfer_fixed;
680 const void * src_data,
681 std::size_t data_size) = p8est_transfer_fixed_begin;
684 p8est_transfer_fixed_end;
691 const int * dest_sizes,
692 const void * src_data,
693 const int * src_sizes) =
694 p8est_transfer_custom;
702 const int * dest_sizes,
703 const void * src_data,
704 const int * src_sizes) = p8est_transfer_custom_begin;
707 p8est_transfer_custom_end;
718 for (
unsigned int c = 0;
719 c < ::GeometryInfo<dim>::max_children_per_cell;
724 P4EST_QUADRANT_INIT(&p4est_children[c]);
727 P8EST_QUADRANT_INIT(&p4est_children[c]);
744 P4EST_QUADRANT_INIT(&quad);
747 P8EST_QUADRANT_INIT(&quad);
780 Assert(coarse_grid_cell < parallel_forest->connectivity->num_trees,
782 return ((coarse_grid_cell >= parallel_forest->first_local_tree) &&
783 (coarse_grid_cell <= parallel_forest->last_local_tree));
803 const int level_1 = (q1 << types<1>::max_n_child_indices_bits) >>
805 const int level_2 = (q2 << types<1>::max_n_child_indices_bits) >>
809 if (level_1 >= level_2)
819 return truncated_id_1 == truncated_id_2;
832 const int level_parent = (q << types<1>::max_n_child_indices_bits) >>
834 const int level_child = level_parent + 1;
837 p4est_children[0] = (q + 1);
855 #endif // DEAL_II_WITH_P4EST
858 #include "p4est_wrappers.inst"