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>
67 const typename ::parallel::distributed::Triangulation<dim,
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;
394 const int8_t * ctc) = p4est_connectivity_new_copy;
402 p4est_connectivity_join_faces;
405 p4est_connectivity_t *connectivity) = p4est_connectivity_destroy;
413 std::size_t data_size,
414 p4est_init_t init_fn,
415 void * user_pointer) = p4est_new_ext;
418 int copy_data) = p4est_copy;
423 int refine_recursive,
424 p4est_refine_t refine_fn,
425 p4est_init_t init_fn) = p4est_refine;
428 int coarsen_recursive,
429 p4est_coarsen_t coarsen_fn,
430 p4est_init_t init_fn) = p4est_coarsen;
434 p4est_init_t init_fn) = p4est_balance;
437 int partition_for_coarsening,
438 p4est_weight_t weight_fn) =
443 int save_data) = p4est_save;
446 const char * filename,
448 std::size_t data_size,
456 const char * filename,
463 const char * filename,
464 std::size_t *length) = p4est_connectivity_load;
471 const char *baseName) =
472 p4est_vtk_write_file;
482 std::size_t data_size,
483 p4est_init_t init_fn,
484 void *user_pointer) = p4est_reset_data;
499 const void * src_data,
500 std::size_t data_size) =
501 p4est_transfer_fixed;
509 const void * src_data,
510 std::size_t data_size) = p4est_transfer_fixed_begin;
513 p4est_transfer_fixed_end;
520 const int * dest_sizes,
521 const void * src_data,
522 const int * src_sizes) =
523 p4est_transfer_custom;
531 const int * dest_sizes,
532 const void * src_data,
533 const int * src_sizes) = p4est_transfer_custom_begin;
536 p4est_transfer_custom_end;
538# ifdef P4EST_SEARCH_LOCAL
544 sc_array_t * points) = p4est_search_partition;
552 double vxyz[3]) = p4est_qcoord_to_vertex;
555 p8est_quadrant_compare;
559 p8est_quadrant_childrenv;
563 p8est_quadrant_overlaps_tree;
568 p8est_quadrant_set_morton;
572 p8est_quadrant_is_equal;
576 p8est_quadrant_is_sibling;
580 p8est_quadrant_is_ancestor;
584 p8est_quadrant_ancestor_id;
590 p8est_comm_find_owner;
616 const int8_t * ctc) = p8est_connectivity_new_copy;
619 p8est_connectivity_t *connectivity) = p8est_connectivity_destroy;
627 p8est_connectivity_join_faces;
635 std::size_t data_size,
636 p8est_init_t init_fn,
637 void * user_pointer) = p8est_new_ext;
640 int copy_data) = p8est_copy;
645 int refine_recursive,
646 p8est_refine_t refine_fn,
647 p8est_init_t init_fn) = p8est_refine;
650 int coarsen_recursive,
651 p8est_coarsen_t coarsen_fn,
652 p8est_init_t init_fn) = p8est_coarsen;
656 p8est_init_t init_fn) = p8est_balance;
659 int partition_for_coarsening,
660 p8est_weight_t weight_fn) =
665 int save_data) = p8est_save;
668 const char * filename,
670 std::size_t data_size,
678 const char * filename,
685 const char * filename,
686 std::size_t *length) = p8est_connectivity_load;
693 const char *baseName) =
694 p8est_vtk_write_file;
704 std::size_t data_size,
705 p8est_init_t init_fn,
706 void *user_pointer) = p8est_reset_data;
721 const void * src_data,
722 std::size_t data_size) =
723 p8est_transfer_fixed;
731 const void * src_data,
732 std::size_t data_size) = p8est_transfer_fixed_begin;
735 p8est_transfer_fixed_end;
742 const int * dest_sizes,
743 const void * src_data,
744 const int * src_sizes) =
745 p8est_transfer_custom;
753 const int * dest_sizes,
754 const void * src_data,
755 const int * src_sizes) = p8est_transfer_custom_begin;
758 p8est_transfer_custom_end;
760# ifdef P4EST_SEARCH_LOCAL
766 sc_array_t * points) = p8est_search_partition;
775 double vxyz[3]) = p8est_qcoord_to_vertex;
784 for (
unsigned int c = 0;
785 c < ::GeometryInfo<dim>::max_children_per_cell;
790 P4EST_QUADRANT_INIT(&p4est_children[c]);
793 P8EST_QUADRANT_INIT(&p4est_children[c]);
810 P4EST_QUADRANT_INIT(&quad);
813 P8EST_QUADRANT_INIT(&quad);
846 Assert(coarse_grid_cell < parallel_forest->connectivity->num_trees,
848 return ((coarse_grid_cell >= parallel_forest->first_local_tree) &&
849 (coarse_grid_cell <= parallel_forest->last_local_tree));
861 connectivity->num_vertices,
862 connectivity->num_trees,
863 connectivity->num_corners,
864 connectivity->vertices,
865 connectivity->tree_to_vertex,
866 connectivity->tree_to_tree,
867 connectivity->tree_to_face,
868 connectivity->tree_to_corner,
869 connectivity->ctt_offset,
870 connectivity->corner_to_tree,
871 connectivity->corner_to_corner);
879 connectivity->num_vertices,
880 connectivity->num_trees,
881 connectivity->num_edges,
882 connectivity->num_corners,
883 connectivity->vertices,
884 connectivity->tree_to_vertex,
885 connectivity->tree_to_tree,
886 connectivity->tree_to_face,
887 connectivity->tree_to_edge,
888 connectivity->ett_offset,
889 connectivity->edge_to_tree,
890 connectivity->edge_to_edge,
891 connectivity->tree_to_corner,
892 connectivity->ctt_offset,
893 connectivity->corner_to_tree,
894 connectivity->corner_to_corner);
915 const int level_1 = (q1 << types<1>::max_n_child_indices_bits) >>
917 const int level_2 = (q2 << types<1>::max_n_child_indices_bits) >>
921 if (level_1 >= level_2)
931 return truncated_id_1 == truncated_id_2;
944 const int level_parent = (q << types<1>::max_n_child_indices_bits) >>
946 const int level_child = level_parent + 1;
949 p4est_children[0] = (q + 1);
971#include "p4est_wrappers.inst"
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
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)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
std::map< unsigned int, std::set<::types::subdomain_id > > * vertices_with_ghost_neighbors