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;
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;
541 p8est_quadrant_compare;
545 p8est_quadrant_childrenv;
549 p8est_quadrant_overlaps_tree;
554 p8est_quadrant_set_morton;
558 p8est_quadrant_is_equal;
562 p8est_quadrant_is_sibling;
566 p8est_quadrant_is_ancestor;
570 p8est_quadrant_ancestor_id;
576 p8est_comm_find_owner;
602 const int8_t * ctc) = p8est_connectivity_new_copy;
605 p8est_connectivity_t *connectivity) = p8est_connectivity_destroy;
613 p8est_connectivity_join_faces;
621 std::size_t data_size,
622 p8est_init_t init_fn,
623 void * user_pointer) = p8est_new_ext;
626 int copy_data) = p8est_copy;
631 int refine_recursive,
632 p8est_refine_t refine_fn,
633 p8est_init_t init_fn) = p8est_refine;
636 int coarsen_recursive,
637 p8est_coarsen_t coarsen_fn,
638 p8est_init_t init_fn) = p8est_coarsen;
642 p8est_init_t init_fn) = p8est_balance;
645 int partition_for_coarsening,
646 p8est_weight_t weight_fn) =
651 int save_data) = p8est_save;
654 const char * filename,
656 std::size_t data_size,
664 const char * filename,
671 const char * filename,
672 std::size_t *length) = p8est_connectivity_load;
679 const char *baseName) =
680 p8est_vtk_write_file;
690 std::size_t data_size,
691 p8est_init_t init_fn,
692 void *user_pointer) = p8est_reset_data;
707 const void * src_data,
708 std::size_t data_size) =
709 p8est_transfer_fixed;
717 const void * src_data,
718 std::size_t data_size) = p8est_transfer_fixed_begin;
721 p8est_transfer_fixed_end;
728 const int * dest_sizes,
729 const void * src_data,
730 const int * src_sizes) =
731 p8est_transfer_custom;
739 const int * dest_sizes,
740 const void * src_data,
741 const int * src_sizes) = p8est_transfer_custom_begin;
744 p8est_transfer_custom_end;
755 for (
unsigned int c = 0;
756 c < ::GeometryInfo<dim>::max_children_per_cell;
761 P4EST_QUADRANT_INIT(&p4est_children[c]);
764 P8EST_QUADRANT_INIT(&p4est_children[c]);
781 P4EST_QUADRANT_INIT(&quad);
784 P8EST_QUADRANT_INIT(&quad);
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));
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);
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);
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) >>
891 if (level_1 >= level_2)
901 return truncated_id_1 == truncated_id_2;
914 const int level_parent = (q << types<1>::max_n_child_indices_bits) >>
916 const int level_child = level_parent + 1;
919 p4est_children[0] = (q + 1);
940#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)
TriangulationBase< dim, spacedim > Triangulation
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
std::map< unsigned int, std::set<::types::subdomain_id > > * vertices_with_ghost_neighbors