21#ifdef DEAL_II_WITH_P4EST
29 template <
int dim,
int spacedim>
30 typename ::Triangulation<dim, spacedim>::cell_iterator
32 const ::parallel::distributed::Triangulation<dim, spacedim>
34 const typename ::internal::p4est::types<dim>::topidx treeidx,
35 const typename ::internal::p4est::types<dim>::quadrant &quad)
37 int i,
l = quad.level;
39 triangulation->get_p4est_tree_to_coarse_cell_permutation()[treeidx];
41 for (i = 0; i <
l; ++i)
43 typename ::Triangulation<dim, spacedim>::cell_iterator cell(
48 Assert(cell->has_children(),
49 ExcMessage(
"p4est quadrant does not correspond to a cell!"));
50 dealii_index = cell->child_index(child_id);
53 typename ::Triangulation<dim, spacedim>::cell_iterator out_cell(
63 template <
int dim,
int spacedim>
66 const typename ::parallel::distributed::Triangulation<dim,
70 std::map<unsigned int, std::set<::types::subdomain_id>>
80 template <
int dim,
int spacedim>
83 typename ::internal::p4est::iter<dim>::corner_info *info,
87 int nsides = info->sides.elem_count;
88 auto *sides =
reinterpret_cast<
89 typename ::internal::p4est::iter<dim>::corner_side *
>(
91 FindGhosts<dim, spacedim> *fg =
92 static_cast<FindGhosts<dim, spacedim> *
>(user_data);
93 sc_array_t *
subids = fg->subids;
94 const ::parallel::distributed::Triangulation<dim, spacedim>
98 std::map<unsigned int, std::set<::types::subdomain_id>>
102 for (i = 0; i < nsides; ++i)
104 if (sides[i].is_ghost)
106 typename ::parallel::distributed::
107 Triangulation<dim, spacedim>::cell_iterator cell =
112 ExcMessage(
"ghost quad did not find ghost cell"));
116 *subid = cell->subdomain_id();
125 nsubs =
static_cast<int>(
subids->elem_count);
129 for (i = 0; i < nsides; ++i)
131 if (!sides[i].is_ghost)
133 typename ::parallel::distributed::
134 Triangulation<dim, spacedim>::cell_iterator cell =
142 for (j = 0; j < nsubs; ++j)
144 (*vertices_with_ghost_neighbors)[cell->vertex_index(
146 .insert(subdomain_ids[j]);
157 template <
int dim,
int spacedim>
160 typename ::internal::p4est::iter<dim>::edge_info *info,
164 int nsides = info->sides.elem_count;
165 auto *sides =
reinterpret_cast<
166 typename ::internal::p4est::iter<dim>::edge_side *
>(
168 auto *fg =
static_cast<FindGhosts<dim, spacedim> *
>(user_data);
169 sc_array_t *
subids = fg->subids;
170 const ::parallel::distributed::Triangulation<dim, spacedim>
174 std::map<unsigned int, std::set<::types::subdomain_id>>
178 for (i = 0; i < nsides; ++i)
180 if (sides[i].is_hanging)
182 for (j = 0; j < 2; ++j)
184 if (sides[i].is.hanging.is_ghost[j])
186 typename ::parallel::distributed::
187 Triangulation<dim, spacedim>::cell_iterator cell =
190 *(sides[i].is.hanging.quad[j]));
194 *subid = cell->subdomain_id();
205 nsubs =
static_cast<int>(
subids->elem_count);
209 for (i = 0; i < nsides; ++i)
211 if (sides[i].is_hanging)
213 for (j = 0; j < 2; ++j)
215 if (!sides[i].is.hanging.is_ghost[j])
217 typename ::parallel::distributed::
218 Triangulation<dim, spacedim>::cell_iterator cell =
221 *(sides[i].is.hanging.quad[j]));
223 for (k = 0; k < nsubs; ++k)
225 (*vertices_with_ghost_neighbors)
227 p8est_edge_corners[sides[i].edge][1 ^ j])]
228 .insert(subdomain_ids[k]);
241 template <
int dim,
int spacedim>
244 typename ::internal::p4est::iter<dim>::face_info *info,
248 int nsides = info->sides.elem_count;
249 auto *sides =
reinterpret_cast<
250 typename ::internal::p4est::iter<dim>::face_side *
>(
252 FindGhosts<dim, spacedim> *fg =
253 static_cast<FindGhosts<dim, spacedim> *
>(user_data);
254 sc_array_t *
subids = fg->subids;
255 const ::parallel::distributed::Triangulation<dim, spacedim>
259 std::map<unsigned int, std::set<::types::subdomain_id>>
261 int limit = (dim == 2) ? 2 : 4;
264 for (i = 0; i < nsides; ++i)
266 if (sides[i].is_hanging)
268 for (j = 0; j < limit; ++j)
270 if (sides[i].is.hanging.is_ghost[j])
272 typename ::parallel::distributed::
273 Triangulation<dim, spacedim>::cell_iterator cell =
276 *(sides[i].is.hanging.quad[j]));
280 *subid = cell->subdomain_id();
291 nsubs =
static_cast<int>(
subids->elem_count);
295 for (i = 0; i < nsides; ++i)
297 if (sides[i].is_hanging)
299 for (j = 0; j < limit; ++j)
301 if (!sides[i].is.hanging.is_ghost[j])
303 typename ::parallel::distributed::
304 Triangulation<dim, spacedim>::cell_iterator cell =
307 *(sides[i].is.hanging.quad[j]));
309 for (k = 0; k < nsubs; ++k)
313 (*vertices_with_ghost_neighbors)
315 p4est_face_corners[sides[i].face]
317 .insert(subdomain_ids[k]);
321 (*vertices_with_ghost_neighbors)
323 p8est_face_corners[sides[i].face]
325 .insert(subdomain_ids[k]);
339 p4est_quadrant_compare;
343 p4est_quadrant_childrenv;
347 p4est_quadrant_overlaps_tree;
352 p4est_quadrant_set_morton;
356 p4est_quadrant_is_equal;
360 p4est_quadrant_is_sibling;
364 p4est_quadrant_is_ancestor;
368 p4est_quadrant_ancestor_id;
374 p4est_comm_find_owner;
393 const int8_t *ctc) = p4est_connectivity_new_copy;
401 p4est_connectivity_join_faces;
404 p4est_connectivity_t *connectivity) = p4est_connectivity_destroy;
412 std::size_t data_size,
413 p4est_init_t init_fn,
414 void *user_pointer) = p4est_new_ext;
417 int copy_data) = p4est_copy;
422 int refine_recursive,
423 p4est_refine_t refine_fn,
424 p4est_init_t init_fn) = p4est_refine;
427 int coarsen_recursive,
428 p4est_coarsen_t coarsen_fn,
429 p4est_init_t init_fn) = p4est_coarsen;
433 p4est_init_t init_fn) = p4est_balance;
436 int partition_for_coarsening,
437 p4est_weight_t weight_fn) =
442 int save_data) = p4est_save;
445 const char *filename,
447 std::size_t data_size,
455 const char *filename,
462 const char *filename,
463 std::size_t *length) = p4est_connectivity_load;
470 const char *baseName) =
471 p4est_vtk_write_file;
481 std::size_t data_size,
482 p4est_init_t init_fn,
483 void *user_pointer) = p4est_reset_data;
498 const void *src_data,
499 std::size_t data_size) =
500 p4est_transfer_fixed;
508 const void *src_data,
509 std::size_t data_size) = p4est_transfer_fixed_begin;
512 p4est_transfer_fixed_end;
519 const int *dest_sizes,
520 const void *src_data,
521 const int *src_sizes) =
522 p4est_transfer_custom;
530 const int *dest_sizes,
531 const void *src_data,
532 const int *src_sizes) = p4est_transfer_custom_begin;
535 p4est_transfer_custom_end;
537# ifdef P4EST_SEARCH_LOCAL
543 sc_array_t *points) = p4est_search_partition;
551 double vxyz[3]) = p4est_qcoord_to_vertex;
554 p8est_quadrant_compare;
558 p8est_quadrant_childrenv;
562 p8est_quadrant_overlaps_tree;
567 p8est_quadrant_set_morton;
571 p8est_quadrant_is_equal;
575 p8est_quadrant_is_sibling;
579 p8est_quadrant_is_ancestor;
583 p8est_quadrant_ancestor_id;
589 p8est_comm_find_owner;
615 const int8_t *ctc) = p8est_connectivity_new_copy;
618 p8est_connectivity_t *connectivity) = p8est_connectivity_destroy;
626 p8est_connectivity_join_faces;
634 std::size_t data_size,
635 p8est_init_t init_fn,
636 void *user_pointer) = p8est_new_ext;
639 int copy_data) = p8est_copy;
644 int refine_recursive,
645 p8est_refine_t refine_fn,
646 p8est_init_t init_fn) = p8est_refine;
649 int coarsen_recursive,
650 p8est_coarsen_t coarsen_fn,
651 p8est_init_t init_fn) = p8est_coarsen;
655 p8est_init_t init_fn) = p8est_balance;
658 int partition_for_coarsening,
659 p8est_weight_t weight_fn) =
664 int save_data) = p8est_save;
667 const char *filename,
669 std::size_t data_size,
677 const char *filename,
684 const char *filename,
685 std::size_t *length) = p8est_connectivity_load;
692 const char *baseName) =
693 p8est_vtk_write_file;
703 std::size_t data_size,
704 p8est_init_t init_fn,
705 void *user_pointer) = p8est_reset_data;
720 const void *src_data,
721 std::size_t data_size) =
722 p8est_transfer_fixed;
730 const void *src_data,
731 std::size_t data_size) = p8est_transfer_fixed_begin;
734 p8est_transfer_fixed_end;
741 const int *dest_sizes,
742 const void *src_data,
743 const int *src_sizes) =
744 p8est_transfer_custom;
752 const int *dest_sizes,
753 const void *src_data,
754 const int *src_sizes) = p8est_transfer_custom_begin;
757 p8est_transfer_custom_end;
759# ifdef P4EST_SEARCH_LOCAL
765 sc_array_t *points) = p8est_search_partition;
774 double vxyz[3]) = p8est_qcoord_to_vertex;
783 for (
unsigned int c = 0;
784 c < ::GeometryInfo<dim>::max_children_per_cell;
789 P4EST_QUADRANT_INIT(&p4est_children[c]);
792 P8EST_QUADRANT_INIT(&p4est_children[c]);
809 P4EST_QUADRANT_INIT(&quad);
812 P8EST_QUADRANT_INIT(&quad);
845 Assert(coarse_grid_cell < parallel_forest->connectivity->num_trees,
847 return ((coarse_grid_cell >= parallel_forest->first_local_tree) &&
848 (coarse_grid_cell <= parallel_forest->last_local_tree));
860 connectivity->num_vertices,
861 connectivity->num_trees,
862 connectivity->num_corners,
863 connectivity->vertices,
864 connectivity->tree_to_vertex,
865 connectivity->tree_to_tree,
866 connectivity->tree_to_face,
867 connectivity->tree_to_corner,
868 connectivity->ctt_offset,
869 connectivity->corner_to_tree,
870 connectivity->corner_to_corner);
878 connectivity->num_vertices,
879 connectivity->num_trees,
880 connectivity->num_edges,
881 connectivity->num_corners,
882 connectivity->vertices,
883 connectivity->tree_to_vertex,
884 connectivity->tree_to_tree,
885 connectivity->tree_to_face,
886 connectivity->tree_to_edge,
887 connectivity->ett_offset,
888 connectivity->edge_to_tree,
889 connectivity->edge_to_edge,
890 connectivity->tree_to_corner,
891 connectivity->ctt_offset,
892 connectivity->corner_to_tree,
893 connectivity->corner_to_corner);
914 const int level_1 = (q1 << types<1>::max_n_child_indices_bits) >>
916 const int level_2 = (q2 << types<1>::max_n_child_indices_bits) >>
920 if (level_1 >= level_2)
930 return truncated_id_1 == truncated_id_2;
943 const int level_parent = (q << types<1>::max_n_child_indices_bits) >>
945 const int level_child = level_parent + 1;
948 p4est_children[0] = (q + 1);
970#include "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)
#define DEAL_II_NOT_IMPLEMENTED()
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)
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)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
std::map< unsigned int, std::set<::types::subdomain_id > > * vertices_with_ghost_neighbors