Reference documentation for deal.II version GIT 29f9da0a34 2023-12-07 10:00:01+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.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2022 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 
19 
21 
22 #ifdef DEAL_II_WITH_P4EST
23 
24 namespace internal
25 {
26  namespace p4est
27  {
28  namespace
29  {
30  template <int dim, int spacedim>
31  typename ::Triangulation<dim, spacedim>::cell_iterator
32  cell_from_quad(
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)
37  {
38  int i, l = quad.level;
39  ::types::global_dof_index dealii_index =
40  triangulation->get_p4est_tree_to_coarse_cell_permutation()[treeidx];
41 
42  for (i = 0; i < l; ++i)
43  {
44  typename ::Triangulation<dim, spacedim>::cell_iterator cell(
45  triangulation, i, dealii_index);
46  const int child_id =
48  &quad, i + 1);
49  Assert(cell->has_children(),
50  ExcMessage("p4est quadrant does not correspond to a cell!"));
51  dealii_index = cell->child_index(child_id);
52  }
53 
54  typename ::Triangulation<dim, spacedim>::cell_iterator out_cell(
55  triangulation, l, dealii_index);
56 
57  return out_cell;
58  }
59 
64  template <int dim, int spacedim>
65  struct FindGhosts
66  {
67  const typename ::parallel::distributed::Triangulation<dim,
68  spacedim>
70  sc_array_t *subids;
71  std::map<unsigned int, std::set<::types::subdomain_id>>
73  };
74 
75 
81  template <int dim, int spacedim>
82  void
83  find_ghosts_corner(
84  typename ::internal::p4est::iter<dim>::corner_info *info,
85  void *user_data)
86  {
87  int i, j;
88  int nsides = info->sides.elem_count;
89  auto *sides = reinterpret_cast<
90  typename ::internal::p4est::iter<dim>::corner_side *>(
91  info->sides.array);
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>
96  *triangulation = fg->triangulation;
97  int nsubs;
98  ::types::subdomain_id *subdomain_ids;
99  std::map<unsigned int, std::set<::types::subdomain_id>>
100  *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
101 
102  subids->elem_count = 0;
103  for (i = 0; i < nsides; ++i)
104  {
105  if (sides[i].is_ghost)
106  {
107  typename ::parallel::distributed::
108  Triangulation<dim, spacedim>::cell_iterator cell =
109  cell_from_quad(triangulation,
110  sides[i].treeid,
111  *(sides[i].quad));
112  Assert(cell->is_ghost(),
113  ExcMessage("ghost quad did not find ghost cell"));
114  ::types::subdomain_id *subid =
115  static_cast<::types::subdomain_id *>(
116  sc_array_push(subids));
117  *subid = cell->subdomain_id();
118  }
119  }
120 
121  if (!subids->elem_count)
122  {
123  return;
124  }
125 
126  nsubs = static_cast<int>(subids->elem_count);
127  subdomain_ids =
128  reinterpret_cast<::types::subdomain_id *>(subids->array);
129 
130  for (i = 0; i < nsides; ++i)
131  {
132  if (!sides[i].is_ghost)
133  {
134  typename ::parallel::distributed::
135  Triangulation<dim, spacedim>::cell_iterator cell =
136  cell_from_quad(triangulation,
137  sides[i].treeid,
138  *(sides[i].quad));
139 
140  Assert(!cell->is_ghost(),
141  ExcMessage("local quad found ghost cell"));
142 
143  for (j = 0; j < nsubs; ++j)
144  {
145  (*vertices_with_ghost_neighbors)[cell->vertex_index(
146  sides[i].corner)]
147  .insert(subdomain_ids[j]);
148  }
149  }
150  }
151 
152  subids->elem_count = 0;
153  }
154 
158  template <int dim, int spacedim>
159  void
160  find_ghosts_edge(
161  typename ::internal::p4est::iter<dim>::edge_info *info,
162  void *user_data)
163  {
164  int i, j, k;
165  int nsides = info->sides.elem_count;
166  auto *sides = reinterpret_cast<
167  typename ::internal::p4est::iter<dim>::edge_side *>(
168  info->sides.array);
169  auto *fg = static_cast<FindGhosts<dim, spacedim> *>(user_data);
170  sc_array_t *subids = fg->subids;
171  const ::parallel::distributed::Triangulation<dim, spacedim>
172  *triangulation = fg->triangulation;
173  int nsubs;
174  ::types::subdomain_id *subdomain_ids;
175  std::map<unsigned int, std::set<::types::subdomain_id>>
176  *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
177 
178  subids->elem_count = 0;
179  for (i = 0; i < nsides; ++i)
180  {
181  if (sides[i].is_hanging)
182  {
183  for (j = 0; j < 2; ++j)
184  {
185  if (sides[i].is.hanging.is_ghost[j])
186  {
187  typename ::parallel::distributed::
188  Triangulation<dim, spacedim>::cell_iterator cell =
189  cell_from_quad(triangulation,
190  sides[i].treeid,
191  *(sides[i].is.hanging.quad[j]));
192  ::types::subdomain_id *subid =
193  static_cast<::types::subdomain_id *>(
194  sc_array_push(subids));
195  *subid = cell->subdomain_id();
196  }
197  }
198  }
199  }
200 
201  if (!subids->elem_count)
202  {
203  return;
204  }
205 
206  nsubs = static_cast<int>(subids->elem_count);
207  subdomain_ids =
208  reinterpret_cast<::types::subdomain_id *>(subids->array);
209 
210  for (i = 0; i < nsides; ++i)
211  {
212  if (sides[i].is_hanging)
213  {
214  for (j = 0; j < 2; ++j)
215  {
216  if (!sides[i].is.hanging.is_ghost[j])
217  {
218  typename ::parallel::distributed::
219  Triangulation<dim, spacedim>::cell_iterator cell =
220  cell_from_quad(triangulation,
221  sides[i].treeid,
222  *(sides[i].is.hanging.quad[j]));
223 
224  for (k = 0; k < nsubs; ++k)
225  {
226  (*vertices_with_ghost_neighbors)
227  [cell->vertex_index(
228  p8est_edge_corners[sides[i].edge][1 ^ j])]
229  .insert(subdomain_ids[k]);
230  }
231  }
232  }
233  }
234  }
235 
236  subids->elem_count = 0;
237  }
238 
242  template <int dim, int spacedim>
243  void
244  find_ghosts_face(
245  typename ::internal::p4est::iter<dim>::face_info *info,
246  void *user_data)
247  {
248  int i, j, k;
249  int nsides = info->sides.elem_count;
250  auto *sides = reinterpret_cast<
251  typename ::internal::p4est::iter<dim>::face_side *>(
252  info->sides.array);
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>
257  *triangulation = fg->triangulation;
258  int nsubs;
259  ::types::subdomain_id *subdomain_ids;
260  std::map<unsigned int, std::set<::types::subdomain_id>>
261  *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
262  int limit = (dim == 2) ? 2 : 4;
263 
264  subids->elem_count = 0;
265  for (i = 0; i < nsides; ++i)
266  {
267  if (sides[i].is_hanging)
268  {
269  for (j = 0; j < limit; ++j)
270  {
271  if (sides[i].is.hanging.is_ghost[j])
272  {
273  typename ::parallel::distributed::
274  Triangulation<dim, spacedim>::cell_iterator cell =
275  cell_from_quad(triangulation,
276  sides[i].treeid,
277  *(sides[i].is.hanging.quad[j]));
278  ::types::subdomain_id *subid =
279  static_cast<::types::subdomain_id *>(
280  sc_array_push(subids));
281  *subid = cell->subdomain_id();
282  }
283  }
284  }
285  }
286 
287  if (!subids->elem_count)
288  {
289  return;
290  }
291 
292  nsubs = static_cast<int>(subids->elem_count);
293  subdomain_ids =
294  reinterpret_cast<::types::subdomain_id *>(subids->array);
295 
296  for (i = 0; i < nsides; ++i)
297  {
298  if (sides[i].is_hanging)
299  {
300  for (j = 0; j < limit; ++j)
301  {
302  if (!sides[i].is.hanging.is_ghost[j])
303  {
304  typename ::parallel::distributed::
305  Triangulation<dim, spacedim>::cell_iterator cell =
306  cell_from_quad(triangulation,
307  sides[i].treeid,
308  *(sides[i].is.hanging.quad[j]));
309 
310  for (k = 0; k < nsubs; ++k)
311  {
312  if (dim == 2)
313  {
314  (*vertices_with_ghost_neighbors)
315  [cell->vertex_index(
316  p4est_face_corners[sides[i].face]
317  [(limit - 1) ^ j])]
318  .insert(subdomain_ids[k]);
319  }
320  else
321  {
322  (*vertices_with_ghost_neighbors)
323  [cell->vertex_index(
324  p8est_face_corners[sides[i].face]
325  [(limit - 1) ^ j])]
326  .insert(subdomain_ids[k]);
327  }
328  }
329  }
330  }
331  }
332  }
333 
334  subids->elem_count = 0;
335  }
336  } // namespace
337 
338 
339  int (&functions<2>::quadrant_compare)(const void *v1, const void *v2) =
340  p4est_quadrant_compare;
341 
343  types<2>::quadrant c[]) =
344  p4est_quadrant_childrenv;
345 
347  const types<2>::quadrant *q) =
348  p4est_quadrant_overlaps_tree;
349 
351  int level,
352  std::uint64_t id) =
353  p4est_quadrant_set_morton;
354 
356  const types<2>::quadrant *q2) =
357  p4est_quadrant_is_equal;
358 
360  const types<2>::quadrant *q2) =
361  p4est_quadrant_is_sibling;
362 
364  const types<2>::quadrant *q2) =
365  p4est_quadrant_is_ancestor;
366 
368  int level) =
369  p4est_quadrant_ancestor_id;
370 
372  const types<2>::locidx which_tree,
373  const types<2>::quadrant *q,
374  const int guess) =
375  p4est_comm_find_owner;
376 
378  types<2>::topidx num_vertices,
379  types<2>::topidx num_trees,
380  types<2>::topidx num_corners,
381  types<2>::topidx num_vtt) = p4est_connectivity_new;
382 
384  types<2>::topidx num_vertices,
385  types<2>::topidx num_trees,
386  types<2>::topidx num_corners,
387  const double *vertices,
388  const types<2>::topidx *ttv,
389  const types<2>::topidx *ttt,
390  const int8_t *ttf,
391  const types<2>::topidx *ttc,
392  const types<2>::topidx *coff,
393  const types<2>::topidx *ctt,
394  const int8_t *ctc) = p4est_connectivity_new_copy;
395 
397  types<2>::topidx tree_left,
398  types<2>::topidx tree_right,
399  int face_left,
400  int face_right,
401  int orientation) =
402  p4est_connectivity_join_faces;
403 
405  p4est_connectivity_t *connectivity) = p4est_connectivity_destroy;
406 
408  MPI_Comm mpicomm,
409  types<2>::connectivity *connectivity,
410  types<2>::locidx min_quadrants,
411  int min_level,
412  int fill_uniform,
413  std::size_t data_size,
414  p4est_init_t init_fn,
415  void *user_pointer) = p4est_new_ext;
416 
418  int copy_data) = p4est_copy;
419 
420  void (&functions<2>::destroy)(types<2>::forest *p4est) = p4est_destroy;
421 
422  void (&functions<2>::refine)(types<2>::forest *p4est,
423  int refine_recursive,
424  p4est_refine_t refine_fn,
425  p4est_init_t init_fn) = p4est_refine;
426 
428  int coarsen_recursive,
429  p4est_coarsen_t coarsen_fn,
430  p4est_init_t init_fn) = p4est_coarsen;
431 
434  p4est_init_t init_fn) = p4est_balance;
435 
437  int partition_for_coarsening,
438  p4est_weight_t weight_fn) =
439  p4est_partition_ext;
440 
441  void (&functions<2>::save)(const char *filename,
442  types<2>::forest *p4est,
443  int save_data) = p4est_save;
444 
446  const char *filename,
447  MPI_Comm mpicomm,
448  std::size_t data_size,
449  int load_data,
450  int autopartition,
451  int broadcasthead,
452  void *user_pointer,
453  types<2>::connectivity **p4est) = p4est_load_ext;
454 
456  const char *filename,
457  types<2>::connectivity *connectivity) = p4est_connectivity_save;
458 
460  types<2>::connectivity *connectivity) = p4est_connectivity_is_valid;
461 
463  const char *filename,
464  std::size_t *length) = p4est_connectivity_load;
465 
466  unsigned int (&functions<2>::checksum)(types<2>::forest *p4est) =
467  p4est_checksum;
468 
470  p4est_geometry_t *,
471  const char *baseName) =
472  p4est_vtk_write_file;
473 
475  types<2>::balance_type btype) =
476  p4est_ghost_new;
477 
479  p4est_ghost_destroy;
480 
482  std::size_t data_size,
483  p4est_init_t init_fn,
484  void *user_pointer) = p4est_reset_data;
485 
486  std::size_t (&functions<2>::forest_memory_used)(types<2>::forest *p4est) =
487  p4est_memory_used;
488 
490  types<2>::connectivity *p4est) = p4est_connectivity_memory_used;
491 
492  constexpr unsigned int functions<2>::max_level;
493 
494  void (&functions<2>::transfer_fixed)(const types<2>::gloidx *dest_gfq,
495  const types<2>::gloidx *src_gfq,
496  MPI_Comm mpicomm,
497  int tag,
498  void *dest_data,
499  const void *src_data,
500  std::size_t data_size) =
501  p4est_transfer_fixed;
502 
504  const types<2>::gloidx *dest_gfq,
505  const types<2>::gloidx *src_gfq,
506  MPI_Comm mpicomm,
507  int tag,
508  void *dest_data,
509  const void *src_data,
510  std::size_t data_size) = p4est_transfer_fixed_begin;
511 
513  p4est_transfer_fixed_end;
514 
515  void (&functions<2>::transfer_custom)(const types<2>::gloidx *dest_gfq,
516  const types<2>::gloidx *src_gfq,
517  MPI_Comm mpicomm,
518  int tag,
519  void *dest_data,
520  const int *dest_sizes,
521  const void *src_data,
522  const int *src_sizes) =
523  p4est_transfer_custom;
524 
526  const types<2>::gloidx *dest_gfq,
527  const types<2>::gloidx *src_gfq,
528  MPI_Comm mpicomm,
529  int tag,
530  void *dest_data,
531  const int *dest_sizes,
532  const void *src_data,
533  const int *src_sizes) = p4est_transfer_custom_begin;
534 
536  p4est_transfer_custom_end;
537 
538 # ifdef P4EST_SEARCH_LOCAL
540  types<2>::forest *p4est,
541  int call_post,
544  sc_array_t *points) = p4est_search_partition;
545 # endif
546 
548  types<2>::connectivity *connectivity,
549  types<2>::topidx treeid,
552  double vxyz[3]) = p4est_qcoord_to_vertex;
553 
554  int (&functions<3>::quadrant_compare)(const void *v1, const void *v2) =
555  p8est_quadrant_compare;
556 
558  types<3>::quadrant c[]) =
559  p8est_quadrant_childrenv;
560 
562  const types<3>::quadrant *q) =
563  p8est_quadrant_overlaps_tree;
564 
566  int level,
567  std::uint64_t id) =
568  p8est_quadrant_set_morton;
569 
571  const types<3>::quadrant *q2) =
572  p8est_quadrant_is_equal;
573 
575  const types<3>::quadrant *q2) =
576  p8est_quadrant_is_sibling;
577 
579  const types<3>::quadrant *q2) =
580  p8est_quadrant_is_ancestor;
581 
583  int level) =
584  p8est_quadrant_ancestor_id;
585 
587  const types<3>::locidx which_tree,
588  const types<3>::quadrant *q,
589  const int guess) =
590  p8est_comm_find_owner;
591 
593  types<3>::topidx num_vertices,
594  types<3>::topidx num_trees,
595  types<3>::topidx num_edges,
596  types<3>::topidx num_ett,
597  types<3>::topidx num_corners,
598  types<3>::topidx num_ctt) = p8est_connectivity_new;
599 
601  types<3>::topidx num_vertices,
602  types<3>::topidx num_trees,
603  types<3>::topidx num_edges,
604  types<3>::topidx num_corners,
605  const double *vertices,
606  const types<3>::topidx *ttv,
607  const types<3>::topidx *ttt,
608  const int8_t *ttf,
609  const types<3>::topidx *tte,
610  const types<3>::topidx *eoff,
611  const types<3>::topidx *ett,
612  const int8_t *ete,
613  const types<3>::topidx *ttc,
614  const types<3>::topidx *coff,
615  const types<3>::topidx *ctt,
616  const int8_t *ctc) = p8est_connectivity_new_copy;
617 
619  p8est_connectivity_t *connectivity) = p8est_connectivity_destroy;
620 
622  types<3>::topidx tree_left,
623  types<3>::topidx tree_right,
624  int face_left,
625  int face_right,
626  int orientation) =
627  p8est_connectivity_join_faces;
628 
630  MPI_Comm mpicomm,
631  types<3>::connectivity *connectivity,
632  types<3>::locidx min_quadrants,
633  int min_level,
634  int fill_uniform,
635  std::size_t data_size,
636  p8est_init_t init_fn,
637  void *user_pointer) = p8est_new_ext;
638 
640  int copy_data) = p8est_copy;
641 
642  void (&functions<3>::destroy)(types<3>::forest *p8est) = p8est_destroy;
643 
644  void (&functions<3>::refine)(types<3>::forest *p8est,
645  int refine_recursive,
646  p8est_refine_t refine_fn,
647  p8est_init_t init_fn) = p8est_refine;
648 
650  int coarsen_recursive,
651  p8est_coarsen_t coarsen_fn,
652  p8est_init_t init_fn) = p8est_coarsen;
653 
656  p8est_init_t init_fn) = p8est_balance;
657 
659  int partition_for_coarsening,
660  p8est_weight_t weight_fn) =
661  p8est_partition_ext;
662 
663  void (&functions<3>::save)(const char *filename,
664  types<3>::forest *p4est,
665  int save_data) = p8est_save;
666 
668  const char *filename,
669  MPI_Comm mpicomm,
670  std::size_t data_size,
671  int load_data,
672  int autopartition,
673  int broadcasthead,
674  void *user_pointer,
675  types<3>::connectivity **p4est) = p8est_load_ext;
676 
678  const char *filename,
679  types<3>::connectivity *connectivity) = p8est_connectivity_save;
680 
682  types<3>::connectivity *connectivity) = p8est_connectivity_is_valid;
683 
685  const char *filename,
686  std::size_t *length) = p8est_connectivity_load;
687 
688  unsigned int (&functions<3>::checksum)(types<3>::forest *p8est) =
689  p8est_checksum;
690 
692  p8est_geometry_t *,
693  const char *baseName) =
694  p8est_vtk_write_file;
695 
697  types<3>::balance_type btype) =
698  p8est_ghost_new;
699 
701  p8est_ghost_destroy;
702 
704  std::size_t data_size,
705  p8est_init_t init_fn,
706  void *user_pointer) = p8est_reset_data;
707 
708  std::size_t (&functions<3>::forest_memory_used)(types<3>::forest *p4est) =
709  p8est_memory_used;
710 
712  types<3>::connectivity *p4est) = p8est_connectivity_memory_used;
713 
714  constexpr unsigned int functions<3>::max_level;
715 
716  void (&functions<3>::transfer_fixed)(const types<3>::gloidx *dest_gfq,
717  const types<3>::gloidx *src_gfq,
718  MPI_Comm mpicomm,
719  int tag,
720  void *dest_data,
721  const void *src_data,
722  std::size_t data_size) =
723  p8est_transfer_fixed;
724 
726  const types<3>::gloidx *dest_gfq,
727  const types<3>::gloidx *src_gfq,
728  MPI_Comm mpicomm,
729  int tag,
730  void *dest_data,
731  const void *src_data,
732  std::size_t data_size) = p8est_transfer_fixed_begin;
733 
735  p8est_transfer_fixed_end;
736 
737  void (&functions<3>::transfer_custom)(const types<3>::gloidx *dest_gfq,
738  const types<3>::gloidx *src_gfq,
739  MPI_Comm mpicomm,
740  int tag,
741  void *dest_data,
742  const int *dest_sizes,
743  const void *src_data,
744  const int *src_sizes) =
745  p8est_transfer_custom;
746 
748  const types<3>::gloidx *dest_gfq,
749  const types<3>::gloidx *src_gfq,
750  MPI_Comm mpicomm,
751  int tag,
752  void *dest_data,
753  const int *dest_sizes,
754  const void *src_data,
755  const int *src_sizes) = p8est_transfer_custom_begin;
756 
758  p8est_transfer_custom_end;
759 
760 # ifdef P4EST_SEARCH_LOCAL
762  types<3>::forest *p4est,
763  int call_post,
766  sc_array_t *points) = p8est_search_partition;
767 # endif
768 
770  types<3>::connectivity *connectivity,
771  types<3>::topidx treeid,
775  double vxyz[3]) = p8est_qcoord_to_vertex;
776 
777  template <int dim>
778  void
780  const typename types<dim>::quadrant &p4est_cell,
781  typename types<dim>::quadrant (
782  &p4est_children)[::GeometryInfo<dim>::max_children_per_cell])
783  {
784  for (unsigned int c = 0;
785  c < ::GeometryInfo<dim>::max_children_per_cell;
786  ++c)
787  switch (dim)
788  {
789  case 2:
790  P4EST_QUADRANT_INIT(&p4est_children[c]);
791  break;
792  case 3:
793  P8EST_QUADRANT_INIT(&p4est_children[c]);
794  break;
795  default:
796  Assert(false, ExcNotImplemented());
797  }
798 
799 
800  functions<dim>::quadrant_childrenv(&p4est_cell, p4est_children);
801  }
802 
803  template <int dim>
804  void
806  {
807  switch (dim)
808  {
809  case 2:
810  P4EST_QUADRANT_INIT(&quad);
811  break;
812  case 3:
813  P8EST_QUADRANT_INIT(&quad);
814  break;
815  default:
816  Assert(false, ExcNotImplemented());
817  }
819  /*level=*/0,
820  /*index=*/0);
821  }
822 
823  template <int dim>
824  bool
826  const typename types<dim>::quadrant &q2)
827  {
828  return functions<dim>::quadrant_is_equal(&q1, &q2);
829  }
830 
831 
832 
833  template <int dim>
834  bool
836  const typename types<dim>::quadrant &q2)
837  {
838  return functions<dim>::quadrant_is_ancestor(&q1, &q2);
839  }
840 
841  template <int dim>
842  bool
843  tree_exists_locally(const typename types<dim>::forest *parallel_forest,
844  const typename types<dim>::topidx coarse_grid_cell)
845  {
846  Assert(coarse_grid_cell < parallel_forest->connectivity->num_trees,
847  ExcInternalError());
848  return ((coarse_grid_cell >= parallel_forest->first_local_tree) &&
849  (coarse_grid_cell <= parallel_forest->last_local_tree));
850  }
851 
852 
853 
854  // template specializations
855 
856  template <>
858  copy_connectivity<2>(const typename types<2>::connectivity *connectivity)
859  {
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);
872  }
873 
874  template <>
876  copy_connectivity<3>(const typename types<3>::connectivity *connectivity)
877  {
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);
895  }
896 
897 
898 
899  template <>
900  bool
901  quadrant_is_equal<1>(const typename types<1>::quadrant &q1,
902  const typename types<1>::quadrant &q2)
903  {
904  return q1 == q2;
905  }
906 
907 
908 
909  template <>
910  bool
912  const types<1>::quadrant &q2)
913  {
914  // determine level of quadrants
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) >>
919 
920  // q1 can be an ancestor of q2 if q1's level is smaller
921  if (level_1 >= level_2)
922  return false;
923 
924  // extract path of quadrants up to level of possible ancestor q1
925  const int truncated_id_1 = (q1 >> (types<1>::n_bits - 1 - level_1))
926  << (types<1>::n_bits - 1 - level_1);
927  const int truncated_id_2 = (q2 >> (types<1>::n_bits - 1 - level_1))
928  << (types<1>::n_bits - 1 - level_1);
929 
930  // compare paths
931  return truncated_id_1 == truncated_id_2;
932  }
933 
934 
935 
936  template <>
937  void
939  const typename types<1>::quadrant &q,
940  typename types<1>::quadrant (
941  &p4est_children)[::GeometryInfo<1>::max_children_per_cell])
942  {
943  // determine the current level of quadrant
944  const int level_parent = (q << types<1>::max_n_child_indices_bits) >>
946  const int level_child = level_parent + 1;
947 
948  // left child: only n_child_indices has to be incremented
949  p4est_children[0] = (q + 1);
950 
951  // right child: increment and set a bit to 1 indicating that it is a right
952  // child
953  p4est_children[1] = (q + 1) | (1 << (types<1>::n_bits - 1 - level_child));
954  }
955 
956 
957 
958  template <>
959  void
961  {
962  quad = 0;
963  }
964 
965  } // namespace p4est
966 } // namespace internal
967 
968 #endif // DEAL_II_WITH_P4EST
969 
970 /*-------------- Explicit Instantiations -------------------------------*/
971 #include "p4est_wrappers.inst"
972 
973 
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:477
#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
static ::ExceptionBase & ExcInternalError()
#define Assert(cond, exc)
Definition: exceptions.h:1631
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcMessage(std::string arg1)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
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_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)
bool quadrant_is_ancestor< 1 >(const types< 1 >::quadrant &q1, const types< 1 >::quadrant &q2)
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)
unsigned int global_dof_index
Definition: types.h:82
unsigned int subdomain_id
Definition: types.h:44
sc_array_t * subids
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
std::map< unsigned int, std::set<::types::subdomain_id > > * vertices_with_ghost_neighbors