19#ifdef DEAL_II_WITH_P4EST
22# include <sc_containers.h>
31# ifndef P4EST_QUADRANT_INIT
32# define P4EST_QUADRANT_INIT(q) \
33 ((void)std::memset((q), -1, sizeof(p4est_quadrant_t)))
36# ifndef P8EST_QUADRANT_INIT
37# define P8EST_QUADRANT_INIT(q) \
38 ((void)std::memset((q), -1, sizeof(p8est_quadrant_t)))
46#ifdef DEAL_II_WITH_P4EST
54 template <
int dim,
int spacedim>
55 typename ::Triangulation<dim, spacedim>::cell_iterator
57 const ::parallel::distributed::Triangulation<dim, spacedim>
59 const typename ::internal::p4est::types<dim>::topidx treeidx,
60 const typename ::internal::p4est::types<dim>::quadrant &quad)
62 int i, l = quad.level;
64 triangulation->get_p4est_tree_to_coarse_cell_permutation()[treeidx];
66 for (i = 0; i < l; ++i)
68 typename ::Triangulation<dim, spacedim>::cell_iterator cell(
69 triangulation, i, dealii_index);
73 Assert(cell->has_children(),
74 ExcMessage(
"p4est quadrant does not correspond to a cell!"));
75 dealii_index = cell->child_index(child_id);
78 typename ::Triangulation<dim, spacedim>::cell_iterator out_cell(
79 triangulation, l, dealii_index);
88 template <
int dim,
int spacedim>
91 const typename ::parallel::distributed::Triangulation<dim,
95 std::map<unsigned int, std::set<::types::subdomain_id>>
96 *vertices_with_ghost_neighbors;
105 template <
int dim,
int spacedim>
108 typename ::internal::p4est::iter<dim>::corner_info *info,
112 int nsides = info->sides.elem_count;
113 auto *sides =
reinterpret_cast<
114 typename ::internal::p4est::iter<dim>::corner_side *
>(
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;
123 std::map<unsigned int, std::set<::types::subdomain_id>>
124 *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
126 subids->elem_count = 0;
127 for (i = 0; i < nsides; ++i)
129 if (sides[i].is_ghost)
131 typename ::parallel::distributed::
132 Triangulation<dim, spacedim>::cell_iterator cell =
133 cell_from_quad(triangulation,
137 ExcMessage(
"ghost quad did not find ghost cell"));
140 sc_array_push(subids));
141 *subid = cell->subdomain_id();
145 if (!subids->elem_count)
150 nsubs =
static_cast<int>(subids->elem_count);
154 for (i = 0; i < nsides; ++i)
156 if (!sides[i].is_ghost)
158 typename ::parallel::distributed::
159 Triangulation<dim, spacedim>::cell_iterator cell =
160 cell_from_quad(triangulation,
167 for (j = 0; j < nsubs; ++j)
169 (*vertices_with_ghost_neighbors)[cell->vertex_index(
171 .insert(subdomain_ids[j]);
176 subids->elem_count = 0;
182 template <
int dim,
int spacedim>
185 typename ::internal::p4est::iter<dim>::edge_info *info,
189 int nsides = info->sides.elem_count;
190 auto *sides =
reinterpret_cast<
191 typename ::internal::p4est::iter<dim>::edge_side *
>(
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;
199 std::map<unsigned int, std::set<::types::subdomain_id>>
200 *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
202 subids->elem_count = 0;
203 for (i = 0; i < nsides; ++i)
205 if (sides[i].is_hanging)
207 for (j = 0; j < 2; ++j)
209 if (sides[i].is.hanging.is_ghost[j])
211 typename ::parallel::distributed::
212 Triangulation<dim, spacedim>::cell_iterator cell =
213 cell_from_quad(triangulation,
215 *(sides[i].is.hanging.quad[j]));
218 sc_array_push(subids));
219 *subid = cell->subdomain_id();
225 if (!subids->elem_count)
230 nsubs =
static_cast<int>(subids->elem_count);
234 for (i = 0; i < nsides; ++i)
236 if (sides[i].is_hanging)
238 for (j = 0; j < 2; ++j)
240 if (!sides[i].is.hanging.is_ghost[j])
242 typename ::parallel::distributed::
243 Triangulation<dim, spacedim>::cell_iterator cell =
244 cell_from_quad(triangulation,
246 *(sides[i].is.hanging.quad[j]));
248 for (k = 0; k < nsubs; ++k)
250 (*vertices_with_ghost_neighbors)
252 p8est_edge_corners[sides[i].edge][1 ^ j])]
253 .insert(subdomain_ids[k]);
260 subids->elem_count = 0;
266 template <
int dim,
int spacedim>
269 typename ::internal::p4est::iter<dim>::face_info *info,
273 int nsides = info->sides.elem_count;
274 auto *sides =
reinterpret_cast<
275 typename ::internal::p4est::iter<dim>::face_side *
>(
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;
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;
288 subids->elem_count = 0;
289 for (i = 0; i < nsides; ++i)
291 if (sides[i].is_hanging)
293 for (j = 0; j < limit; ++j)
295 if (sides[i].is.hanging.is_ghost[j])
297 typename ::parallel::distributed::
298 Triangulation<dim, spacedim>::cell_iterator cell =
299 cell_from_quad(triangulation,
301 *(sides[i].is.hanging.quad[j]));
304 sc_array_push(subids));
305 *subid = cell->subdomain_id();
311 if (!subids->elem_count)
316 nsubs =
static_cast<int>(subids->elem_count);
320 for (i = 0; i < nsides; ++i)
322 if (sides[i].is_hanging)
324 for (j = 0; j < limit; ++j)
326 if (!sides[i].is.hanging.is_ghost[j])
328 typename ::parallel::distributed::
329 Triangulation<dim, spacedim>::cell_iterator cell =
330 cell_from_quad(triangulation,
332 *(sides[i].is.hanging.quad[j]));
334 for (k = 0; k < nsubs; ++k)
338 (*vertices_with_ghost_neighbors)
340 p4est_face_corners[sides[i].face]
342 .insert(subdomain_ids[k]);
346 (*vertices_with_ghost_neighbors)
348 p8est_face_corners[sides[i].face]
350 .insert(subdomain_ids[k]);
358 subids->elem_count = 0;
364 p4est_quadrant_compare;
368 p4est_quadrant_childrenv;
372 p4est_quadrant_overlaps_tree;
377 p4est_quadrant_set_morton;
387 p4est_quadrant_is_equal;
391 p4est_quadrant_is_sibling;
395 p4est_quadrant_is_ancestor;
399 p4est_quadrant_ancestor_id;
405 p4est_comm_find_owner;
417 const double *vertices,
424 const int8_t *ctc) = p4est_connectivity_new_copy;
432 p4est_connectivity_join_faces;
435 p4est_connectivity_t *connectivity) = p4est_connectivity_destroy;
443 std::size_t data_size,
444 p4est_init_t init_fn,
445 void *user_pointer) = p4est_new_ext;
448 int copy_data) = p4est_copy;
453 int refine_recursive,
454 p4est_refine_t refine_fn,
455 p4est_init_t init_fn) = p4est_refine;
458 int coarsen_recursive,
459 p4est_coarsen_t coarsen_fn,
460 p4est_init_t init_fn) = p4est_coarsen;
464 p4est_init_t init_fn) = p4est_balance;
467 int partition_for_coarsening,
468 p4est_weight_t weight_fn) =
473 int save_data) = p4est_save;
476 const char *filename,
478 std::size_t data_size,
486 const char *filename,
493 const char *filename,
494 std::size_t *length) = p4est_connectivity_load;
501 const char *baseName) =
502 p4est_vtk_write_file;
512 std::size_t data_size,
513 p4est_init_t init_fn,
514 void *user_pointer) = p4est_reset_data;
529 const void *src_data,
530 std::size_t data_size) =
531 p4est_transfer_fixed;
539 const void *src_data,
540 std::size_t data_size) = p4est_transfer_fixed_begin;
543 p4est_transfer_fixed_end;
550 const int *dest_sizes,
551 const void *src_data,
552 const int *src_sizes) =
553 p4est_transfer_custom;
561 const int *dest_sizes,
562 const void *src_data,
563 const int *src_sizes) = p4est_transfer_custom_begin;
566 p4est_transfer_custom_end;
573 sc_array_t *points) = p4est_search_partition;
580 double vxyz[3]) = p4est_qcoord_to_vertex;
583 p8est_quadrant_compare;
587 p8est_quadrant_childrenv;
591 p8est_quadrant_overlaps_tree;
596 p8est_quadrant_set_morton;
606 p8est_quadrant_is_equal;
610 p8est_quadrant_is_sibling;
614 p8est_quadrant_is_ancestor;
618 p8est_quadrant_ancestor_id;
624 p8est_comm_find_owner;
639 const double *vertices,
650 const int8_t *ctc) = p8est_connectivity_new_copy;
653 p8est_connectivity_t *connectivity) = p8est_connectivity_destroy;
661 p8est_connectivity_join_faces;
669 std::size_t data_size,
670 p8est_init_t init_fn,
671 void *user_pointer) = p8est_new_ext;
674 int copy_data) = p8est_copy;
679 int refine_recursive,
680 p8est_refine_t refine_fn,
681 p8est_init_t init_fn) = p8est_refine;
684 int coarsen_recursive,
685 p8est_coarsen_t coarsen_fn,
686 p8est_init_t init_fn) = p8est_coarsen;
690 p8est_init_t init_fn) = p8est_balance;
693 int partition_for_coarsening,
694 p8est_weight_t weight_fn) =
699 int save_data) = p8est_save;
702 const char *filename,
704 std::size_t data_size,
712 const char *filename,
719 const char *filename,
720 std::size_t *length) = p8est_connectivity_load;
727 const char *baseName) =
728 p8est_vtk_write_file;
738 std::size_t data_size,
739 p8est_init_t init_fn,
740 void *user_pointer) = p8est_reset_data;
755 const void *src_data,
756 std::size_t data_size) =
757 p8est_transfer_fixed;
765 const void *src_data,
766 std::size_t data_size) = p8est_transfer_fixed_begin;
769 p8est_transfer_fixed_end;
776 const int *dest_sizes,
777 const void *src_data,
778 const int *src_sizes) =
779 p8est_transfer_custom;
787 const int *dest_sizes,
788 const void *src_data,
789 const int *src_sizes) = p8est_transfer_custom_begin;
792 p8est_transfer_custom_end;
799 sc_array_t *points) = p8est_search_partition;
807 double vxyz[3]) = p8est_qcoord_to_vertex;
816 for (
unsigned int c = 0;
817 c < ::GeometryInfo<dim>::max_children_per_cell;
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));
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);
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);
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) >>
932 if (level_1 >= level_2)
942 return truncated_id_1 == truncated_id_2;
955 const int level_parent = (q << types<1>::max_n_child_indices_bits) >>
957 const int level_child = level_parent + 1;
960 p4est_children[0] = (q + 1);
982#include "distributed/p4est_wrappers.inst"
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#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
unsigned int global_dof_index
#define P8EST_QUADRANT_INIT(q)
#define P4EST_QUADRANT_INIT(q)
static constexpr unsigned int max_children_per_cell