Reference documentation for deal.II version 9.1.0-pre
p4est_wrappers.cc
1 #include <deal.II/distributed/p4est_wrappers.h>
2 #include <deal.II/distributed/tria.h>
3 
4 DEAL_II_NAMESPACE_OPEN
5 
6 #ifdef DEAL_II_WITH_P4EST
7 
8 namespace internal
9 {
10  namespace p4est
11  {
12  namespace
13  {
14  template <int dim, int spacedim>
15  typename ::Triangulation<dim, spacedim>::cell_iterator
16  cell_from_quad(
17  const ::parallel::distributed::Triangulation<dim, spacedim>
18  *triangulation,
19  const typename ::internal::p4est::types<dim>::topidx treeidx,
20  const typename ::internal::p4est::types<dim>::quadrant &quad)
21  {
22  int i, l = quad.level;
23  ::types::global_dof_index dealii_index =
24  triangulation->get_p4est_tree_to_coarse_cell_permutation()[treeidx];
25 
26  for (i = 0; i < l; i++)
27  {
28  typename ::Triangulation<dim, spacedim>::cell_iterator cell(
29  triangulation, i, dealii_index);
30  const int child_id =
32  &quad, i + 1);
33  Assert(cell->has_children(),
34  ExcMessage("p4est quadrant does not correspond to a cell!"));
35  dealii_index = cell->child_index(child_id);
36  }
37 
38  typename ::Triangulation<dim, spacedim>::cell_iterator out_cell(
39  triangulation, l, dealii_index);
40 
41  return out_cell;
42  }
43 
48  template <int dim, int spacedim>
49  struct FindGhosts
50  {
51  const typename ::parallel::distributed::Triangulation<dim,
52  spacedim>
53  * triangulation;
54  sc_array_t *subids;
55  std::map<unsigned int, std::set<::types::subdomain_id>>
56  *vertices_with_ghost_neighbors;
57  };
58 
59 
65  template <int dim, int spacedim>
66  void
67  find_ghosts_corner(
68  typename ::internal::p4est::iter<dim>::corner_info *info,
69  void * user_data)
70  {
71  int i, j;
72  int nsides = info->sides.elem_count;
73  typename ::internal::p4est::iter<dim>::corner_side *sides =
74  (typename ::internal::p4est::iter<dim>::corner_side
75  *)(info->sides.array);
76  FindGhosts<dim, spacedim> *fg =
77  static_cast<FindGhosts<dim, spacedim> *>(user_data);
78  sc_array_t *subids = fg->subids;
79  const ::parallel::distributed::Triangulation<dim, spacedim>
80  * triangulation = fg->triangulation;
81  int nsubs;
82  ::types::subdomain_id *subdomain_ids;
83  std::map<unsigned int, std::set<::types::subdomain_id>>
84  *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
85 
86  subids->elem_count = 0;
87  for (i = 0; i < nsides; i++)
88  {
89  if (sides[i].is_ghost)
90  {
91  typename ::parallel::distributed::
92  Triangulation<dim, spacedim>::cell_iterator cell =
93  cell_from_quad(triangulation,
94  sides[i].treeid,
95  *(sides[i].quad));
96  Assert(cell->is_ghost(),
97  ExcMessage("ghost quad did not find ghost cell"));
98  ::types::subdomain_id *subid =
99  static_cast<::types::subdomain_id *>(
100  sc_array_push(subids));
101  *subid = cell->subdomain_id();
102  }
103  }
104 
105  if (!subids->elem_count)
106  {
107  return;
108  }
109 
110  nsubs = (int)subids->elem_count;
111  subdomain_ids = (::types::subdomain_id *)(subids->array);
112 
113  for (i = 0; i < nsides; i++)
114  {
115  if (!sides[i].is_ghost)
116  {
117  typename ::parallel::distributed::
118  Triangulation<dim, spacedim>::cell_iterator cell =
119  cell_from_quad(triangulation,
120  sides[i].treeid,
121  *(sides[i].quad));
122 
123  Assert(!cell->is_ghost(),
124  ExcMessage("local quad found ghost cell"));
125 
126  for (j = 0; j < nsubs; j++)
127  {
128  (*vertices_with_ghost_neighbors)[cell->vertex_index(
129  sides[i].corner)]
130  .insert(subdomain_ids[j]);
131  }
132  }
133  }
134 
135  subids->elem_count = 0;
136  }
137 
141  template <int dim, int spacedim>
142  void
143  find_ghosts_edge(
144  typename ::internal::p4est::iter<dim>::edge_info *info,
145  void * user_data)
146  {
147  int i, j, k;
148  int nsides = info->sides.elem_count;
149  typename ::internal::p4est::iter<dim>::edge_side *sides =
150  (typename ::internal::p4est::iter<dim>::edge_side *)(info->sides
151  .array);
152  FindGhosts<dim, spacedim> *fg =
153  static_cast<FindGhosts<dim, spacedim> *>(user_data);
154  sc_array_t *subids = fg->subids;
155  const ::parallel::distributed::Triangulation<dim, spacedim>
156  * triangulation = fg->triangulation;
157  int nsubs;
158  ::types::subdomain_id *subdomain_ids;
159  std::map<unsigned int, std::set<::types::subdomain_id>>
160  *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
161 
162  subids->elem_count = 0;
163  for (i = 0; i < nsides; i++)
164  {
165  if (sides[i].is_hanging)
166  {
167  for (j = 0; j < 2; j++)
168  {
169  if (sides[i].is.hanging.is_ghost[j])
170  {
171  typename ::parallel::distributed::
172  Triangulation<dim, spacedim>::cell_iterator cell =
173  cell_from_quad(triangulation,
174  sides[i].treeid,
175  *(sides[i].is.hanging.quad[j]));
176  ::types::subdomain_id *subid =
177  static_cast<::types::subdomain_id *>(
178  sc_array_push(subids));
179  *subid = cell->subdomain_id();
180  }
181  }
182  }
183  }
184 
185  if (!subids->elem_count)
186  {
187  return;
188  }
189 
190  nsubs = (int)subids->elem_count;
191  subdomain_ids = (::types::subdomain_id *)(subids->array);
192 
193  for (i = 0; i < nsides; i++)
194  {
195  if (sides[i].is_hanging)
196  {
197  for (j = 0; j < 2; j++)
198  {
199  if (!sides[i].is.hanging.is_ghost[j])
200  {
201  typename ::parallel::distributed::
202  Triangulation<dim, spacedim>::cell_iterator cell =
203  cell_from_quad(triangulation,
204  sides[i].treeid,
205  *(sides[i].is.hanging.quad[j]));
206 
207  for (k = 0; k < nsubs; k++)
208  {
209  (*vertices_with_ghost_neighbors)
210  [cell->vertex_index(
211  p8est_edge_corners[sides[i].edge][1 ^ j])]
212  .insert(subdomain_ids[k]);
213  }
214  }
215  }
216  }
217  }
218 
219  subids->elem_count = 0;
220  }
221 
225  template <int dim, int spacedim>
226  void
227  find_ghosts_face(
228  typename ::internal::p4est::iter<dim>::face_info *info,
229  void * user_data)
230  {
231  int i, j, k;
232  int nsides = info->sides.elem_count;
233  typename ::internal::p4est::iter<dim>::face_side *sides =
234  (typename ::internal::p4est::iter<dim>::face_side *)(info->sides
235  .array);
236  FindGhosts<dim, spacedim> *fg =
237  static_cast<FindGhosts<dim, spacedim> *>(user_data);
238  sc_array_t *subids = fg->subids;
239  const ::parallel::distributed::Triangulation<dim, spacedim>
240  * triangulation = fg->triangulation;
241  int nsubs;
242  ::types::subdomain_id *subdomain_ids;
243  std::map<unsigned int, std::set<::types::subdomain_id>>
244  * vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
245  int limit = (dim == 2) ? 2 : 4;
246 
247  subids->elem_count = 0;
248  for (i = 0; i < nsides; i++)
249  {
250  if (sides[i].is_hanging)
251  {
252  for (j = 0; j < limit; j++)
253  {
254  if (sides[i].is.hanging.is_ghost[j])
255  {
256  typename ::parallel::distributed::
257  Triangulation<dim, spacedim>::cell_iterator cell =
258  cell_from_quad(triangulation,
259  sides[i].treeid,
260  *(sides[i].is.hanging.quad[j]));
261  ::types::subdomain_id *subid =
262  static_cast<::types::subdomain_id *>(
263  sc_array_push(subids));
264  *subid = cell->subdomain_id();
265  }
266  }
267  }
268  }
269 
270  if (!subids->elem_count)
271  {
272  return;
273  }
274 
275  nsubs = (int)subids->elem_count;
276  subdomain_ids = (::types::subdomain_id *)(subids->array);
277 
278  for (i = 0; i < nsides; i++)
279  {
280  if (sides[i].is_hanging)
281  {
282  for (j = 0; j < limit; j++)
283  {
284  if (!sides[i].is.hanging.is_ghost[j])
285  {
286  typename ::parallel::distributed::
287  Triangulation<dim, spacedim>::cell_iterator cell =
288  cell_from_quad(triangulation,
289  sides[i].treeid,
290  *(sides[i].is.hanging.quad[j]));
291 
292  for (k = 0; k < nsubs; k++)
293  {
294  if (dim == 2)
295  {
296  (*vertices_with_ghost_neighbors)
297  [cell->vertex_index(
298  p4est_face_corners[sides[i].face]
299  [(limit - 1) ^ j])]
300  .insert(subdomain_ids[k]);
301  }
302  else
303  {
304  (*vertices_with_ghost_neighbors)
305  [cell->vertex_index(
306  p8est_face_corners[sides[i].face]
307  [(limit - 1) ^ j])]
308  .insert(subdomain_ids[k]);
309  }
310  }
311  }
312  }
313  }
314  }
315 
316  subids->elem_count = 0;
317  }
318  } // namespace
319 
320 
321  int (&functions<2>::quadrant_compare)(const void *v1, const void *v2) =
322  p4est_quadrant_compare;
323 
324  void (&functions<2>::quadrant_childrenv)(const types<2>::quadrant *q,
325  types<2>::quadrant c[]) =
326  p4est_quadrant_childrenv;
327 
328  int (&functions<2>::quadrant_overlaps_tree)(types<2>::tree * tree,
329  const types<2>::quadrant *q) =
330  p4est_quadrant_overlaps_tree;
331 
332  void (&functions<2>::quadrant_set_morton)(types<2>::quadrant *quadrant,
333  int level,
334  uint64_t id) =
335  p4est_quadrant_set_morton;
336 
337  int (&functions<2>::quadrant_is_equal)(const types<2>::quadrant *q1,
338  const types<2>::quadrant *q2) =
339  p4est_quadrant_is_equal;
340 
341  int (&functions<2>::quadrant_is_sibling)(const types<2>::quadrant *q1,
342  const types<2>::quadrant *q2) =
343  p4est_quadrant_is_sibling;
344 
345  int (&functions<2>::quadrant_is_ancestor)(const types<2>::quadrant *q1,
346  const types<2>::quadrant *q2) =
347  p4est_quadrant_is_ancestor;
348 
349  int (&functions<2>::quadrant_ancestor_id)(const types<2>::quadrant *q,
350  int level) =
351  p4est_quadrant_ancestor_id;
352 
353  int (&functions<2>::comm_find_owner)(types<2>::forest * p4est,
354  const types<2>::locidx which_tree,
355  const types<2>::quadrant *q,
356  const int guess) =
357  p4est_comm_find_owner;
358 
359  types<2>::connectivity *(&functions<2>::connectivity_new)(
360  types<2>::topidx num_vertices,
361  types<2>::topidx num_trees,
362  types<2>::topidx num_corners,
363  types<2>::topidx num_vtt) = p4est_connectivity_new;
364 
365  void (&functions<2>::connectivity_join_faces)(types<2>::connectivity *conn,
366  types<2>::topidx tree_left,
367  types<2>::topidx tree_right,
368  int face_left,
369  int face_right,
370  int orientation) =
371  p4est_connectivity_join_faces;
372 
373  void (&functions<2>::connectivity_destroy)(
374  p4est_connectivity_t *connectivity) = p4est_connectivity_destroy;
375 
376  types<2>::forest *(&functions<2>::new_forest)(
377  MPI_Comm mpicomm,
378  types<2>::connectivity *connectivity,
379  types<2>::locidx min_quadrants,
380  int min_level,
381  int fill_uniform,
382  size_t data_size,
383  p4est_init_t init_fn,
384  void * user_pointer) = p4est_new_ext;
385 
386  void (&functions<2>::destroy)(types<2>::forest *p4est) = p4est_destroy;
387 
388  void (&functions<2>::refine)(types<2>::forest *p4est,
389  int refine_recursive,
390  p4est_refine_t refine_fn,
391  p4est_init_t init_fn) = p4est_refine;
392 
393  void (&functions<2>::coarsen)(types<2>::forest *p4est,
394  int coarsen_recursive,
395  p4est_coarsen_t coarsen_fn,
396  p4est_init_t init_fn) = p4est_coarsen;
397 
398  void (&functions<2>::balance)(types<2>::forest * p4est,
399  types<2>::balance_type btype,
400  p4est_init_t init_fn) = p4est_balance;
401 
402  types<2>::gloidx (&functions<2>::partition)(types<2>::forest *p4est,
403  int partition_for_coarsening,
404  p4est_weight_t weight_fn) =
405  p4est_partition_ext;
406 
407  void (&functions<2>::save)(const char * filename,
408  types<2>::forest *p4est,
409  int save_data) = p4est_save;
410 
411  types<2>::forest *(&functions<2>::load_ext)(
412  const char * filename,
413  MPI_Comm mpicomm,
414  std::size_t data_size,
415  int load_data,
416  int autopartition,
417  int broadcasthead,
418  void * user_pointer,
419  types<2>::connectivity **p4est) = p4est_load_ext;
420 
421  int (&functions<2>::connectivity_save)(
422  const char * filename,
423  types<2>::connectivity *connectivity) = p4est_connectivity_save;
424 
425  int (&functions<2>::connectivity_is_valid)(
426  types<2>::connectivity *connectivity) = p4est_connectivity_is_valid;
427 
428  types<2>::connectivity *(&functions<2>::connectivity_load)(
429  const char *filename,
430  size_t * length) = p4est_connectivity_load;
431 
432  unsigned int (&functions<2>::checksum)(types<2>::forest *p4est) =
433  p4est_checksum;
434 
435  void (&functions<2>::vtk_write_file)(types<2>::forest *p4est,
436  p4est_geometry_t *,
437  const char *baseName) =
438  p4est_vtk_write_file;
439 
440  types<2>::ghost *(&functions<2>::ghost_new)(types<2>::forest * p4est,
441  types<2>::balance_type btype) =
442  p4est_ghost_new;
443 
444  void (&functions<2>::ghost_destroy)(types<2>::ghost *ghost) =
445  p4est_ghost_destroy;
446 
447  void (&functions<2>::reset_data)(types<2>::forest *p4est,
448  size_t data_size,
449  p4est_init_t init_fn,
450  void *user_pointer) = p4est_reset_data;
451 
452  size_t (&functions<2>::forest_memory_used)(types<2>::forest *p4est) =
453  p4est_memory_used;
454 
455  size_t (&functions<2>::connectivity_memory_used)(
456  types<2>::connectivity *p4est) = p4est_connectivity_memory_used;
457 
458  template <int dim, int spacedim>
459  std::map<unsigned int, std::set<::types::subdomain_id>>
460  compute_vertices_with_ghost_neighbors(
461  const typename ::parallel::distributed::Triangulation<dim, spacedim>
462  & tria,
463  typename ::internal::p4est::types<dim>::forest *parallel_forest,
464  typename ::internal::p4est::types<dim>::ghost * parallel_ghost)
465  {
466  std::map<unsigned int, std::set<::types::subdomain_id>>
467  vertices_with_ghost_neighbors;
468 
469  ::internal::p4est::FindGhosts<dim, spacedim> fg;
470  fg.subids = sc_array_new(sizeof(::types::subdomain_id));
471  fg.triangulation = &tria;
472  fg.vertices_with_ghost_neighbors = &vertices_with_ghost_neighbors;
473 
474  switch (dim)
475  {
476  case 2:
477  p4est_iterate(
478  reinterpret_cast<::internal::p4est::types<2>::forest *>(
479  parallel_forest),
480  reinterpret_cast<::internal::p4est::types<2>::ghost *>(
481  parallel_ghost),
482  static_cast<void *>(&fg),
483  nullptr,
484  find_ghosts_face<2, spacedim>,
485  find_ghosts_corner<2, spacedim>);
486  break;
487 
488  case 3:
489  p8est_iterate(
490  reinterpret_cast<::internal::p4est::types<3>::forest *>(
491  parallel_forest),
492  reinterpret_cast<::internal::p4est::types<3>::ghost *>(
493  parallel_ghost),
494  static_cast<void *>(&fg),
495  nullptr,
496  find_ghosts_face<3, 3>,
497  find_ghosts_edge<3, 3>,
498  find_ghosts_corner<3, 3>);
499  break;
500 
501  default:
502  Assert(false, ExcNotImplemented());
503  }
504 
505  sc_array_destroy(fg.subids);
506 
507  return vertices_with_ghost_neighbors;
508  }
509 
510  constexpr unsigned int functions<2>::max_level;
511 
512  void (&functions<2>::transfer_fixed)(const types<2>::gloidx *dest_gfq,
513  const types<2>::gloidx *src_gfq,
514  MPI_Comm mpicomm,
515  int tag,
516  void * dest_data,
517  const void * src_data,
518  size_t data_size) =
519  p4est_transfer_fixed;
520 
521  types<2>::transfer_context *(&functions<2>::transfer_fixed_begin)(
522  const types<2>::gloidx *dest_gfq,
523  const types<2>::gloidx *src_gfq,
524  MPI_Comm mpicomm,
525  int tag,
526  void * dest_data,
527  const void * src_data,
528  size_t data_size) = p4est_transfer_fixed_begin;
529 
530  void (&functions<2>::transfer_fixed_end)(types<2>::transfer_context *tc) =
531  p4est_transfer_fixed_end;
532 
533  void (&functions<2>::transfer_custom)(const types<2>::gloidx *dest_gfq,
534  const types<2>::gloidx *src_gfq,
535  MPI_Comm mpicomm,
536  int tag,
537  void * dest_data,
538  const int * dest_sizes,
539  const void * src_data,
540  const int * src_sizes) =
541  p4est_transfer_custom;
542 
543  types<2>::transfer_context *(&functions<2>::transfer_custom_begin)(
544  const types<2>::gloidx *dest_gfq,
545  const types<2>::gloidx *src_gfq,
546  MPI_Comm mpicomm,
547  int tag,
548  void * dest_data,
549  const int * dest_sizes,
550  const void * src_data,
551  const int * src_sizes) = p4est_transfer_custom_begin;
552 
553  void (&functions<2>::transfer_custom_end)(types<2>::transfer_context *tc) =
554  p4est_transfer_custom_end;
555 
556 
557 
558  int (&functions<3>::quadrant_compare)(const void *v1, const void *v2) =
559  p8est_quadrant_compare;
560 
561  void (&functions<3>::quadrant_childrenv)(const types<3>::quadrant *q,
562  types<3>::quadrant c[]) =
563  p8est_quadrant_childrenv;
564 
565  int (&functions<3>::quadrant_overlaps_tree)(types<3>::tree * tree,
566  const types<3>::quadrant *q) =
567  p8est_quadrant_overlaps_tree;
568 
569  void (&functions<3>::quadrant_set_morton)(types<3>::quadrant *quadrant,
570  int level,
571  uint64_t id) =
572  p8est_quadrant_set_morton;
573 
574  int (&functions<3>::quadrant_is_equal)(const types<3>::quadrant *q1,
575  const types<3>::quadrant *q2) =
576  p8est_quadrant_is_equal;
577 
578  int (&functions<3>::quadrant_is_sibling)(const types<3>::quadrant *q1,
579  const types<3>::quadrant *q2) =
580  p8est_quadrant_is_sibling;
581 
582  int (&functions<3>::quadrant_is_ancestor)(const types<3>::quadrant *q1,
583  const types<3>::quadrant *q2) =
584  p8est_quadrant_is_ancestor;
585 
586  int (&functions<3>::quadrant_ancestor_id)(const types<3>::quadrant *q,
587  int level) =
588  p8est_quadrant_ancestor_id;
589 
590  int (&functions<3>::comm_find_owner)(types<3>::forest * p4est,
591  const types<3>::locidx which_tree,
592  const types<3>::quadrant *q,
593  const int guess) =
594  p8est_comm_find_owner;
595 
596  types<3>::connectivity *(&functions<3>::connectivity_new)(
597  types<3>::topidx num_vertices,
598  types<3>::topidx num_trees,
599  types<3>::topidx num_edges,
600  types<3>::topidx num_ett,
601  types<3>::topidx num_corners,
602  types<3>::topidx num_ctt) = p8est_connectivity_new;
603 
604  void (&functions<3>::connectivity_destroy)(
605  p8est_connectivity_t *connectivity) = p8est_connectivity_destroy;
606 
607  void (&functions<3>::connectivity_join_faces)(types<3>::connectivity *conn,
608  types<3>::topidx tree_left,
609  types<3>::topidx tree_right,
610  int face_left,
611  int face_right,
612  int orientation) =
613  p8est_connectivity_join_faces;
614 
615  types<3>::forest *(&functions<3>::new_forest)(
616  MPI_Comm mpicomm,
617  types<3>::connectivity *connectivity,
618  types<3>::locidx min_quadrants,
619  int min_level,
620  int fill_uniform,
621  size_t data_size,
622  p8est_init_t init_fn,
623  void * user_pointer) = p8est_new_ext;
624 
625  void (&functions<3>::destroy)(types<3>::forest *p8est) = p8est_destroy;
626 
627  void (&functions<3>::refine)(types<3>::forest *p8est,
628  int refine_recursive,
629  p8est_refine_t refine_fn,
630  p8est_init_t init_fn) = p8est_refine;
631 
632  void (&functions<3>::coarsen)(types<3>::forest *p8est,
633  int coarsen_recursive,
634  p8est_coarsen_t coarsen_fn,
635  p8est_init_t init_fn) = p8est_coarsen;
636 
637  void (&functions<3>::balance)(types<3>::forest * p8est,
638  types<3>::balance_type btype,
639  p8est_init_t init_fn) = p8est_balance;
640 
641  types<3>::gloidx (&functions<3>::partition)(types<3>::forest *p8est,
642  int partition_for_coarsening,
643  p8est_weight_t weight_fn) =
644  p8est_partition_ext;
645 
646  void (&functions<3>::save)(const char * filename,
647  types<3>::forest *p4est,
648  int save_data) = p8est_save;
649 
650  types<3>::forest *(&functions<3>::load_ext)(
651  const char * filename,
652  MPI_Comm mpicomm,
653  std::size_t data_size,
654  int load_data,
655  int autopartition,
656  int broadcasthead,
657  void * user_pointer,
658  types<3>::connectivity **p4est) = p8est_load_ext;
659 
660  int (&functions<3>::connectivity_save)(
661  const char * filename,
662  types<3>::connectivity *connectivity) = p8est_connectivity_save;
663 
664  int (&functions<3>::connectivity_is_valid)(
665  types<3>::connectivity *connectivity) = p8est_connectivity_is_valid;
666 
667  types<3>::connectivity *(&functions<3>::connectivity_load)(
668  const char *filename,
669  size_t * length) = p8est_connectivity_load;
670 
671  unsigned int (&functions<3>::checksum)(types<3>::forest *p8est) =
672  p8est_checksum;
673 
674  void (&functions<3>::vtk_write_file)(types<3>::forest *p8est,
675  p8est_geometry_t *,
676  const char *baseName) =
677  p8est_vtk_write_file;
678 
679  types<3>::ghost *(&functions<3>::ghost_new)(types<3>::forest * p4est,
680  types<3>::balance_type btype) =
681  p8est_ghost_new;
682 
683  void (&functions<3>::ghost_destroy)(types<3>::ghost *ghost) =
684  p8est_ghost_destroy;
685 
686  void (&functions<3>::reset_data)(types<3>::forest *p4est,
687  size_t data_size,
688  p8est_init_t init_fn,
689  void *user_pointer) = p8est_reset_data;
690 
691  size_t (&functions<3>::forest_memory_used)(types<3>::forest *p4est) =
692  p8est_memory_used;
693 
694  size_t (&functions<3>::connectivity_memory_used)(
695  types<3>::connectivity *p4est) = p8est_connectivity_memory_used;
696 
697  constexpr unsigned int functions<3>::max_level;
698 
699  void (&functions<3>::transfer_fixed)(const types<3>::gloidx *dest_gfq,
700  const types<3>::gloidx *src_gfq,
701  MPI_Comm mpicomm,
702  int tag,
703  void * dest_data,
704  const void * src_data,
705  size_t data_size) =
706  p8est_transfer_fixed;
707 
708  types<3>::transfer_context *(&functions<3>::transfer_fixed_begin)(
709  const types<3>::gloidx *dest_gfq,
710  const types<3>::gloidx *src_gfq,
711  MPI_Comm mpicomm,
712  int tag,
713  void * dest_data,
714  const void * src_data,
715  size_t data_size) = p8est_transfer_fixed_begin;
716 
717  void (&functions<3>::transfer_fixed_end)(types<3>::transfer_context *tc) =
718  p8est_transfer_fixed_end;
719 
720  void (&functions<3>::transfer_custom)(const types<3>::gloidx *dest_gfq,
721  const types<3>::gloidx *src_gfq,
722  MPI_Comm mpicomm,
723  int tag,
724  void * dest_data,
725  const int * dest_sizes,
726  const void * src_data,
727  const int * src_sizes) =
728  p8est_transfer_custom;
729 
730  types<3>::transfer_context *(&functions<3>::transfer_custom_begin)(
731  const types<3>::gloidx *dest_gfq,
732  const types<3>::gloidx *src_gfq,
733  MPI_Comm mpicomm,
734  int tag,
735  void * dest_data,
736  const int * dest_sizes,
737  const void * src_data,
738  const int * src_sizes) = p8est_transfer_custom_begin;
739 
740  void (&functions<3>::transfer_custom_end)(types<3>::transfer_context *tc) =
741  p8est_transfer_custom_end;
742 
743 
744 
745  template <int dim>
746  void
747  init_quadrant_children(
748  const typename types<dim>::quadrant &p4est_cell,
749  typename types<dim>::quadrant (
750  &p4est_children)[::GeometryInfo<dim>::max_children_per_cell])
751  {
752  for (unsigned int c = 0;
753  c < ::GeometryInfo<dim>::max_children_per_cell;
754  ++c)
755  switch (dim)
756  {
757  case 2:
758  P4EST_QUADRANT_INIT(&p4est_children[c]);
759  break;
760  case 3:
761  P8EST_QUADRANT_INIT(&p4est_children[c]);
762  break;
763  default:
764  Assert(false, ExcNotImplemented());
765  }
766 
767 
768  functions<dim>::quadrant_childrenv(&p4est_cell, p4est_children);
769  }
770 
771  template <int dim>
772  void
773  init_coarse_quadrant(typename types<dim>::quadrant &quad)
774  {
775  switch (dim)
776  {
777  case 2:
778  P4EST_QUADRANT_INIT(&quad);
779  break;
780  case 3:
781  P8EST_QUADRANT_INIT(&quad);
782  break;
783  default:
784  Assert(false, ExcNotImplemented());
785  }
786  functions<dim>::quadrant_set_morton(&quad,
787  /*level=*/0,
788  /*index=*/0);
789  }
790 
791  template <int dim>
792  bool
793  quadrant_is_equal(const typename types<dim>::quadrant &q1,
794  const typename types<dim>::quadrant &q2)
795  {
796  return functions<dim>::quadrant_is_equal(&q1, &q2);
797  }
798 
799 
800 
801  template <int dim>
802  bool
803  quadrant_is_ancestor(const typename types<dim>::quadrant &q1,
804  const typename types<dim>::quadrant &q2)
805  {
806  return functions<dim>::quadrant_is_ancestor(&q1, &q2);
807  }
808 
809  template <int dim>
810  bool
811  tree_exists_locally(const typename types<dim>::forest *parallel_forest,
812  const typename types<dim>::topidx coarse_grid_cell)
813  {
814  Assert(coarse_grid_cell < parallel_forest->connectivity->num_trees,
815  ExcInternalError());
816  return ((coarse_grid_cell >= parallel_forest->first_local_tree) &&
817  (coarse_grid_cell <= parallel_forest->last_local_tree));
818  }
819  } // namespace p4est
820 } // namespace internal
821 
822 #endif // DEAL_II_WITH_P4EST
823 
824 /*-------------- Explicit Instantiations -------------------------------*/
825 #include "p4est_wrappers.inst"
826 
827 
828 DEAL_II_NAMESPACE_CLOSE
unsigned long long int global_dof_index
Definition: types.h:72
static::ExceptionBase & ExcMessage(std::string arg1)
unsigned int subdomain_id
Definition: types.h:43
#define Assert(cond, exc)
Definition: exceptions.h:1227
static::ExceptionBase & ExcNotImplemented()
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
static::ExceptionBase & ExcInternalError()