Reference documentation for deal.II version GIT d6cf33b98c 2023-09-22 19:50:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
p4est_wrappers.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2023 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_p4est_wrappers_h
17 #define dealii_p4est_wrappers_h
18 
19 #include <deal.II/base/config.h>
20 
22 
23 #ifdef DEAL_II_WITH_P4EST
24 # include <p4est_bits.h>
25 # include <p4est_communication.h>
26 # include <p4est_extended.h>
27 # include <p4est_ghost.h>
28 # include <p4est_iterate.h>
29 # include <p4est_search.h>
30 # include <p4est_vtk.h>
31 # include <p8est_bits.h>
32 # include <p8est_communication.h>
33 # include <p8est_extended.h>
34 # include <p8est_ghost.h>
35 # include <p8est_iterate.h>
36 # include <p8est_search.h>
37 # include <p8est_vtk.h>
38 
39 # include <limits>
40 
42 
43 // Forward declaration
44 # ifndef DOXYGEN
45 namespace parallel
46 {
47  namespace distributed
48  {
49  template <int dim, int spacedim>
50  DEAL_II_CXX20_REQUIRES((concepts::is_valid_dim_spacedim<dim, spacedim>))
51  class Triangulation;
52  }
53 } // namespace parallel
54 # endif
55 
56 namespace internal
57 {
58  namespace p4est
59  {
67  template <int>
68  struct types;
69 
70  // these struct mimics p4est for 1d
71  template <>
72  struct types<1>
73  {
74  // id of a quadrant is an integeger
75  using quadrant = int;
76 
77  // maximum number of children
78  static const int max_n_child_indices_bits = 27;
79 
80  // number of bits the data type of id has
81  static const int n_bits = std::numeric_limits<quadrant>::digits;
82  };
83 
84  template <>
85  struct types<2>
86  {
87  using connectivity = p4est_connectivity_t;
88  using forest = p4est_t;
89  using tree = p4est_tree_t;
90  using quadrant = p4est_quadrant_t;
91  using quadrant_coord = p4est_qcoord_t;
92  using topidx = p4est_topidx_t;
93  using locidx = p4est_locidx_t;
94  using gloidx = p4est_gloidx_t;
95  using balance_type = p4est_connect_type_t;
96  using ghost = p4est_ghost_t;
97  using transfer_context = p4est_transfer_context_t;
98 # ifdef P4EST_SEARCH_LOCAL
99  using search_partition_callback = p4est_search_partition_t;
100 # endif
101  };
102 
103  template <>
104  struct types<3>
105  {
106  using connectivity = p8est_connectivity_t;
107  using forest = p8est_t;
108  using tree = p8est_tree_t;
109  using quadrant = p8est_quadrant_t;
110  using quadrant_coord = p4est_qcoord_t;
111  using topidx = p4est_topidx_t;
112  using locidx = p4est_locidx_t;
113  using gloidx = p4est_gloidx_t;
114  using balance_type = p8est_connect_type_t;
115  using ghost = p8est_ghost_t;
116  using transfer_context = p8est_transfer_context_t;
117 # ifdef P4EST_SEARCH_LOCAL
118  using search_partition_callback = p8est_search_partition_t;
119 # endif
120  };
121 
122 
123 
131  template <int dim>
132  struct functions;
133 
134  template <>
135  struct functions<2>
136  {
137  static int (&quadrant_compare)(const void *v1, const void *v2);
138 
139  static void (&quadrant_childrenv)(const types<2>::quadrant *q,
140  types<2>::quadrant c[]);
141 
142  static int (&quadrant_overlaps_tree)(types<2>::tree *tree,
143  const types<2>::quadrant *q);
144 
145  static void (&quadrant_set_morton)(types<2>::quadrant *quadrant,
146  int level,
147  std::uint64_t id);
148 
149  static int (&quadrant_is_equal)(const types<2>::quadrant *q1,
150  const types<2>::quadrant *q2);
151 
152  static int (&quadrant_is_sibling)(const types<2>::quadrant *q1,
153  const types<2>::quadrant *q2);
154 
155  static int (&quadrant_is_ancestor)(const types<2>::quadrant *q1,
156  const types<2>::quadrant *q2);
157 
158  static int (&quadrant_ancestor_id)(const types<2>::quadrant *q,
159  int level);
160 
161  static int (&comm_find_owner)(types<2>::forest *p4est,
162  const types<2>::locidx which_tree,
163  const types<2>::quadrant *q,
164  const int guess);
165 
167  types<2>::topidx num_vertices,
168  types<2>::topidx num_trees,
169  types<2>::topidx num_corners,
170  types<2>::topidx num_vtt);
171 
173  types<2>::topidx num_vertices,
174  types<2>::topidx num_trees,
175  types<2>::topidx num_corners,
176  const double *vertices,
177  const types<2>::topidx *ttv,
178  const types<2>::topidx *ttt,
179  const int8_t *ttf,
180  const types<2>::topidx *ttc,
181  const types<2>::topidx *coff,
182  const types<2>::topidx *ctt,
183  const int8_t *ctc);
184 
185  static void (&connectivity_join_faces)(types<2>::connectivity *conn,
186  types<2>::topidx tree_left,
187  types<2>::topidx tree_right,
188  int face_left,
189  int face_right,
190  int orientation);
191 
192 
193 
194  static void (&connectivity_destroy)(p4est_connectivity_t *connectivity);
195 
197  MPI_Comm mpicomm,
198  types<2>::connectivity *connectivity,
199  types<2>::locidx min_quadrants,
200  int min_level,
201  int fill_uniform,
202  std::size_t data_size,
203  p4est_init_t init_fn,
204  void *user_pointer);
205 
207  int copy_data);
208 
209  static void (&destroy)(types<2>::forest *p4est);
210 
211  static void (&refine)(types<2>::forest *p4est,
212  int refine_recursive,
213  p4est_refine_t refine_fn,
214  p4est_init_t init_fn);
215 
216  static void (&coarsen)(types<2>::forest *p4est,
217  int coarsen_recursive,
218  p4est_coarsen_t coarsen_fn,
219  p4est_init_t init_fn);
220 
221  static void (&balance)(types<2>::forest *p4est,
223  p4est_init_t init_fn);
224 
226  int partition_for_coarsening,
227  p4est_weight_t weight_fn);
228 
229  static void (&save)(const char *filename,
230  types<2>::forest *p4est,
231  int save_data);
232 
233  static types<2>::forest *(&load_ext)(const char *filename,
234  MPI_Comm mpicomm,
235  std::size_t data_size,
236  int load_data,
237  int autopartition,
238  int broadcasthead,
239  void *user_pointer,
240  types<2>::connectivity **p4est);
241 
242  static int (&connectivity_save)(const char *filename,
243  types<2>::connectivity *connectivity);
244 
245  static int (&connectivity_is_valid)(types<2>::connectivity *connectivity);
246 
247  static types<2>::connectivity *(&connectivity_load)(const char *filename,
248  std::size_t *length);
249 
250  static unsigned int (&checksum)(types<2>::forest *p4est);
251 
252  static void (&vtk_write_file)(types<2>::forest *p4est,
253  p4est_geometry_t *,
254  const char *baseName);
255 
257  types<2>::balance_type btype);
258 
259  static void (&ghost_destroy)(types<2>::ghost *ghost);
260 
261  static void (&reset_data)(types<2>::forest *p4est,
262  std::size_t data_size,
263  p4est_init_t init_fn,
264  void *user_pointer);
265 
266  static std::size_t (&forest_memory_used)(types<2>::forest *p4est);
267 
268  static std::size_t (&connectivity_memory_used)(
269  types<2>::connectivity *p4est);
270 
271  template <int spacedim>
272  static void
274  ::internal::p4est::types<2>::ghost *parallel_ghost,
275  void *user_data);
276 
277  static constexpr unsigned int max_level = P4EST_MAXLEVEL;
278 
279  static void (&transfer_fixed)(const types<2>::gloidx *dest_gfq,
280  const types<2>::gloidx *src_gfq,
281  MPI_Comm mpicomm,
282  int tag,
283  void *dest_data,
284  const void *src_data,
285  std::size_t data_size);
286 
288  const types<2>::gloidx *dest_gfq,
289  const types<2>::gloidx *src_gfq,
290  MPI_Comm mpicomm,
291  int tag,
292  void *dest_data,
293  const void *src_data,
294  std::size_t data_size);
295 
296  static void (&transfer_fixed_end)(types<2>::transfer_context *tc);
297 
298  static void (&transfer_custom)(const types<2>::gloidx *dest_gfq,
299  const types<2>::gloidx *src_gfq,
300  MPI_Comm mpicomm,
301  int tag,
302  void *dest_data,
303  const int *dest_sizes,
304  const void *src_data,
305  const int *src_sizes);
306 
308  const types<2>::gloidx *dest_gfq,
309  const types<2>::gloidx *src_gfq,
310  MPI_Comm mpicomm,
311  int tag,
312  void *dest_data,
313  const int *dest_sizes,
314  const void *src_data,
315  const int *src_sizes);
316 
317  static void (&transfer_custom_end)(types<2>::transfer_context *tc);
318 
319 # ifdef P4EST_SEARCH_LOCAL
320  static void (&search_partition)(
321  types<2>::forest *forest,
322  int call_post,
325  sc_array_t *points);
326 # endif
327 
328  static void (&quadrant_coord_to_vertex)(
329  types<2>::connectivity *connectivity,
330  types<2>::topidx treeid,
333  double vxyz[3]);
334  };
335 
336 
337  template <>
338  struct functions<3>
339  {
340  static int (&quadrant_compare)(const void *v1, const void *v2);
341 
342  static void (&quadrant_childrenv)(const types<3>::quadrant *q,
343  types<3>::quadrant c[]);
344 
345  static int (&quadrant_overlaps_tree)(types<3>::tree *tree,
346  const types<3>::quadrant *q);
347 
348  static void (&quadrant_set_morton)(types<3>::quadrant *quadrant,
349  int level,
350  std::uint64_t id);
351 
352  static int (&quadrant_is_equal)(const types<3>::quadrant *q1,
353  const types<3>::quadrant *q2);
354 
355  static int (&quadrant_is_sibling)(const types<3>::quadrant *q1,
356  const types<3>::quadrant *q2);
357 
358  static int (&quadrant_is_ancestor)(const types<3>::quadrant *q1,
359  const types<3>::quadrant *q2);
360  static int (&quadrant_ancestor_id)(const types<3>::quadrant *q,
361  int level);
362 
363  static int (&comm_find_owner)(types<3>::forest *p4est,
364  const types<3>::locidx which_tree,
365  const types<3>::quadrant *q,
366  const int guess);
367 
369  types<3>::topidx num_vertices,
370  types<3>::topidx num_trees,
371  types<3>::topidx num_edges,
372  types<3>::topidx num_ett,
373  types<3>::topidx num_corners,
374  types<3>::topidx num_ctt);
375 
377  types<3>::topidx num_vertices,
378  types<3>::topidx num_trees,
379  types<3>::topidx num_edges,
380  types<3>::topidx num_corners,
381  const double *vertices,
382  const types<3>::topidx *ttv,
383  const types<3>::topidx *ttt,
384  const int8_t *ttf,
385  const types<3>::topidx *tte,
386  const types<3>::topidx *eoff,
387  const types<3>::topidx *ett,
388  const int8_t *ete,
389  const types<3>::topidx *ttc,
390  const types<3>::topidx *coff,
391  const types<3>::topidx *ctt,
392  const int8_t *ctc);
393 
394  static void (&connectivity_join_faces)(types<3>::connectivity *conn,
395  types<3>::topidx tree_left,
396  types<3>::topidx tree_right,
397  int face_left,
398  int face_right,
399  int orientation);
400 
401  static void (&connectivity_destroy)(p8est_connectivity_t *connectivity);
402 
404  MPI_Comm mpicomm,
405  types<3>::connectivity *connectivity,
406  types<3>::locidx min_quadrants,
407  int min_level,
408  int fill_uniform,
409  std::size_t data_size,
410  p8est_init_t init_fn,
411  void *user_pointer);
412 
414  int copy_data);
415 
416  static void (&destroy)(types<3>::forest *p8est);
417 
418  static void (&refine)(types<3>::forest *p8est,
419  int refine_recursive,
420  p8est_refine_t refine_fn,
421  p8est_init_t init_fn);
422 
423  static void (&coarsen)(types<3>::forest *p8est,
424  int coarsen_recursive,
425  p8est_coarsen_t coarsen_fn,
426  p8est_init_t init_fn);
427 
428  static void (&balance)(types<3>::forest *p8est,
430  p8est_init_t init_fn);
431 
433  int partition_for_coarsening,
434  p8est_weight_t weight_fn);
435 
436  static void (&save)(const char *filename,
437  types<3>::forest *p4est,
438  int save_data);
439 
440  static types<3>::forest *(&load_ext)(const char *filename,
441  MPI_Comm mpicomm,
442  std::size_t data_size,
443  int load_data,
444  int autopartition,
445  int broadcasthead,
446  void *user_pointer,
447  types<3>::connectivity **p4est);
448 
449  static int (&connectivity_save)(const char *filename,
450  types<3>::connectivity *connectivity);
451 
452  static int (&connectivity_is_valid)(types<3>::connectivity *connectivity);
453 
454  static types<3>::connectivity *(&connectivity_load)(const char *filename,
455  std::size_t *length);
456 
457  static unsigned int (&checksum)(types<3>::forest *p8est);
458 
459  static void (&vtk_write_file)(types<3>::forest *p8est,
460  p8est_geometry_t *,
461  const char *baseName);
463  types<3>::balance_type btype);
464 
465  static void (&ghost_destroy)(types<3>::ghost *ghost);
466 
467  static void (&reset_data)(types<3>::forest *p4est,
468  std::size_t data_size,
469  p8est_init_t init_fn,
470  void *user_pointer);
471 
472  static std::size_t (&forest_memory_used)(types<3>::forest *p4est);
473 
474  static std::size_t (&connectivity_memory_used)(
475  types<3>::connectivity *p4est);
476 
477  static constexpr unsigned int max_level = P8EST_MAXLEVEL;
478 
479  static void (&transfer_fixed)(const types<3>::gloidx *dest_gfq,
480  const types<3>::gloidx *src_gfq,
481  MPI_Comm mpicomm,
482  int tag,
483  void *dest_data,
484  const void *src_data,
485  std::size_t data_size);
486 
488  const types<3>::gloidx *dest_gfq,
489  const types<3>::gloidx *src_gfq,
490  MPI_Comm mpicomm,
491  int tag,
492  void *dest_data,
493  const void *src_data,
494  std::size_t data_size);
495 
496  static void (&transfer_fixed_end)(types<3>::transfer_context *tc);
497 
498  static void (&transfer_custom)(const types<3>::gloidx *dest_gfq,
499  const types<3>::gloidx *src_gfq,
500  MPI_Comm mpicomm,
501  int tag,
502  void *dest_data,
503  const int *dest_sizes,
504  const void *src_data,
505  const int *src_sizes);
506 
508  const types<3>::gloidx *dest_gfq,
509  const types<3>::gloidx *src_gfq,
510  MPI_Comm mpicomm,
511  int tag,
512  void *dest_data,
513  const int *dest_sizes,
514  const void *src_data,
515  const int *src_sizes);
516 
517  static void (&transfer_custom_end)(types<3>::transfer_context *tc);
518 
519 # ifdef P4EST_SEARCH_LOCAL
520  static void (&search_partition)(
521  types<3>::forest *forest,
522  int call_post,
525  sc_array_t *points);
526 # endif
527 
528  static void (&quadrant_coord_to_vertex)(
529  types<3>::connectivity *connectivity,
530  types<3>::topidx treeid,
534  double vxyz[3]);
535  };
536 
537 
538 
545  template <int dim>
546  struct iter;
547 
548  template <>
549  struct iter<2>
550  {
551  using corner_info = p4est_iter_corner_info_t;
552  using corner_side = p4est_iter_corner_side_t;
553  using corner_iter = p4est_iter_corner_t;
554  using face_info = p4est_iter_face_info_t;
555  using face_side = p4est_iter_face_side_t;
556  using face_iter = p4est_iter_face_t;
557  };
558 
559  template <>
560  struct iter<3>
561  {
562  using corner_info = p8est_iter_corner_info_t;
563  using corner_side = p8est_iter_corner_side_t;
564  using corner_iter = p8est_iter_corner_t;
565  using edge_info = p8est_iter_edge_info_t;
566  using edge_side = p8est_iter_edge_side_t;
567  using edge_iter = p8est_iter_edge_t;
568  using face_info = p8est_iter_face_info_t;
569  using face_side = p8est_iter_face_side_t;
570  using face_iter = p8est_iter_face_t;
571  };
572 
573 
574 
579  template <int dim>
580  void
582  const typename types<dim>::quadrant &p4est_cell,
583  typename types<dim>::quadrant (
584  &p4est_children)[::GeometryInfo<dim>::max_children_per_cell]);
585 
586 
587 
591  template <int dim>
592  void
594 
595 
596 
600  template <int dim>
601  bool
602  quadrant_is_equal(const typename types<dim>::quadrant &q1,
603  const typename types<dim>::quadrant &q2);
604 
605 
606 
610  template <int dim>
611  bool
612  quadrant_is_ancestor(const typename types<dim>::quadrant &q1,
613  const typename types<dim>::quadrant &q2);
614 
615 
616 
620  template <int dim>
621  bool
622  tree_exists_locally(const typename types<dim>::forest *parallel_forest,
623  const typename types<dim>::topidx coarse_grid_cell);
624 
625 
629  template <int dim>
630  typename types<dim>::connectivity *
631  copy_connectivity(const typename types<dim>::connectivity *connectivity);
632 
633 # ifndef DOXYGEN
634  template <>
635  typename types<2>::connectivity *
636  copy_connectivity<2>(const typename types<2>::connectivity *connectivity);
637 
638  template <>
639  typename types<3>::connectivity *
640  copy_connectivity<3>(const typename types<3>::connectivity *connectivity);
641 # endif
642  } // namespace p4est
643 } // namespace internal
644 
646 
647 #endif // DEAL_II_WITH_P4EST
648 
649 #endif // dealii_p4est_wrappers_h
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:477
#define DEAL_II_CXX20_REQUIRES(condition)
Definition: config.h:166
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:478
Point< 3 > vertices[4]
unsigned int level
Definition: grid_out.cc:4617
const unsigned int v1
Definition: grid_tools.cc:1063
void refine(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
void coarsen(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold)
void partition(const SparsityPattern &sparsity_pattern, const unsigned int n_partitions, std::vector< unsigned int > &partition_indices, const Partitioner partitioner=Partitioner::metis)
types< 3 >::connectivity * copy_connectivity< 3 >(const typename types< 3 >::connectivity *connectivity)
types< 2 >::connectivity * copy_connectivity< 2 >(const typename types< 2 >::connectivity *connectivity)
void init_coarse_quadrant(typename types< dim >::quadrant &quad)
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)
types< dim >::connectivity * copy_connectivity(const typename types< dim >::connectivity *connectivity)
bool tree_exists_locally(const typename types< dim >::forest *parallel_forest, const typename types< dim >::topidx coarse_grid_cell)
Definition: types.h:33
static types< 2 >::transfer_context *(&) transfer_custom_begin(const types< 2 >::gloidx *dest_gfq, const types< 2 >::gloidx *src_gfq, MPI_Comm mpicomm, int tag, void *dest_data, const int *dest_sizes, const void *src_data, const int *src_sizes)
static types< 2 >::ghost *(&) ghost_new(types< 2 >::forest *p4est, types< 2 >::balance_type btype)
static types< 2 >::forest *(&) load_ext(const char *filename, MPI_Comm mpicomm, std::size_t data_size, int load_data, int autopartition, int broadcasthead, void *user_pointer, types< 2 >::connectivity **p4est)
static types< 2 >::connectivity *(&) connectivity_load(const char *filename, std::size_t *length)
static types< 2 >::forest *(&) new_forest(MPI_Comm mpicomm, types< 2 >::connectivity *connectivity, types< 2 >::locidx min_quadrants, int min_level, int fill_uniform, std::size_t data_size, p4est_init_t init_fn, void *user_pointer)
static types< 2 >::connectivity *(&) connectivity_new(types< 2 >::topidx num_vertices, types< 2 >::topidx num_trees, types< 2 >::topidx num_corners, types< 2 >::topidx num_vtt)
static void iterate(::internal::p4est::types< 2 >::forest *parallel_forest, ::internal::p4est::types< 2 >::ghost *parallel_ghost, void *user_data)
static types< 2 >::connectivity *(&) connectivity_new_copy(types< 2 >::topidx num_vertices, types< 2 >::topidx num_trees, types< 2 >::topidx num_corners, const double *vertices, const types< 2 >::topidx *ttv, const types< 2 >::topidx *ttt, const int8_t *ttf, const types< 2 >::topidx *ttc, const types< 2 >::topidx *coff, const types< 2 >::topidx *ctt, const int8_t *ctc)
static types< 2 >::transfer_context *(&) transfer_fixed_begin(const types< 2 >::gloidx *dest_gfq, const types< 2 >::gloidx *src_gfq, MPI_Comm mpicomm, int tag, void *dest_data, const void *src_data, std::size_t data_size)
static types< 2 >::forest *(&) copy_forest(types< 2 >::forest *input, int copy_data)
static types< 3 >::forest *(&) copy_forest(types< 3 >::forest *input, int copy_data)
static types< 3 >::transfer_context *(&) transfer_custom_begin(const types< 3 >::gloidx *dest_gfq, const types< 3 >::gloidx *src_gfq, MPI_Comm mpicomm, int tag, void *dest_data, const int *dest_sizes, const void *src_data, const int *src_sizes)
static types< 3 >::transfer_context *(&) transfer_fixed_begin(const types< 3 >::gloidx *dest_gfq, const types< 3 >::gloidx *src_gfq, MPI_Comm mpicomm, int tag, void *dest_data, const void *src_data, std::size_t data_size)
static types< 3 >::ghost *(&) ghost_new(types< 3 >::forest *p4est, types< 3 >::balance_type btype)
static types< 3 >::connectivity *(&) connectivity_load(const char *filename, std::size_t *length)
static types< 3 >::connectivity *(&) connectivity_new_copy(types< 3 >::topidx num_vertices, types< 3 >::topidx num_trees, types< 3 >::topidx num_edges, types< 3 >::topidx num_corners, const double *vertices, const types< 3 >::topidx *ttv, const types< 3 >::topidx *ttt, const int8_t *ttf, const types< 3 >::topidx *tte, const types< 3 >::topidx *eoff, const types< 3 >::topidx *ett, const int8_t *ete, const types< 3 >::topidx *ttc, const types< 3 >::topidx *coff, const types< 3 >::topidx *ctt, const int8_t *ctc)
static types< 3 >::forest *(&) load_ext(const char *filename, MPI_Comm mpicomm, std::size_t data_size, int load_data, int autopartition, int broadcasthead, void *user_pointer, types< 3 >::connectivity **p4est)
static types< 3 >::forest *(&) new_forest(MPI_Comm mpicomm, types< 3 >::connectivity *connectivity, types< 3 >::locidx min_quadrants, int min_level, int fill_uniform, std::size_t data_size, p8est_init_t init_fn, void *user_pointer)
static types< 3 >::connectivity *(&) connectivity_new(types< 3 >::topidx num_vertices, types< 3 >::topidx num_trees, types< 3 >::topidx num_edges, types< 3 >::topidx num_ett, types< 3 >::topidx num_corners, types< 3 >::topidx num_ctt)
p4est_iter_corner_t corner_iter
p4est_iter_corner_info_t corner_info
p4est_iter_face_info_t face_info
p4est_iter_face_side_t face_side
p4est_iter_corner_side_t corner_side
p4est_iter_face_t face_iter
p8est_iter_corner_t corner_iter
p8est_iter_corner_info_t corner_info
p8est_iter_face_info_t face_info
p8est_iter_edge_side_t edge_side
p8est_iter_edge_t edge_iter
p8est_iter_edge_info_t edge_info
p8est_iter_face_side_t face_side
p8est_iter_face_t face_iter
p8est_iter_corner_side_t corner_side
p4est_connectivity_t connectivity
p4est_connect_type_t balance_type
p4est_transfer_context_t transfer_context
p8est_connectivity_t connectivity
p8est_transfer_context_t transfer_context
p8est_connect_type_t balance_type