Reference documentation for deal.II version 8.4.2
tria.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2016 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/base/utilities.h>
18 #include <deal.II/base/memory_consumption.h>
19 #include <deal.II/base/logstream.h>
20 #include <deal.II/lac/sparsity_tools.h>
21 #include <deal.II/lac/dynamic_sparsity_pattern.h>
22 #include <deal.II/grid/tria.h>
23 #include <deal.II/grid/tria_accessor.h>
24 #include <deal.II/grid/tria_iterator.h>
25 #include <deal.II/grid/grid_tools.h>
26 #include <deal.II/distributed/tria.h>
27 
28 #ifdef DEAL_II_WITH_P4EST
29 # include <p4est_bits.h>
30 # include <p4est_extended.h>
31 # include <p4est_vtk.h>
32 # include <p4est_ghost.h>
33 # include <p4est_communication.h>
34 # include <p4est_iterate.h>
35 
36 # include <p8est_bits.h>
37 # include <p8est_extended.h>
38 # include <p8est_vtk.h>
39 # include <p8est_ghost.h>
40 # include <p8est_communication.h>
41 # include <p8est_iterate.h>
42 #endif
43 
44 #include <algorithm>
45 #include <numeric>
46 #include <iostream>
47 #include <fstream>
48 
49 
50 DEAL_II_NAMESPACE_OPEN
51 
52 
53 #ifdef DEAL_II_WITH_P4EST
54 
55 namespace internal
56 {
57  namespace p4est
58  {
66  template <int dim> struct functions;
67 
68  template <> struct functions<2>
69  {
70  static
71  int (&quadrant_compare) (const void *v1, const void *v2);
72 
73  static
74  void (&quadrant_childrenv) (const types<2>::quadrant *q,
75  types<2>::quadrant c[]);
76 
77  static
78  int (&quadrant_overlaps_tree) (types<2>::tree *tree,
79  const types<2>::quadrant *q);
80 
81  static
82  void (&quadrant_set_morton) (types<2>::quadrant *quadrant,
83  int level,
84  uint64_t id);
85 
86  static
87  int (&quadrant_is_equal) (const types<2>::quadrant *q1,
88  const types<2>::quadrant *q2);
89 
90  static
91  int (&quadrant_is_sibling) (const types<2>::quadrant *q1,
92  const types<2>::quadrant *q2);
93 
94  static
95  int (&quadrant_is_ancestor) (const types<2>::quadrant *q1,
96  const types<2>::quadrant *q2);
97 
98  static
99  int (&quadrant_ancestor_id) (const types<2>::quadrant *q,
100  int level);
101 
102  static
103  int (&comm_find_owner) (types<2>::forest *p4est,
104  const types<2>::locidx which_tree,
105  const types<2>::quadrant *q,
106  const int guess);
107 
108  static
109  types<2>::connectivity *(&connectivity_new) (types<2>::topidx num_vertices,
110  types<2>::topidx num_trees,
111  types<2>::topidx num_corners,
112  types<2>::topidx num_vtt);
113 
114  static
115  void (&connectivity_join_faces) (types<2>::connectivity *conn,
116  types<2>::topidx tree_left,
117  types<2>::topidx tree_right,
118  int face_left,
119  int face_right,
120  int orientation);
121 
122 
123 
124  static
125  void (&connectivity_destroy) (p4est_connectivity_t *connectivity);
126 
127  static
128  types<2>::forest *(&new_forest) (MPI_Comm mpicomm,
129  types<2>::connectivity *connectivity,
130  types<2>::locidx min_quadrants,
131  int min_level,
132  int fill_uniform,
133  size_t data_size,
134  p4est_init_t init_fn,
135  void *user_pointer);
136 
137  static
138  void (&destroy) (types<2>::forest *p4est);
139 
140  static
141  void (&refine) (types<2>::forest *p4est,
142  int refine_recursive,
143  p4est_refine_t refine_fn,
144  p4est_init_t init_fn);
145 
146  static
147  void (&coarsen) (types<2>::forest *p4est,
148  int coarsen_recursive,
149  p4est_coarsen_t coarsen_fn,
150  p4est_init_t init_fn);
151 
152  static
153  void (&balance) (types<2>::forest *p4est,
155  p4est_init_t init_fn);
156 
157 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
158  static
159  p4est_gloidx_t (&partition) (types<2>::forest *p4est,
160  int partition_for_coarsening,
161  p4est_weight_t weight_fn);
162 #else
163  static
164  void (&partition) (types<2>::forest *p4est,
165  int partition_for_coarsening,
166  p4est_weight_t weight_fn);
167 #endif
168 
169  static
170  void (&save) (const char *filename,
171  types<2>::forest *p4est,
172  int save_data);
173 
174 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
175  static
176  types<2>::forest *(&load_ext) (const char *filename,
177  MPI_Comm mpicomm,
178  size_t data_size,
179  int load_data,
180  int autopartition,
181  int broadcasthead,
182  void *user_pointer,
183  types<2>::connectivity **p4est);
184 #else
185  static
186  types<2>::forest *(&load) (const char *filename,
187  MPI_Comm mpicomm,
188  size_t data_size,
189  int load_data,
190  void *user_pointer,
191  types<2>::connectivity **p4est);
192 #endif
193 
194 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
195  static
196  int (&connectivity_save) (const char *filename,
197  types<2>::connectivity *connectivity);
198 #else
199  static
200  void (&connectivity_save) (const char *filename,
201  types<2>::connectivity *connectivity);
202 #endif
203 
204  static
205  int (&connectivity_is_valid) (types<2>::connectivity *connectivity);
206 
207 #if DEAL_II_P4EST_VERSION_GTE(1,0,0,0)
208  static
209  types<2>::connectivity *(&connectivity_load) (const char *filename,
210  size_t *length);
211 #elif DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
212  static
213  types<2>::connectivity *(&connectivity_load) (const char *filename,
214  long unsigned *length);
215 #else
216  static
217  types<2>::connectivity *(&connectivity_load) (const char *filename,
218  long *length);
219 #endif
220 
221  static
222  unsigned int (&checksum) (types<2>::forest *p4est);
223 
224  static
225  void (&vtk_write_file) (types<2>::forest *p4est,
226  p4est_geometry_t *,
227  const char *baseName);
228 
229  static
230  types<2>::ghost *(&ghost_new) (types<2>::forest *p4est,
231  types<2>::balance_type btype);
232 
233  static
234  void (&ghost_destroy) (types<2>::ghost *ghost);
235 
236  static
237  void (&reset_data) (types<2>::forest *p4est,
238  size_t data_size,
239  p4est_init_t init_fn,
240  void *user_pointer);
241 
242  static
243  size_t (&forest_memory_used) (types<2>::forest *p4est);
244 
245  static
246  size_t (&connectivity_memory_used) (types<2>::connectivity *p4est);
247 
248  static const unsigned max_level;
249  };
250 
251  int (&functions<2>::quadrant_compare) (const void *v1, const void *v2)
252  = p4est_quadrant_compare;
253 
255  types<2>::quadrant c[])
256  = p4est_quadrant_childrenv;
257 
259  const types<2>::quadrant *q)
260  = p4est_quadrant_overlaps_tree;
261 
263  int level,
264  uint64_t id)
265  = p4est_quadrant_set_morton;
266 
268  const types<2>::quadrant *q2)
269  = p4est_quadrant_is_equal;
270 
272  const types<2>::quadrant *q2)
273  = p4est_quadrant_is_sibling;
274 
276  const types<2>::quadrant *q2)
277  = p4est_quadrant_is_ancestor;
278 
280  int level)
281  = p4est_quadrant_ancestor_id;
282 
284  const types<2>::locidx which_tree,
285  const types<2>::quadrant *q,
286  const int guess)
287  = p4est_comm_find_owner;
288 
290  types<2>::topidx num_trees,
291  types<2>::topidx num_corners,
292  types<2>::topidx num_vtt)
293  = p4est_connectivity_new;
294 
295 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,1)
297  types<2>::topidx tree_left,
298  types<2>::topidx tree_right,
299  int face_left,
300  int face_right,
301  int orientation)
302  = p4est_connectivity_join_faces;
303 #endif
304 
305  void (&functions<2>::connectivity_destroy) (p4est_connectivity_t *connectivity)
306  = p4est_connectivity_destroy;
307 
308  types<2>::forest *(&functions<2>::new_forest) (MPI_Comm mpicomm,
309  types<2>::connectivity *connectivity,
310  types<2>::locidx min_quadrants,
311  int min_level,
312  int fill_uniform,
313  size_t data_size,
314  p4est_init_t init_fn,
315  void *user_pointer)
316  = p4est_new_ext;
317 
318  void (&functions<2>::destroy) (types<2>::forest *p4est)
319  = p4est_destroy;
320 
321  void (&functions<2>::refine) (types<2>::forest *p4est,
322  int refine_recursive,
323  p4est_refine_t refine_fn,
324  p4est_init_t init_fn)
325  = p4est_refine;
326 
327  void (&functions<2>::coarsen) (types<2>::forest *p4est,
328  int coarsen_recursive,
329  p4est_coarsen_t coarsen_fn,
330  p4est_init_t init_fn)
331  = p4est_coarsen;
332 
333  void (&functions<2>::balance) (types<2>::forest *p4est,
335  p4est_init_t init_fn)
336  = p4est_balance;
337 
338 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
339  p4est_gloidx_t (&functions<2>::partition) (types<2>::forest *p4est,
340  int partition_for_coarsening,
341  p4est_weight_t weight_fn)
342  = p4est_partition_ext;
343 #else
345  int partition_for_coarsening,
346  p4est_weight_t weight_fn)
347  = p4est_partition_ext;
348 #endif
349 
350  void (&functions<2>::save) (const char *filename,
351  types<2>::forest *p4est,
352  int save_data)
353  = p4est_save;
354 
355 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
357  (&functions<2>::load_ext) (const char *filename,
358  MPI_Comm mpicomm,
359  std::size_t data_size,
360  int load_data,
361  int autopartition,
362  int broadcasthead,
363  void *user_pointer,
364  types<2>::connectivity **p4est)
365  = p4est_load_ext;
366 #else
368  (&functions<2>::load) (const char *filename,
369  MPI_Comm mpicomm,
370  std::size_t data_size,
371  int load_data,
372  void *user_pointer,
373  types<2>::connectivity **p4est)
374  = p4est_load;
375 #endif
376 
377 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
378  int (&functions<2>::connectivity_save) (const char *filename,
379  types<2>::connectivity *connectivity)
380  = p4est_connectivity_save;
381 #else
382  void (&functions<2>::connectivity_save) (const char *filename,
383  types<2>::connectivity *connectivity)
384  = p4est_connectivity_save;
385 #endif
386 
388  *connectivity)
389  = p4est_connectivity_is_valid;
390 
391 #if DEAL_II_P4EST_VERSION_GTE(1,0,0,0)
393  (&functions<2>::connectivity_load) (const char *filename,
394  size_t *length)
395  = p4est_connectivity_load;
396 #elif DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
398  (&functions<2>::connectivity_load) (const char *filename,
399  long unsigned *length)
400  = p4est_connectivity_load;
401 #else
403  (&functions<2>::connectivity_load) (const char *filename,
404  long *length)
405  = p4est_connectivity_load;
406 #endif
407 
408  unsigned int (&functions<2>::checksum) (types<2>::forest *p4est)
409  = p4est_checksum;
410 
412  p4est_geometry_t *,
413  const char *baseName)
414  = p4est_vtk_write_file;
415 
418  = p4est_ghost_new;
419 
421  = p4est_ghost_destroy;
422 
424  size_t data_size,
425  p4est_init_t init_fn,
426  void *user_pointer)
427  = p4est_reset_data;
428 
430  = p4est_memory_used;
431 
433  = p4est_connectivity_memory_used;
434 
435  const unsigned int functions<2>::max_level = P4EST_MAXLEVEL;
436 
437  template <> struct functions<3>
438  {
439  static
440  int (&quadrant_compare) (const void *v1, const void *v2);
441 
442  static
443  void (&quadrant_childrenv) (const types<3>::quadrant *q,
444  types<3>::quadrant c[]);
445 
446  static
447  int (&quadrant_overlaps_tree) (types<3>::tree *tree,
448  const types<3>::quadrant *q);
449 
450  static
451  void (&quadrant_set_morton) (types<3>::quadrant *quadrant,
452  int level,
453  uint64_t id);
454 
455  static
456  int (&quadrant_is_equal) (const types<3>::quadrant *q1,
457  const types<3>::quadrant *q2);
458 
459  static
460  int (&quadrant_is_sibling) (const types<3>::quadrant *q1,
461  const types<3>::quadrant *q2);
462 
463  static
464  int (&quadrant_is_ancestor) (const types<3>::quadrant *q1,
465  const types<3>::quadrant *q2);
466 
467  static
468  int (&quadrant_ancestor_id) (const types<3>::quadrant *q,
469  int level);
470 
471  static
472  int (&comm_find_owner) (types<3>::forest *p4est,
473  const types<3>::locidx which_tree,
474  const types<3>::quadrant *q,
475  const int guess);
476 
477  static
478  types<3>::connectivity *(&connectivity_new) (types<3>::topidx num_vertices,
479  types<3>::topidx num_trees,
480  types<3>::topidx num_edges,
481  types<3>::topidx num_ett,
482  types<3>::topidx num_corners,
483  types<3>::topidx num_ctt);
484 
485  static
486  void (&connectivity_join_faces) (types<3>::connectivity *conn,
487  types<3>::topidx tree_left,
488  types<3>::topidx tree_right,
489  int face_left,
490  int face_right,
491  int orientation);
492 
493  static
494  void (&connectivity_destroy) (p8est_connectivity_t *connectivity);
495 
496  static
497  types<3>::forest *(&new_forest) (MPI_Comm mpicomm,
498  types<3>::connectivity *connectivity,
499  types<3>::locidx min_quadrants,
500  int min_level,
501  int fill_uniform,
502  size_t data_size,
503  p8est_init_t init_fn,
504  void *user_pointer);
505 
506  static
507  void (&destroy) (types<3>::forest *p8est);
508 
509  static
510  void (&refine) (types<3>::forest *p8est,
511  int refine_recursive,
512  p8est_refine_t refine_fn,
513  p8est_init_t init_fn);
514 
515  static
516  void (&coarsen) (types<3>::forest *p8est,
517  int coarsen_recursive,
518  p8est_coarsen_t coarsen_fn,
519  p8est_init_t init_fn);
520 
521  static
522  void (&balance) (types<3>::forest *p8est,
524  p8est_init_t init_fn);
525 
526 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
527  static
528  p4est_gloidx_t (&partition) (types<3>::forest *p8est,
529  int partition_for_coarsening,
530  p8est_weight_t weight_fn);
531 #else
532  static
533  void (&partition) (types<3>::forest *p8est,
534  int partition_for_coarsening,
535  p8est_weight_t weight_fn);
536 #endif
537 
538  static
539  void (&save) (const char *filename,
540  types<3>::forest *p4est,
541  int save_data);
542 
543 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
544  static
545  types<3>::forest *(&load_ext) (const char *filename,
546  MPI_Comm mpicomm,
547  std::size_t data_size,
548  int load_data,
549  int autopartition,
550  int broadcasthead,
551  void *user_pointer,
552  types<3>::connectivity **p4est);
553 #else
554  static
555  types<3>::forest *(&load) (const char *filename,
556  MPI_Comm mpicomm,
557  std::size_t data_size,
558  int load_data,
559  void *user_pointer,
560  types<3>::connectivity **p4est);
561 #endif
562 
563 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
564  static
565  int (&connectivity_save) (const char *filename,
566  types<3>::connectivity *connectivity);
567 #else
568  static
569  void (&connectivity_save) (const char *filename,
570  types<3>::connectivity *connectivity);
571 #endif
572 
573  static
574  int (&connectivity_is_valid) (types<3>::connectivity *connectivity);
575 
576 #if DEAL_II_P4EST_VERSION_GTE(1,0,0,0)
577  static
578  types<3>::connectivity *(&connectivity_load) (const char *filename,
579  size_t *length);
580 #elif DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
581  static
582  types<3>::connectivity *(&connectivity_load) (const char *filename,
583  long unsigned *length);
584 #else
585  static
586  types<3>::connectivity *(&connectivity_load) (const char *filename,
587  long *length);
588 #endif
589 
590  static
591  unsigned int (&checksum) (types<3>::forest *p8est);
592 
593  static
594  void (&vtk_write_file) (types<3>::forest *p8est,
595  p8est_geometry_t *,
596  const char *baseName);
597 
598  static
599  types<3>::ghost *(&ghost_new) (types<3>::forest *p4est,
600  types<3>::balance_type btype);
601 
602  static
603  void (&ghost_destroy) (types<3>::ghost *ghost);
604 
605  static
606  void (&reset_data) (types<3>::forest *p4est,
607  size_t data_size,
608  p8est_init_t init_fn,
609  void *user_pointer);
610 
611  static
612  size_t (&forest_memory_used) (types<3>::forest *p4est);
613 
614  static
615  size_t (&connectivity_memory_used) (types<3>::connectivity *p4est);
616 
617  static const unsigned max_level;
618  };
619 
620 
621  int (&functions<3>::quadrant_compare) (const void *v1, const void *v2)
622  = p8est_quadrant_compare;
623 
625  types<3>::quadrant c[])
626  = p8est_quadrant_childrenv;
627 
629  const types<3>::quadrant *q)
630  = p8est_quadrant_overlaps_tree;
631 
633  int level,
634  uint64_t id)
635  = p8est_quadrant_set_morton;
636 
638  const types<3>::quadrant *q2)
639  = p8est_quadrant_is_equal;
640 
642  const types<3>::quadrant *q2)
643  = p8est_quadrant_is_sibling;
644 
646  const types<3>::quadrant *q2)
647  = p8est_quadrant_is_ancestor;
648 
650  int level)
651  = p8est_quadrant_ancestor_id;
652 
654  const types<3>::locidx which_tree,
655  const types<3>::quadrant *q,
656  const int guess)
657  = p8est_comm_find_owner;
658 
660  types<3>::topidx num_trees,
661  types<3>::topidx num_edges,
662  types<3>::topidx num_ett,
663  types<3>::topidx num_corners,
664  types<3>::topidx num_ctt)
665  = p8est_connectivity_new;
666 
667  void (&functions<3>::connectivity_destroy) (p8est_connectivity_t *connectivity)
668  = p8est_connectivity_destroy;
669 
670 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,1)
672  types<3>::topidx tree_left,
673  types<3>::topidx tree_right,
674  int face_left,
675  int face_right,
676  int orientation)
677  = p8est_connectivity_join_faces;
678 #endif
679 
680  types<3>::forest *(&functions<3>::new_forest) (MPI_Comm mpicomm,
681  types<3>::connectivity *connectivity,
682  types<3>::locidx min_quadrants,
683  int min_level,
684  int fill_uniform,
685  size_t data_size,
686  p8est_init_t init_fn,
687  void *user_pointer)
688  = p8est_new_ext;
689 
690  void (&functions<3>::destroy) (types<3>::forest *p8est)
691  = p8est_destroy;
692 
693  void (&functions<3>::refine) (types<3>::forest *p8est,
694  int refine_recursive,
695  p8est_refine_t refine_fn,
696  p8est_init_t init_fn)
697  = p8est_refine;
698 
699  void (&functions<3>::coarsen) (types<3>::forest *p8est,
700  int coarsen_recursive,
701  p8est_coarsen_t coarsen_fn,
702  p8est_init_t init_fn)
703  = p8est_coarsen;
704 
705  void (&functions<3>::balance) (types<3>::forest *p8est,
707  p8est_init_t init_fn)
708  = p8est_balance;
709 
710 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
711  p4est_gloidx_t (&functions<3>::partition) (types<3>::forest *p8est,
712  int partition_for_coarsening,
713  p8est_weight_t weight_fn)
714  = p8est_partition_ext;
715 #else
717  int partition_for_coarsening,
718  p8est_weight_t weight_fn)
719  = p8est_partition_ext;
720 #endif
721 
722  void (&functions<3>::save) (const char *filename,
723  types<3>::forest *p4est,
724  int save_data)
725  = p8est_save;
726 
727 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
729  (&functions<3>::load_ext) (const char *filename,
730  MPI_Comm mpicomm,
731  std::size_t data_size,
732  int load_data,
733  int autopartition,
734  int broadcasthead,
735  void *user_pointer,
736  types<3>::connectivity **p4est)
737  = p8est_load_ext;
738 #else
740  (&functions<3>::load) (const char *filename,
741  MPI_Comm mpicomm,
742  std::size_t data_size,
743  int load_data,
744  void *user_pointer,
745  types<3>::connectivity **p4est)
746  = p8est_load;
747 #endif
748 
749 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
750  int (&functions<3>::connectivity_save) (const char *filename,
751  types<3>::connectivity *connectivity)
752  = p8est_connectivity_save;
753 #else
754  void (&functions<3>::connectivity_save) (const char *filename,
755  types<3>::connectivity *connectivity)
756  = p8est_connectivity_save;
757 #endif
758 
760  *connectivity)
761  = p8est_connectivity_is_valid;
762 
763 #if DEAL_II_P4EST_VERSION_GTE(1,0,0,0)
765  (&functions<3>::connectivity_load) (const char *filename,
766  size_t *length)
767  = p8est_connectivity_load;
768 #elif DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
770  (&functions<3>::connectivity_load) (const char *filename,
771  long unsigned *length)
772  = p8est_connectivity_load;
773 #else
775  (&functions<3>::connectivity_load) (const char *filename,
776  long *length)
777  = p8est_connectivity_load;
778 #endif
779 
780  unsigned int (&functions<3>::checksum) (types<3>::forest *p8est)
781  = p8est_checksum;
782 
784  p8est_geometry_t *,
785  const char *baseName)
786  = p8est_vtk_write_file;
787 
790  = p8est_ghost_new;
791 
793  = p8est_ghost_destroy;
794 
796  size_t data_size,
797  p8est_init_t init_fn,
798  void *user_pointer)
799  = p8est_reset_data;
800 
802  = p8est_memory_used;
803 
805  = p8est_connectivity_memory_used;
806 
807  const unsigned int functions<3>::max_level = P8EST_MAXLEVEL;
808 
809 
810  template <int dim>
811  void
812  init_quadrant_children
813  (const typename types<dim>::quadrant &p4est_cell,
815  {
816 
817  for (unsigned int c=0;
818  c<GeometryInfo<dim>::max_children_per_cell; ++c)
819  switch (dim)
820  {
821  case 2:
822  P4EST_QUADRANT_INIT(&p4est_children[c]);
823  break;
824  case 3:
825  P8EST_QUADRANT_INIT(&p4est_children[c]);
826  break;
827  default:
828  Assert (false, ExcNotImplemented());
829  }
830 
831 
833  p4est_children);
834 
835  }
836 
837 
838  template <int dim>
839  void
840  init_coarse_quadrant(typename types<dim>::quadrant &quad)
841  {
842  switch (dim)
843  {
844  case 2:
845  P4EST_QUADRANT_INIT(&quad);
846  break;
847  case 3:
848  P8EST_QUADRANT_INIT(&quad);
849  break;
850  default:
851  Assert (false, ExcNotImplemented());
852  }
854  /*level=*/0,
855  /*index=*/0);
856  }
857 
858 
859  template <int dim>
860  bool
861  quadrant_is_equal (const typename types<dim>::quadrant &q1,
862  const typename types<dim>::quadrant &q2)
863  {
864  return functions<dim>::quadrant_is_equal(&q1, &q2);
865  }
866 
867 
868 
869  template <int dim>
870  bool
871  quadrant_is_ancestor (const typename types<dim>::quadrant &q1,
872  const typename types<dim>::quadrant &q2)
873  {
874  return functions<dim>::quadrant_is_ancestor(&q1, &q2);
875  }
876 
882  template <int dim> struct iter;
883 
884  template <> struct iter<2>
885  {
886  typedef p4est_iter_corner_info_t corner_info;
887  typedef p4est_iter_corner_side_t corner_side;
888  typedef p4est_iter_corner_t corner_iter;
889  typedef p4est_iter_face_info_t face_info;
890  typedef p4est_iter_face_side_t face_side;
891  typedef p4est_iter_face_t face_iter;
892  };
893 
894  template <> struct iter<3>
895  {
896  typedef p8est_iter_corner_info_t corner_info;
897  typedef p8est_iter_corner_side_t corner_side;
898  typedef p8est_iter_corner_t corner_iter;
899  typedef p8est_iter_edge_info_t edge_info;
900  typedef p8est_iter_edge_side_t edge_side;
901  typedef p8est_iter_edge_t edge_iter;
902  typedef p8est_iter_face_info_t face_info;
903  typedef p8est_iter_face_side_t face_side;
904  typedef p8est_iter_face_t face_iter;
905  };
906 
907  }
908 }
909 
910 
911 namespace
912 {
913  template <int dim, int spacedim>
914  void
915  get_vertex_to_cell_mappings (const Triangulation<dim,spacedim> &triangulation,
916  std::vector<unsigned int> &vertex_touch_count,
917  std::vector<std::list<
918  std::pair<typename Triangulation<dim,spacedim>::active_cell_iterator,unsigned int> > >
919  &vertex_to_cell)
920  {
921  vertex_touch_count.resize (triangulation.n_vertices());
922  vertex_to_cell.resize (triangulation.n_vertices());
923 
925  cell = triangulation.begin_active();
926  cell != triangulation.end(); ++cell)
927  for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
928  {
929  ++vertex_touch_count[cell->vertex_index(v)];
930  vertex_to_cell[cell->vertex_index(v)]
931  .push_back (std::make_pair (cell, v));
932  }
933  }
934 
935 
936 
937  template <int dim, int spacedim>
938  void
939  get_edge_to_cell_mappings (const Triangulation<dim,spacedim> &triangulation,
940  std::vector<unsigned int> &edge_touch_count,
941  std::vector<std::list<
942  std::pair<typename Triangulation<dim,spacedim>::active_cell_iterator,unsigned int> > >
943  &edge_to_cell)
944  {
945  Assert (triangulation.n_levels() == 1, ExcInternalError());
946 
947  edge_touch_count.resize (triangulation.n_active_lines());
948  edge_to_cell.resize (triangulation.n_active_lines());
949 
951  cell = triangulation.begin_active();
952  cell != triangulation.end(); ++cell)
953  for (unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++l)
954  {
955  ++edge_touch_count[cell->line(l)->index()];
956  edge_to_cell[cell->line(l)->index()]
957  .push_back (std::make_pair (cell, l));
958  }
959  }
960 
961 
962 
967  template <int dim, int spacedim>
968  void
969  set_vertex_and_cell_info (const Triangulation<dim,spacedim> &triangulation,
970  const std::vector<unsigned int> &vertex_touch_count,
971  const std::vector<std::list<
972  std::pair<typename Triangulation<dim,spacedim>::active_cell_iterator,unsigned int> > >
973  &vertex_to_cell,
974  const std::vector<types::global_dof_index> &coarse_cell_to_p4est_tree_permutation,
975  const bool set_vertex_info,
976  typename internal::p4est::types<dim>::connectivity *connectivity)
977  {
978  // copy the vertices into the connectivity structure. the triangulation
979  // exports the array of vertices, but some of the entries are sometimes
980  // unused; this shouldn't be the case for a newly created triangulation,
981  // but make sure
982  //
983  // note that p4est stores coordinates as a triplet of values even in 2d
984  Assert (triangulation.get_used_vertices().size() ==
985  triangulation.get_vertices().size(),
986  ExcInternalError());
987  Assert (std::find (triangulation.get_used_vertices().begin(),
988  triangulation.get_used_vertices().end(),
989  false)
990  == triangulation.get_used_vertices().end(),
991  ExcInternalError());
992  if (set_vertex_info == true)
993  for (unsigned int v=0; v<triangulation.n_vertices(); ++v)
994  {
995  connectivity->vertices[3*v ] = triangulation.get_vertices()[v][0];
996  connectivity->vertices[3*v+1] = triangulation.get_vertices()[v][1];
997  connectivity->vertices[3*v+2] = (spacedim == 2 ?
998  0
999  :
1000  triangulation.get_vertices()[v][2]);
1001  }
1002 
1003  // next store the tree_to_vertex indices (each tree is here only a single
1004  // cell in the coarse mesh). p4est requires vertex numbering in clockwise
1005  // orientation
1006  //
1007  // while we're at it, also copy the neighborship information between cells
1009  cell = triangulation.begin_active(),
1010  endc = triangulation.end();
1011  for (; cell != endc; ++cell)
1012  {
1013  const unsigned int
1014  index = coarse_cell_to_p4est_tree_permutation[cell->index()];
1015 
1016  for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
1017  {
1018  if (set_vertex_info == true)
1019  connectivity->tree_to_vertex[index*GeometryInfo<dim>::vertices_per_cell+v] = cell->vertex_index(v);
1020  connectivity->tree_to_corner[index*GeometryInfo<dim>::vertices_per_cell+v] = cell->vertex_index(v);
1021  }
1022 
1023  // neighborship information. if a cell is at a boundary, then enter
1024  // the index of the cell itself here
1025  for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
1026  if (cell->face(f)->at_boundary() == false)
1027  connectivity->tree_to_tree[index*GeometryInfo<dim>::faces_per_cell + f]
1028  = coarse_cell_to_p4est_tree_permutation[cell->neighbor(f)->index()];
1029  else
1030  connectivity->tree_to_tree[index*GeometryInfo<dim>::faces_per_cell + f]
1031  = coarse_cell_to_p4est_tree_permutation[cell->index()];
1032 
1033  // fill tree_to_face, which is essentially neighbor_to_neighbor;
1034  // however, we have to remap the resulting face number as well
1035  for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
1036  if (cell->face(f)->at_boundary() == false)
1037  {
1038  switch (dim)
1039  {
1040  case 2:
1041  {
1042  connectivity->tree_to_face[index*GeometryInfo<dim>::faces_per_cell + f]
1043  = cell->neighbor_of_neighbor (f);
1044  break;
1045  }
1046 
1047  case 3:
1048  {
1049  /*
1050  * The values for tree_to_face are in 0..23 where ttf % 6
1051  * gives the face number and ttf / 4 the face orientation
1052  * code. The orientation is determined as follows. Let
1053  * my_face and other_face be the two face numbers of the
1054  * connecting trees in 0..5. Then the first face vertex of
1055  * the lower of my_face and other_face connects to a face
1056  * vertex numbered 0..3 in the higher of my_face and
1057  * other_face. The face orientation is defined as this
1058  * number. If my_face == other_face, treating either of
1059  * both faces as the lower one leads to the same result.
1060  */
1061 
1062  connectivity->tree_to_face[index*6 + f]
1063  = cell->neighbor_of_neighbor (f);
1064 
1065  unsigned int face_idx_list[2] =
1066  {f, cell->neighbor_of_neighbor (f)};
1068  cell_list[2] = {cell, cell->neighbor(f)};
1069  unsigned int smaller_idx = 0;
1070 
1071  if (f>cell->neighbor_of_neighbor (f))
1072  smaller_idx = 1;
1073 
1074  unsigned int larger_idx = (smaller_idx+1) % 2;
1075  //smaller = *_list[smaller_idx]
1076  //larger = *_list[larger_idx]
1077 
1078  unsigned int v = 0;
1079 
1080  // global vertex index of vertex 0 on face of cell with
1081  // smaller local face index
1082  unsigned int g_idx =
1083  cell_list[smaller_idx]->vertex_index(
1085  face_idx_list[smaller_idx],
1086  0,
1087  cell_list[smaller_idx]->face_orientation(face_idx_list[smaller_idx]),
1088  cell_list[smaller_idx]->face_flip(face_idx_list[smaller_idx]),
1089  cell_list[smaller_idx]->face_rotation(face_idx_list[smaller_idx]))
1090  );
1091 
1092  // loop over vertices on face from other cell and compare
1093  // global vertex numbers
1094  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_face; ++i)
1095  {
1096  unsigned int idx
1097  =
1098  cell_list[larger_idx]->vertex_index(
1100  face_idx_list[larger_idx],
1101  i)
1102  );
1103 
1104  if (idx==g_idx)
1105  {
1106  v = i;
1107  break;
1108  }
1109  }
1110 
1111  connectivity->tree_to_face[index*6 + f] += 6*v;
1112  break;
1113  }
1114 
1115  default:
1116  Assert (false, ExcNotImplemented());
1117  }
1118  }
1119  else
1120  connectivity->tree_to_face[index*GeometryInfo<dim>::faces_per_cell + f] = f;
1121  }
1122 
1123  // now fill the vertex information
1124  connectivity->ctt_offset[0] = 0;
1125  std::partial_sum (vertex_touch_count.begin(),
1126  vertex_touch_count.end(),
1127  &connectivity->ctt_offset[1]);
1128 
1130  num_vtt = std::accumulate (vertex_touch_count.begin(),
1131  vertex_touch_count.end(),
1132  0);
1133  (void)num_vtt;
1134  Assert (connectivity->ctt_offset[triangulation.n_vertices()] ==
1135  num_vtt,
1136  ExcInternalError());
1137 
1138  for (unsigned int v=0; v<triangulation.n_vertices(); ++v)
1139  {
1140  Assert (vertex_to_cell[v].size() == vertex_touch_count[v],
1141  ExcInternalError());
1142 
1143  typename std::list<std::pair
1145  unsigned int> >::const_iterator
1146  p = vertex_to_cell[v].begin();
1147  for (unsigned int c=0; c<vertex_touch_count[v]; ++c, ++p)
1148  {
1149  connectivity->corner_to_tree[connectivity->ctt_offset[v]+c]
1150  = coarse_cell_to_p4est_tree_permutation[p->first->index()];
1151  connectivity->corner_to_corner[connectivity->ctt_offset[v]+c]
1152  = p->second;
1153  }
1154  }
1155  }
1156 
1157 
1158 
1159  template <int dim, int spacedim>
1160  bool
1161  tree_exists_locally (const typename internal::p4est::types<dim>::forest *parallel_forest,
1162  const typename internal::p4est::types<dim>::topidx coarse_grid_cell)
1163  {
1164  Assert (coarse_grid_cell < parallel_forest->connectivity->num_trees,
1165  ExcInternalError());
1166  return ((coarse_grid_cell >= parallel_forest->first_local_tree)
1167  &&
1168  (coarse_grid_cell <= parallel_forest->last_local_tree));
1169  }
1170 
1171 
1172  template <int dim, int spacedim>
1173  void
1174  delete_all_children_and_self (const typename Triangulation<dim,spacedim>::cell_iterator &cell)
1175  {
1176  if (cell->has_children())
1177  for (unsigned int c=0; c<cell->n_children(); ++c)
1178  delete_all_children_and_self<dim,spacedim> (cell->child(c));
1179  else
1180  cell->set_coarsen_flag ();
1181  }
1182 
1183 
1184 
1185  template <int dim, int spacedim>
1186  void
1187  delete_all_children (const typename Triangulation<dim,spacedim>::cell_iterator &cell)
1188  {
1189  if (cell->has_children())
1190  for (unsigned int c=0; c<cell->n_children(); ++c)
1191  delete_all_children_and_self<dim,spacedim> (cell->child(c));
1192  }
1193 
1194 
1195  template <int dim, int spacedim>
1196  void
1197  determine_level_subdomain_id_recursively (const typename internal::p4est::types<dim>::tree &tree,
1198  const typename internal::p4est::types<dim>::locidx &tree_index,
1199  const typename Triangulation<dim,spacedim>::cell_iterator &dealii_cell,
1200  const typename internal::p4est::types<dim>::quadrant &p4est_cell,
1201  typename internal::p4est::types<dim>::forest &forest,
1202  const types::subdomain_id my_subdomain,
1203  const std::vector<std::vector<bool> > &marked_vertices)
1204  {
1205  if (dealii_cell->level_subdomain_id()==numbers::artificial_subdomain_id)
1206  {
1207  //important: only assign the level_subdomain_id if it is a ghost cell
1208  // even though we could fill in all.
1209  bool used = false;
1210  for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
1211  {
1212  if (marked_vertices[dealii_cell->level()][dealii_cell->vertex_index(v)])
1213  {
1214  used = true;
1215  break;
1216  }
1217  }
1218 
1219  if (used)
1220  {
1222  tree_index,
1223  &p4est_cell,
1224  my_subdomain);
1225  Assert((owner!=-2) && (owner!=-1), ExcMessage("p4est should know the owner."));
1226  dealii_cell->set_level_subdomain_id(owner);
1227  }
1228 
1229  }
1230 
1231  if (dealii_cell->has_children ())
1232  {
1234  p4est_child[GeometryInfo<dim>::max_children_per_cell];
1235  for (unsigned int c=0;
1236  c<GeometryInfo<dim>::max_children_per_cell; ++c)
1237  switch (dim)
1238  {
1239  case 2:
1240  P4EST_QUADRANT_INIT(&p4est_child[c]);
1241  break;
1242  case 3:
1243  P8EST_QUADRANT_INIT(&p4est_child[c]);
1244  break;
1245  default:
1246  Assert (false, ExcNotImplemented());
1247  }
1248 
1249 
1251  quadrant_childrenv (&p4est_cell,
1252  p4est_child);
1253 
1254  for (unsigned int c=0;
1255  c<GeometryInfo<dim>::max_children_per_cell; ++c)
1256  {
1257  determine_level_subdomain_id_recursively <dim,spacedim> (tree,tree_index,
1258  dealii_cell->child(c),
1259  p4est_child[c],
1260  forest,
1261  my_subdomain,
1262  marked_vertices);
1263  }
1264  }
1265  }
1266 
1267 
1268  template <int dim, int spacedim>
1269  void
1270  match_tree_recursively (const typename internal::p4est::types<dim>::tree &tree,
1271  const typename Triangulation<dim,spacedim>::cell_iterator &dealii_cell,
1272  const typename internal::p4est::types<dim>::quadrant &p4est_cell,
1273  const typename internal::p4est::types<dim>::forest &forest,
1274  const types::subdomain_id my_subdomain)
1275  {
1276  // check if this cell exists in the local p4est cell
1277  if (sc_array_bsearch(const_cast<sc_array_t *>(&tree.quadrants),
1278  &p4est_cell,
1280  != -1)
1281  {
1282  // yes, cell found in local part of p4est
1283  delete_all_children<dim,spacedim> (dealii_cell);
1284  if (!dealii_cell->has_children())
1285  dealii_cell->set_subdomain_id(my_subdomain);
1286  }
1287  else
1288  {
1289  // no, cell not found in local part of p4est. this means that the
1290  // local part is more refined than the current cell. if this cell has
1291  // no children of its own, we need to refine it, and if it does
1292  // already have children then loop over all children and see if they
1293  // are locally available as well
1294  if (dealii_cell->has_children () == false)
1295  dealii_cell->set_refine_flag ();
1296  else
1297  {
1299  p4est_child[GeometryInfo<dim>::max_children_per_cell];
1300  for (unsigned int c=0;
1301  c<GeometryInfo<dim>::max_children_per_cell; ++c)
1302  switch (dim)
1303  {
1304  case 2:
1305  P4EST_QUADRANT_INIT(&p4est_child[c]);
1306  break;
1307  case 3:
1308  P8EST_QUADRANT_INIT(&p4est_child[c]);
1309  break;
1310  default:
1311  Assert (false, ExcNotImplemented());
1312  }
1313 
1314 
1316  quadrant_childrenv (&p4est_cell,
1317  p4est_child);
1318 
1319  for (unsigned int c=0;
1320  c<GeometryInfo<dim>::max_children_per_cell; ++c)
1322  quadrant_overlaps_tree (const_cast<typename internal::p4est::types<dim>::tree *>(&tree),
1323  &p4est_child[c])
1324  == false)
1325  {
1326  // no, this child is locally not available in the p4est.
1327  // delete all its children but, because this may not be
1328  // successfull, make sure to mark all children recursively
1329  // as not local.
1330  delete_all_children<dim,spacedim> (dealii_cell->child(c));
1331  dealii_cell->child(c)
1332  ->recursively_set_subdomain_id(numbers::artificial_subdomain_id);
1333  }
1334  else
1335  {
1336  // at least some part of the tree rooted in this child is
1337  // locally available
1338  match_tree_recursively<dim,spacedim> (tree,
1339  dealii_cell->child(c),
1340  p4est_child[c],
1341  forest,
1342  my_subdomain);
1343  }
1344  }
1345  }
1346  }
1347 
1348 
1349  template <int dim, int spacedim>
1350  void
1351  match_quadrant (const ::Triangulation<dim,spacedim> *tria,
1352  unsigned int dealii_index,
1353  typename internal::p4est::types<dim>::quadrant &ghost_quadrant,
1354  unsigned int ghost_owner)
1355  {
1356  int i, child_id;
1357  int l = ghost_quadrant.level;
1358 
1359  for (i = 0; i < l; i++)
1360  {
1361  typename Triangulation<dim,spacedim>::cell_iterator cell (tria, i, dealii_index);
1362  if (cell->has_children () == false)
1363  {
1364  cell->clear_coarsen_flag();
1365  cell->set_refine_flag ();
1366  return;
1367  }
1368 
1369  child_id = internal::p4est::functions<dim>::quadrant_ancestor_id (&ghost_quadrant, i + 1);
1370  dealii_index = cell->child_index(child_id);
1371  }
1372 
1373  typename Triangulation<dim,spacedim>::cell_iterator cell (tria, l, dealii_index);
1374  if (cell->has_children())
1375  delete_all_children<dim,spacedim> (cell);
1376  else
1377  {
1378  cell->clear_coarsen_flag();
1379  cell->set_subdomain_id(ghost_owner);
1380  }
1381  }
1382 
1383 
1384 
1385  template <int dim, int spacedim>
1386  void
1387  attach_mesh_data_recursively (const typename internal::p4est::types<dim>::tree &tree,
1388  const typename Triangulation<dim,spacedim>::cell_iterator &dealii_cell,
1389  const typename internal::p4est::types<dim>::quadrant &p4est_cell,
1390  const typename std::list<std::pair<unsigned int, typename std_cxx11::function<
1392  typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus,
1393  void *)
1394  > > > &attached_data_pack_callbacks)
1395  {
1396  typedef std::list<std::pair<unsigned int, typename std_cxx11::function<
1398  typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus,
1399  void *)
1400  > > > callback_list_t;
1401 
1402  int idx = sc_array_bsearch(const_cast<sc_array_t *>(&tree.quadrants),
1403  &p4est_cell,
1405 
1406  if (idx == -1 && (internal::p4est::functions<dim>::
1407  quadrant_overlaps_tree (const_cast<typename internal::p4est::types<dim>::tree *>(&tree),
1408  &p4est_cell)
1409  == false))
1410  return; //this quadrant and none of its children belongs to us.
1411 
1412  bool p4est_has_children = (idx == -1);
1413 
1414  if (p4est_has_children && dealii_cell->has_children())
1415  {
1416  //recurse further
1418  p4est_child[GeometryInfo<dim>::max_children_per_cell];
1419  for (unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
1420  switch (dim)
1421  {
1422  case 2:
1423  P4EST_QUADRANT_INIT(&p4est_child[c]);
1424  break;
1425  case 3:
1426  P8EST_QUADRANT_INIT(&p4est_child[c]);
1427  break;
1428  default:
1429  Assert (false, ExcNotImplemented());
1430  }
1431 
1433  quadrant_childrenv (&p4est_cell, p4est_child);
1434 
1435  for (unsigned int c=0;
1436  c<GeometryInfo<dim>::max_children_per_cell; ++c)
1437  {
1438  attach_mesh_data_recursively<dim,spacedim> (tree,
1439  dealii_cell->child(c),
1440  p4est_child[c],
1441  attached_data_pack_callbacks);
1442  }
1443  }
1444  else if (!p4est_has_children && !dealii_cell->has_children())
1445  {
1446  //this active cell didn't change
1448  q = static_cast<typename internal::p4est::types<dim>::quadrant *> (
1449  sc_array_index (const_cast<sc_array_t *>(&tree.quadrants), idx)
1450  );
1451  *static_cast<typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus *>(q->p.user_data) = parallel::distributed::Triangulation<dim,spacedim>::CELL_PERSIST;
1452 
1453  for (typename callback_list_t::const_iterator it = attached_data_pack_callbacks.begin();
1454  it != attached_data_pack_callbacks.end();
1455  ++it)
1456  {
1457  void *ptr = static_cast<char *>(q->p.user_data) + (*it).first; //add offset
1458  ((*it).second)(dealii_cell,
1460  ptr);
1461  }
1462  }
1463  else if (p4est_has_children)
1464  {
1465  //this cell got refined
1466 
1467  //attach to the first child, because we can only attach to active
1468  // quadrants
1470  p4est_child[GeometryInfo<dim>::max_children_per_cell];
1471  for (unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
1472  switch (dim)
1473  {
1474  case 2:
1475  P4EST_QUADRANT_INIT(&p4est_child[c]);
1476  break;
1477  case 3:
1478  P8EST_QUADRANT_INIT(&p4est_child[c]);
1479  break;
1480  default:
1481  Assert (false, ExcNotImplemented());
1482  }
1483 
1485  quadrant_childrenv (&p4est_cell, p4est_child);
1486  int child0_idx = sc_array_bsearch(const_cast<sc_array_t *>(&tree.quadrants),
1487  &p4est_child[0],
1489  Assert(child0_idx != -1, ExcMessage("the first child should exist as an active quadrant!"));
1490 
1492  q = static_cast<typename internal::p4est::types<dim>::quadrant *> (
1493  sc_array_index (const_cast<sc_array_t *>(&tree.quadrants), child0_idx)
1494  );
1495  *static_cast<typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus *>(q->p.user_data) = parallel::distributed::Triangulation<dim,spacedim>::CELL_REFINE;
1496 
1497  for (typename callback_list_t::const_iterator it = attached_data_pack_callbacks.begin();
1498  it != attached_data_pack_callbacks.end();
1499  ++it)
1500  {
1501  void *ptr = static_cast<char *>(q->p.user_data) + (*it).first; //add offset
1502 
1503  ((*it).second)(dealii_cell,
1505  ptr);
1506  }
1507 
1508  //mark other children as invalid, so that unpack only happens once
1509  for (unsigned int i=1; i<GeometryInfo<dim>::max_children_per_cell; ++i)
1510  {
1511  int child_idx = sc_array_bsearch(const_cast<sc_array_t *>(&tree.quadrants),
1512  &p4est_child[i],
1514  q = static_cast<typename internal::p4est::types<dim>::quadrant *> (
1515  sc_array_index (const_cast<sc_array_t *>(&tree.quadrants), child_idx)
1516  );
1517  *static_cast<typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus *>(q->p.user_data) = parallel::distributed::Triangulation<dim,spacedim>::CELL_INVALID;
1518  }
1519 
1520 
1521  }
1522  else
1523  {
1524  //its children got coarsened into this cell
1526  q = static_cast<typename internal::p4est::types<dim>::quadrant *> (
1527  sc_array_index (const_cast<sc_array_t *>(&tree.quadrants), idx)
1528  );
1529  *static_cast<typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus *>(q->p.user_data) = parallel::distributed::Triangulation<dim,spacedim>::CELL_COARSEN;
1530 
1531  for (typename callback_list_t::const_iterator it = attached_data_pack_callbacks.begin();
1532  it != attached_data_pack_callbacks.end();
1533  ++it)
1534  {
1535  void *ptr = static_cast<char *>(q->p.user_data) + (*it).first; //add offset
1536  ((*it).second)(dealii_cell,
1538  ptr);
1539  }
1540  }
1541  }
1542 
1543  template <int dim, int spacedim>
1544  void
1545  get_cell_weights_recursively (const typename internal::p4est::types<dim>::tree &tree,
1546  const typename Triangulation<dim,spacedim>::cell_iterator &dealii_cell,
1547  const typename internal::p4est::types<dim>::quadrant &p4est_cell,
1548  const typename Triangulation<dim,spacedim>::Signals &signals,
1549  std::vector<unsigned int> &weight)
1550  {
1551  const int idx = sc_array_bsearch(const_cast<sc_array_t *>(&tree.quadrants),
1552  &p4est_cell,
1554 
1555  if (idx == -1 && (internal::p4est::functions<dim>::
1556  quadrant_overlaps_tree (const_cast<typename internal::p4est::types<dim>::tree *>(&tree),
1557  &p4est_cell)
1558  == false))
1559  return; // This quadrant and none of its children belongs to us.
1560 
1561  const bool p4est_has_children = (idx == -1);
1562 
1563  if (p4est_has_children && dealii_cell->has_children())
1564  {
1565  //recurse further
1567  p4est_child[GeometryInfo<dim>::max_children_per_cell];
1568  for (unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
1569  switch (dim)
1570  {
1571  case 2:
1572  P4EST_QUADRANT_INIT(&p4est_child[c]);
1573  break;
1574  case 3:
1575  P8EST_QUADRANT_INIT(&p4est_child[c]);
1576  break;
1577  default:
1578  Assert (false, ExcNotImplemented());
1579  }
1580 
1582  quadrant_childrenv (&p4est_cell, p4est_child);
1583 
1584  for (unsigned int c=0;
1585  c<GeometryInfo<dim>::max_children_per_cell; ++c)
1586  {
1587  get_cell_weights_recursively<dim,spacedim> (tree,
1588  dealii_cell->child(c),
1589  p4est_child[c],
1590  signals,
1591  weight);
1592  }
1593  }
1594  else if (!p4est_has_children && !dealii_cell->has_children())
1595  {
1596  // This active cell didn't change
1597  weight.push_back(1000);
1598  weight.back() += signals.cell_weight(dealii_cell,
1600  }
1601  else if (p4est_has_children)
1602  {
1603  // This cell will be refined
1604  unsigned int parent_weight(1000);
1605  parent_weight += signals.cell_weight(dealii_cell,
1607 
1608  for (unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
1609  {
1610  // We assign the weight of the parent cell equally to all children
1611  weight.push_back(parent_weight);
1612  }
1613  }
1614  else
1615  {
1616  // This cell's children will be coarsened into this cell
1617  weight.push_back(1000);
1618  weight.back() += signals.cell_weight(dealii_cell,
1620  }
1621  }
1622 
1623 
1624  template <int dim, int spacedim>
1625  void
1626  post_mesh_data_recursively (const typename internal::p4est::types<dim>::tree &tree,
1627  const typename Triangulation<dim,spacedim>::cell_iterator &dealii_cell,
1628  const typename Triangulation<dim,spacedim>::cell_iterator &parent_cell,
1629  const typename internal::p4est::types<dim>::quadrant &p4est_cell,
1630  const unsigned int offset,
1631  const typename std_cxx11::function<
1632  void(typename parallel::distributed::Triangulation<dim,spacedim>::cell_iterator, typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus, void *)
1633  > &unpack_callback)
1634  {
1635  int idx = sc_array_bsearch(const_cast<sc_array_t *>(&tree.quadrants),
1636  &p4est_cell,
1638  if (idx == -1 && (internal::p4est::functions<dim>::
1639  quadrant_overlaps_tree (const_cast<typename internal::p4est::types<dim>::tree *>(&tree),
1640  &p4est_cell)
1641  == false))
1642  // this quadrant and none of its children belong to us.
1643  return;
1644 
1645 
1646  const bool p4est_has_children = (idx == -1);
1647  if (p4est_has_children)
1648  {
1649  Assert(dealii_cell->has_children(), ExcInternalError());
1650 
1651  //recurse further
1653  p4est_child[GeometryInfo<dim>::max_children_per_cell];
1654  for (unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
1655  switch (dim)
1656  {
1657  case 2:
1658  P4EST_QUADRANT_INIT(&p4est_child[c]);
1659  break;
1660  case 3:
1661  P8EST_QUADRANT_INIT(&p4est_child[c]);
1662  break;
1663  default:
1664  Assert (false, ExcNotImplemented());
1665  }
1666 
1668  quadrant_childrenv (&p4est_cell, p4est_child);
1669 
1670  for (unsigned int c=0;
1671  c<GeometryInfo<dim>::max_children_per_cell; ++c)
1672  {
1673  post_mesh_data_recursively<dim,spacedim> (tree,
1674  dealii_cell->child(c),
1675  dealii_cell,
1676  p4est_child[c],
1677  offset,
1678  unpack_callback);
1679  }
1680  }
1681  else
1682  {
1683  Assert(! dealii_cell->has_children(), ExcInternalError());
1684 
1686  q = static_cast<typename internal::p4est::types<dim>::quadrant *> (
1687  sc_array_index (const_cast<sc_array_t *>(&tree.quadrants), idx)
1688  );
1689 
1690  void *ptr = static_cast<char *>(q->p.user_data) + offset;
1691  typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus
1692  status = * static_cast<
1693  typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus *
1694  >(q->p.user_data);
1695  switch (status)
1696  {
1698  {
1699  unpack_callback(dealii_cell, status, ptr);
1700  break;
1701  }
1703  {
1704  unpack_callback(parent_cell, status, ptr);
1705  break;
1706  }
1708  {
1709  unpack_callback(dealii_cell, status, ptr);
1710  break;
1711  }
1713  {
1714  break;
1715  }
1716  default:
1717  AssertThrow (false, ExcInternalError());
1718  }
1719  }
1720  }
1721 
1722 
1723 
1729  template <int dim, int spacedim>
1730  class RefineAndCoarsenList
1731  {
1732  public:
1733  RefineAndCoarsenList (const Triangulation<dim,spacedim> &triangulation,
1734  const std::vector<types::global_dof_index> &p4est_tree_to_coarse_cell_permutation,
1735  const types::subdomain_id my_subdomain);
1736 
1745  static
1746  int
1747  refine_callback (typename internal::p4est::types<dim>::forest *forest,
1748  typename internal::p4est::types<dim>::topidx coarse_cell_index,
1749  typename internal::p4est::types<dim>::quadrant *quadrant);
1750 
1755  static
1756  int
1757  coarsen_callback (typename internal::p4est::types<dim>::forest *forest,
1758  typename internal::p4est::types<dim>::topidx coarse_cell_index,
1759  typename internal::p4est::types<dim>::quadrant *children[]);
1760 
1761  bool pointers_are_at_end () const;
1762 
1763  private:
1764  std::vector<typename internal::p4est::types<dim>::quadrant> refine_list;
1765  typename std::vector<typename internal::p4est::types<dim>::quadrant>::const_iterator current_refine_pointer;
1766 
1767  std::vector<typename internal::p4est::types<dim>::quadrant> coarsen_list;
1768  typename std::vector<typename internal::p4est::types<dim>::quadrant>::const_iterator current_coarsen_pointer;
1769 
1770  void build_lists (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
1771  const typename internal::p4est::types<dim>::quadrant &p4est_cell,
1772  const unsigned int myid);
1773  };
1774 
1775 
1776 
1777  template <int dim, int spacedim>
1778  bool
1779  RefineAndCoarsenList<dim,spacedim>::
1780  pointers_are_at_end () const
1781  {
1782  return ((current_refine_pointer == refine_list.end())
1783  &&
1784  (current_coarsen_pointer == coarsen_list.end()));
1785  }
1786 
1787 
1788 
1789  template <int dim, int spacedim>
1790  RefineAndCoarsenList<dim,spacedim>::
1791  RefineAndCoarsenList (const Triangulation<dim,spacedim> &triangulation,
1792  const std::vector<types::global_dof_index> &p4est_tree_to_coarse_cell_permutation,
1793  const types::subdomain_id my_subdomain)
1794  {
1795  // count how many flags are set and allocate that much memory
1796  unsigned int n_refine_flags = 0,
1797  n_coarsen_flags = 0;
1798  for (typename Triangulation<dim,spacedim>::active_cell_iterator
1799  cell = triangulation.begin_active();
1800  cell != triangulation.end(); ++cell)
1801  {
1802  //skip cells that are not local
1803  if (cell->subdomain_id() != my_subdomain)
1804  continue;
1805 
1806  if (cell->refine_flag_set())
1807  ++n_refine_flags;
1808  else if (cell->coarsen_flag_set())
1809  ++n_coarsen_flags;
1810  }
1811 
1812  refine_list.reserve (n_refine_flags);
1813  coarsen_list.reserve (n_coarsen_flags);
1814 
1815 
1816  // now build the lists of cells that are flagged. note that p4est will
1817  // traverse its cells in the order in which trees appear in the
1818  // forest. this order is not the same as the order of coarse cells in the
1819  // deal.II Triangulation because we have translated everything by the
1820  // coarse_cell_to_p4est_tree_permutation permutation. in order to make
1821  // sure that the output array is already in the correct order, traverse
1822  // our coarse cells in the same order in which p4est will:
1823  for (unsigned int c=0; c<triangulation.n_cells(0); ++c)
1824  {
1825  unsigned int coarse_cell_index =
1826  p4est_tree_to_coarse_cell_permutation[c];
1827 
1829  cell (&triangulation, 0, coarse_cell_index);
1830 
1831  typename internal::p4est::types<dim>::quadrant p4est_cell;
1833  quadrant_set_morton (&p4est_cell,
1834  /*level=*/0,
1835  /*index=*/0);
1836  p4est_cell.p.which_tree = c;
1837  build_lists (cell, p4est_cell, my_subdomain);
1838  }
1839 
1840 
1841  Assert(refine_list.size() == n_refine_flags,
1842  ExcInternalError());
1843  Assert(coarsen_list.size() == n_coarsen_flags,
1844  ExcInternalError());
1845 
1846  // make sure that our ordering in fact worked
1847  for (unsigned int i=1; i<refine_list.size(); ++i)
1848  Assert (refine_list[i].p.which_tree >=
1849  refine_list[i-1].p.which_tree,
1850  ExcInternalError());
1851  for (unsigned int i=1; i<coarsen_list.size(); ++i)
1852  Assert (coarsen_list[i].p.which_tree >=
1853  coarsen_list[i-1].p.which_tree,
1854  ExcInternalError());
1855 
1856  current_refine_pointer = refine_list.begin();
1857  current_coarsen_pointer = coarsen_list.begin();
1858  }
1859 
1860 
1861 
1862  template <int dim, int spacedim>
1863  void
1864  RefineAndCoarsenList<dim,spacedim>::
1865  build_lists (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
1866  const typename internal::p4est::types<dim>::quadrant &p4est_cell,
1867  const types::subdomain_id my_subdomain)
1868  {
1869  if (!cell->has_children())
1870  {
1871  if (cell->subdomain_id() == my_subdomain)
1872  {
1873  if (cell->refine_flag_set())
1874  refine_list.push_back (p4est_cell);
1875  else if (cell->coarsen_flag_set())
1876  coarsen_list.push_back (p4est_cell);
1877  }
1878  }
1879  else
1880  {
1882  p4est_child[GeometryInfo<dim>::max_children_per_cell];
1883  for (unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
1884  switch (dim)
1885  {
1886  case 2:
1887  P4EST_QUADRANT_INIT(&p4est_child[c]);
1888  break;
1889  case 3:
1890  P8EST_QUADRANT_INIT(&p4est_child[c]);
1891  break;
1892  default:
1893  Assert (false, ExcNotImplemented());
1894  }
1896  quadrant_childrenv (&p4est_cell,
1897  p4est_child);
1898  for (unsigned int c=0;
1899  c<GeometryInfo<dim>::max_children_per_cell; ++c)
1900  {
1901  p4est_child[c].p.which_tree = p4est_cell.p.which_tree;
1902  build_lists (cell->child(c),
1903  p4est_child[c],
1904  my_subdomain);
1905  }
1906  }
1907  }
1908 
1909 
1910  template <int dim, int spacedim>
1911  int
1912  RefineAndCoarsenList<dim,spacedim>::
1913  refine_callback (typename internal::p4est::types<dim>::forest *forest,
1914  typename internal::p4est::types<dim>::topidx coarse_cell_index,
1915  typename internal::p4est::types<dim>::quadrant *quadrant)
1916  {
1917  RefineAndCoarsenList<dim,spacedim> *this_object
1918  = reinterpret_cast<RefineAndCoarsenList<dim,spacedim>*>(forest->user_pointer);
1919 
1920  // if there are no more cells in our list the current cell can't be
1921  // flagged for refinement
1922  if (this_object->current_refine_pointer == this_object->refine_list.end())
1923  return false;
1924 
1925  Assert (coarse_cell_index <=
1926  this_object->current_refine_pointer->p.which_tree,
1927  ExcInternalError());
1928 
1929  // if p4est hasn't yet reached the tree of the next flagged cell the
1930  // current cell can't be flagged for refinement
1931  if (coarse_cell_index <
1932  this_object->current_refine_pointer->p.which_tree)
1933  return false;
1934 
1935  // now we're in the right tree in the forest
1936  Assert (coarse_cell_index <=
1937  this_object->current_refine_pointer->p.which_tree,
1938  ExcInternalError());
1939 
1940  // make sure that the p4est loop over cells hasn't gotten ahead of our own
1941  // pointer
1943  quadrant_compare (quadrant,
1944  &*this_object->current_refine_pointer) <= 0,
1945  ExcInternalError());
1946 
1947  // now, if the p4est cell is one in the list, it is supposed to be refined
1949  quadrant_is_equal (quadrant, &*this_object->current_refine_pointer))
1950  {
1951  ++this_object->current_refine_pointer;
1952  return true;
1953  }
1954 
1955  // p4est cell is not in list
1956  return false;
1957  }
1958 
1959 
1960 
1961  template <int dim, int spacedim>
1962  int
1963  RefineAndCoarsenList<dim,spacedim>::
1964  coarsen_callback (typename internal::p4est::types<dim>::forest *forest,
1965  typename internal::p4est::types<dim>::topidx coarse_cell_index,
1966  typename internal::p4est::types<dim>::quadrant *children[])
1967  {
1968  RefineAndCoarsenList<dim,spacedim> *this_object
1969  = reinterpret_cast<RefineAndCoarsenList<dim,spacedim>*>(forest->user_pointer);
1970 
1971  // if there are no more cells in our list the current cell can't be
1972  // flagged for coarsening
1973  if (this_object->current_coarsen_pointer ==
1974  this_object->coarsen_list.end())
1975  return false;
1976 
1977  Assert (coarse_cell_index <=
1978  this_object->current_coarsen_pointer->p.which_tree,
1979  ExcInternalError());
1980 
1981  // if p4est hasn't yet reached the tree of the next flagged cell the
1982  // current cell can't be flagged for coarsening
1983  if (coarse_cell_index <
1984  this_object->current_coarsen_pointer->p.which_tree)
1985  return false;
1986 
1987  // now we're in the right tree in the forest
1988  Assert (coarse_cell_index <=
1989  this_object->current_coarsen_pointer->p.which_tree,
1990  ExcInternalError());
1991 
1992  // make sure that the p4est loop over cells hasn't gotten ahead of our own
1993  // pointer
1995  quadrant_compare (children[0],
1996  &*this_object->current_coarsen_pointer) <= 0,
1997  ExcInternalError());
1998 
1999  // now, if the p4est cell is one in the list, it is supposed to be
2000  // coarsened
2002  quadrant_is_equal (children[0],
2003  &*this_object->current_coarsen_pointer))
2004  {
2005  // move current pointer one up
2006  ++this_object->current_coarsen_pointer;
2007 
2008  // note that the next 3 cells in our list need to correspond to the
2009  // other siblings of the cell we have just found
2010  for (unsigned int c=1; c<GeometryInfo<dim>::max_children_per_cell; ++c)
2011  {
2013  quadrant_is_equal (children[c],
2014  &*this_object->current_coarsen_pointer),
2015  ExcInternalError());
2016  ++this_object->current_coarsen_pointer;
2017  }
2018 
2019  return true;
2020  }
2021 
2022  // p4est cell is not in list
2023  return false;
2024  }
2025 
2026 
2027 
2034  template <int dim, int spacedim>
2035  class PartitionWeights
2036  {
2037  public:
2043  PartitionWeights (const std::vector<unsigned int> &cell_weights);
2044 
2052  static
2053  int
2054  cell_weight (typename internal::p4est::types<dim>::forest *forest,
2055  typename internal::p4est::types<dim>::topidx coarse_cell_index,
2056  typename internal::p4est::types<dim>::quadrant *quadrant);
2057 
2058  private:
2059  std::vector<unsigned int> cell_weights_list;
2060  std::vector<unsigned int>::const_iterator current_pointer;
2061  };
2062 
2063 
2064  template <int dim, int spacedim>
2065  PartitionWeights<dim,spacedim>::
2066  PartitionWeights (const std::vector<unsigned int> &cell_weights)
2067  :
2068  cell_weights_list(cell_weights)
2069  {
2070  // set the current pointer to the first element of the list, given that
2071  // we will walk through it sequentially
2072  current_pointer = cell_weights_list.begin();
2073  }
2074 
2075 
2076  template <int dim, int spacedim>
2077  int
2078  PartitionWeights<dim,spacedim>::
2079  cell_weight (typename internal::p4est::types<dim>::forest *forest,
2082  {
2083  // the function gets two additional arguments, but we don't need them
2084  // since we know in which order p4est will walk through the cells
2085  // and have already built our weight lists in this order
2086 
2087  PartitionWeights<dim,spacedim> *this_object
2088  = reinterpret_cast<PartitionWeights<dim,spacedim>*>(forest->user_pointer);
2089 
2090  Assert (this_object->current_pointer >= this_object->cell_weights_list.begin(),
2091  ExcInternalError());
2092  Assert (this_object->current_pointer < this_object->cell_weights_list.end(),
2093  ExcInternalError());
2094 
2095  // get the weight, increment the pointer, and return the weight
2096  return *this_object->current_pointer++;
2097  }
2098 }
2099 
2100 
2101 namespace parallel
2102 {
2103  namespace distributed
2104  {
2105 
2106  /* ---------------------- class Triangulation<dim,spacedim> ------------------------------ */
2107 
2108 
2109  template <int dim, int spacedim>
2111  Triangulation (MPI_Comm mpi_communicator,
2112  const typename ::Triangulation<dim,spacedim>::MeshSmoothing smooth_grid,
2113  const Settings settings_)
2114  :
2115  // do not check for distorted cells
2116  ::parallel::Triangulation<dim,spacedim>
2117  (mpi_communicator,
2118  smooth_grid,
2119  false),
2120  settings(settings_),
2121  triangulation_has_content (false),
2122  connectivity (0),
2123  parallel_forest (0),
2124  refinement_in_progress (false),
2125  attached_data_size(0),
2126  n_attached_datas(0),
2127  n_attached_deserialize(0)
2128  {
2129  parallel_ghost = 0;
2130  }
2131 
2132 
2133 
2134  template <int dim, int spacedim>
2136  {
2137  clear ();
2138 
2140  ExcInternalError());
2141  Assert (connectivity == 0, ExcInternalError());
2142  Assert (parallel_forest == 0, ExcInternalError());
2143  Assert (refinement_in_progress == false, ExcInternalError());
2144  }
2145 
2146 
2147 
2148 
2149  template <int dim, int spacedim>
2150  void
2153  const std::vector<CellData<dim> > &cells,
2154  const SubCellData &subcelldata)
2155  {
2156  try
2157  {
2159  create_triangulation (vertices, cells, subcelldata);
2160  }
2161  catch (const typename ::Triangulation<dim,spacedim>::DistortedCellList &)
2162  {
2163  // the underlying triangulation should not be checking for distorted
2164  // cells
2165  AssertThrow (false, ExcInternalError());
2166  }
2167 
2168  // note that now we have some content in the p4est objects and call the
2169  // functions that do the actual work (which are dimension dependent, so
2170  // separate)
2172 
2174 
2176 
2177  try
2178  {
2180  }
2181  catch (const typename Triangulation<dim>::DistortedCellList &)
2182  {
2183  // the underlying triangulation should not be checking for distorted
2184  // cells
2185  AssertThrow (false, ExcInternalError());
2186  }
2187 
2188  this->update_number_cache ();
2189  }
2190 
2191 
2192  // This anonymous namespace contains utility for
2193  // the function Triangulation::communicate_locally_moved_vertices
2194  namespace CommunicateLocallyMovedVertices
2195  {
2196  namespace
2197  {
2203  template <int dim, int spacedim>
2204  struct CellInfo
2205  {
2206  // store all the tree_indices we send/receive consecutively (n_cells entries)
2207  std::vector<unsigned int> tree_index;
2208  // store all the quadrants we send/receive consecutively (n_cells entries)
2209  std::vector<typename ::internal::p4est::types<dim>::quadrant> quadrants;
2210  // store for each cell the number of vertices we send/receive
2211  // and then the vertex indices (for each cell: n_vertices+1 entries)
2212  std::vector<unsigned int> vertex_indices;
2213  // store for each cell the vertices we send/receive
2214  // (for each cell n_vertices entries)
2215  std::vector<::Point<spacedim> > vertices;
2216  // for receiving and unpacking data we need to store pointers to the
2217  // first vertex and vertex_index on each cell additionally
2218  // both vectors have as many entries as there are cells
2219  std::vector<unsigned int * > first_vertex_indices;
2220  std::vector<::Point<spacedim>* > first_vertices;
2221 
2222  unsigned int bytes_for_buffer () const
2223  {
2224  return (sizeof(unsigned int) +
2225  tree_index.size() * sizeof(unsigned int) +
2226  quadrants.size() * sizeof(typename ::internal::p4est
2227  ::types<dim>::quadrant) +
2228  vertices.size() * sizeof(::Point<spacedim>)) +
2229  vertex_indices.size() * sizeof(unsigned int);
2230  }
2231 
2232  void pack_data (std::vector<char> &buffer) const
2233  {
2234  buffer.resize(bytes_for_buffer());
2235 
2236  char *ptr = &buffer[0];
2237 
2238  const unsigned int num_cells = tree_index.size();
2239  std::memcpy(ptr, &num_cells, sizeof(unsigned int));
2240  ptr += sizeof(unsigned int);
2241 
2242  std::memcpy(ptr,
2243  &tree_index[0],
2244  num_cells*sizeof(unsigned int));
2245  ptr += num_cells*sizeof(unsigned int);
2246 
2247  std::memcpy(ptr,
2248  &quadrants[0],
2249  num_cells * sizeof(typename ::internal::p4est::
2250  types<dim>::quadrant));
2251  ptr += num_cells*sizeof(typename ::internal::p4est::types<dim>::
2252  quadrant);
2253 
2254  std::memcpy(ptr,
2255  &vertex_indices[0],
2256  vertex_indices.size() * sizeof(unsigned int));
2257  ptr += vertex_indices.size() * sizeof(unsigned int);
2258 
2259  std::memcpy(ptr,
2260  &vertices[0],
2261  vertices.size() * sizeof(::Point<spacedim>));
2262  ptr += vertices.size() * sizeof(::Point<spacedim>);
2263 
2264  Assert (ptr == &buffer[0]+buffer.size(),
2265  ExcInternalError());
2266 
2267  }
2268 
2269  void unpack_data (const std::vector<char> &buffer)
2270  {
2271  const char *ptr = &buffer[0];
2272  unsigned int cells;
2273  memcpy(&cells, ptr, sizeof(unsigned int));
2274  ptr += sizeof(unsigned int);
2275 
2276  tree_index.resize(cells);
2277  memcpy(&tree_index[0],ptr,sizeof(unsigned int)*cells);
2278  ptr += sizeof(unsigned int)*cells;
2279 
2280  quadrants.resize(cells);
2281  memcpy(&quadrants[0],ptr,
2282  sizeof(typename ::internal::p4est::types<dim>::quadrant)*cells);
2283  ptr += sizeof(typename ::internal::p4est::types<dim>::quadrant)*cells;
2284 
2285  vertex_indices.clear();
2286  first_vertex_indices.resize(cells);
2287  std::vector<unsigned int> n_vertices_on_cell(cells);
2288  std::vector<unsigned int> first_indices (cells);
2289  for (unsigned int c=0; c<cells; ++c)
2290  {
2291  // The first 'vertex index' is the number of vertices.
2292  // Additionally, we need to store the pointer to this
2293  // vertex index with respect to the std::vector
2294  const unsigned int *const vertex_index
2295  = reinterpret_cast<const unsigned int *const>(ptr);
2296  first_indices[c] = vertex_indices.size();
2297  vertex_indices.push_back(*vertex_index);
2298  n_vertices_on_cell[c] = *vertex_index;
2299  ptr += sizeof(unsigned int);
2300  // Now copy all the 'real' vertex_indices
2301  vertex_indices.resize(vertex_indices.size() + n_vertices_on_cell[c]);
2302  memcpy(&vertex_indices[vertex_indices.size() - n_vertices_on_cell[c]],
2303  ptr, n_vertices_on_cell[c]*sizeof(unsigned int));
2304  ptr += n_vertices_on_cell[c]*sizeof(unsigned int);
2305  }
2306  for (unsigned int c=0; c<cells; ++c)
2307  first_vertex_indices[c] = &vertex_indices[first_indices[c]];
2308 
2309  vertices.clear();
2310  first_vertices.resize(cells);
2311  for (unsigned int c=0; c<cells; ++c)
2312  {
2313  // We need to store a pointer to the first vertex.
2314  const ::Point<spacedim> *const vertex
2315  = reinterpret_cast<const ::Point<spacedim> * const>(ptr);
2316  first_indices[c] = vertices.size();
2317  vertices.push_back(*vertex);
2318  ptr += sizeof(::Point<spacedim>);
2319  vertices.resize(vertices.size() + n_vertices_on_cell[c]-1);
2320  memcpy(&vertices[vertices.size() - (n_vertices_on_cell[c]-1)],
2321  ptr, (n_vertices_on_cell[c]-1)*sizeof(::Point<spacedim>));
2322  ptr += (n_vertices_on_cell[c]-1)*sizeof(::Point<spacedim>);
2323  }
2324  for (unsigned int c=0; c<cells; ++c)
2325  first_vertices[c] = &vertices[first_indices[c]];
2326 
2327  Assert (ptr == &buffer[0]+buffer.size(),
2328  ExcInternalError());
2329  }
2330  };
2331 
2332 
2333 
2334  // This function is responsible for gathering the information
2335  // we want to send to each process.
2336  // For each dealii cell on the coarsest level the corresponding
2337  // p4est_cell has to be provided when calling this function.
2338  // By recursing through all children we consider each active cell.
2339  // vertices_with_ghost_neighbors tells us which vertices
2340  // are in the ghost layer and for which processes they might
2341  // be interesting.
2342  // Whether a vertex has actually been updated locally is
2343  // stored in vertex_locally_moved. Only those are considered
2344  // for sending.
2345  // The gathered information is saved into needs_to_get_cell.
2346  template <int dim, int spacedim>
2347  void
2348  fill_vertices_recursively (const typename parallel::distributed::Triangulation<dim,spacedim> &tria,
2349  const unsigned int tree_index,
2350  const typename Triangulation<dim,spacedim>::cell_iterator &dealii_cell,
2351  const typename ::internal::p4est::types<dim>::quadrant &p4est_cell,
2352  const std::map<unsigned int, std::set<::types::subdomain_id> > &vertices_with_ghost_neighbors,
2353  const std::vector<bool> &vertex_locally_moved,
2354  std::map<::types::subdomain_id, CellInfo<dim, spacedim> > &needs_to_get_cell)
2355  {
2356  // see if we have to
2357  // recurse...
2358  if (dealii_cell->has_children())
2359  {
2360  typename ::internal::p4est::types<dim>::quadrant
2361  p4est_child[GeometryInfo<dim>::max_children_per_cell];
2362  ::internal::p4est::init_quadrant_children<dim>(p4est_cell, p4est_child);
2363 
2364 
2365  for (unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
2366  fill_vertices_recursively<dim,spacedim>(tria,
2367  tree_index,
2368  dealii_cell->child(c),
2369  p4est_child[c],
2370  vertices_with_ghost_neighbors,
2371  vertex_locally_moved,
2372  needs_to_get_cell);
2373  return;
2374  }
2375 
2376  // We're at a leaf cell. If the cell is locally owned, we may
2377  // have to send its vertices to other processors if any of
2378  // its vertices is adjacent to a ghost cell and has been moved
2379  //
2380  // If one of the vertices of the cell is interesting,
2381  // send all moved vertices of the cell to all processors
2382  // adjacent to all cells adjacent to this vertex
2383  if (dealii_cell->is_locally_owned())
2384  {
2385  std::set<::types::subdomain_id> send_to;
2386  for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
2387  {
2388  const std::map<unsigned int, std::set<::types::subdomain_id> >::const_iterator
2389  neighbor_subdomains_of_vertex
2390  = vertices_with_ghost_neighbors.find (dealii_cell->vertex_index(v));
2391 
2392  if (neighbor_subdomains_of_vertex
2393  != vertices_with_ghost_neighbors.end())
2394  {
2395  Assert(neighbor_subdomains_of_vertex->second.size()!=0,
2396  ExcInternalError());
2397  send_to.insert(neighbor_subdomains_of_vertex->second.begin(),
2398  neighbor_subdomains_of_vertex->second.end());
2399  }
2400  }
2401 
2402  if (send_to.size() > 0)
2403  {
2404  std::vector<unsigned int> vertex_indices;
2405  std::vector<::Point<spacedim> > local_vertices;
2406  for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
2407  if (vertex_locally_moved[dealii_cell->vertex_index(v)])
2408  {
2409  vertex_indices.push_back(v);
2410  local_vertices.push_back(dealii_cell->vertex(v));
2411  }
2412 
2413  if (vertex_indices.size()>0)
2414  for (std::set<::types::subdomain_id>::iterator it=send_to.begin();
2415  it!=send_to.end(); ++it)
2416  {
2417  const ::types::subdomain_id subdomain = *it;
2418 
2419  // get an iterator to what needs to be sent to that
2420  // subdomain (if already exists), or create such an object
2421  const typename std::map<::types::subdomain_id, CellInfo<dim, spacedim> >::iterator
2422  p
2423  = needs_to_get_cell.insert (std::make_pair(subdomain,
2424  CellInfo<dim,spacedim>()))
2425  .first;
2426 
2427  p->second.tree_index.push_back(tree_index);
2428  p->second.quadrants.push_back(p4est_cell);
2429 
2430  p->second.vertex_indices.push_back(vertex_indices.size());
2431  p->second.vertex_indices.insert(p->second.vertex_indices.end(),
2432  vertex_indices.begin(),
2433  vertex_indices.end());
2434 
2435  p->second.vertices.insert(p->second.vertices.end(),
2436  local_vertices.begin(),
2437  local_vertices.end());
2438  }
2439  }
2440  }
2441  }
2442 
2443 
2444 
2445  // After the cell data has been received this function is responsible
2446  // for moving the vertices in the corresponding ghost layer locally.
2447  // As in fill_vertices_recursively for each dealii cell on the
2448  // coarsest level the corresponding p4est_cell has to be provided
2449  // when calling this function. By recursing through through all
2450  // children we consider each active cell.
2451  // Additionally, we need to give a pointer to the first vertex indices
2452  // and vertices. Since the first information saved in vertex_indices
2453  // is the number of vertices this all the information we need.
2454  template <int dim, int spacedim>
2455  void
2456  set_vertices_recursively (
2458  const typename ::internal::p4est::types<dim>::quadrant &p4est_cell,
2459  const typename Triangulation<dim,spacedim>::cell_iterator &dealii_cell,
2460  const typename ::internal::p4est::types<dim>::quadrant &quadrant,
2461  const ::Point<spacedim> *const vertices,
2462  const unsigned int *const vertex_indices)
2463  {
2464  if (::internal::p4est::quadrant_is_equal<dim>(p4est_cell, quadrant))
2465  {
2466  Assert(!dealii_cell->is_artificial(), ExcInternalError());
2467  Assert(!dealii_cell->has_children(), ExcInternalError());
2468  Assert(!dealii_cell->is_locally_owned(), ExcInternalError());
2469 
2470  const unsigned int n_vertices = vertex_indices[0];
2471 
2472  // update dof indices of cell
2473  for (unsigned int i=0; i<n_vertices; ++i)
2474  dealii_cell->vertex(vertex_indices[i+1]) = vertices[i];
2475 
2476  return;
2477  }
2478 
2479  if (! dealii_cell->has_children())
2480  return;
2481 
2482  if (! ::internal::p4est::quadrant_is_ancestor<dim> (p4est_cell, quadrant))
2483  return;
2484 
2485  typename ::internal::p4est::types<dim>::quadrant
2486  p4est_child[GeometryInfo<dim>::max_children_per_cell];
2487  ::internal::p4est::init_quadrant_children<dim>(p4est_cell, p4est_child);
2488 
2489  for (unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
2490  set_vertices_recursively<dim,spacedim> (tria, p4est_child[c],
2491  dealii_cell->child(c),
2492  quadrant, vertices,
2493  vertex_indices);
2494  }
2495  }
2496  }
2497 
2498 
2499 
2500  template <int dim, int spacedim>
2501  void
2503  {
2504  triangulation_has_content = false;
2505 
2506  if (parallel_ghost != 0)
2507  {
2509  parallel_ghost = 0;
2510  }
2511 
2512  if (parallel_forest != 0)
2513  {
2515  parallel_forest = 0;
2516  }
2517 
2518  if (connectivity != 0)
2519  {
2521  connectivity = 0;
2522  }
2523 
2524  coarse_cell_to_p4est_tree_permutation.resize (0);
2525  p4est_tree_to_coarse_cell_permutation.resize (0);
2526 
2528 
2530 
2531  this->update_number_cache ();
2532  }
2533 
2534  template <int dim, int spacedim>
2535  bool
2537  {
2538  if (this->n_global_levels()<=1)
2539  return false; // can not have hanging nodes without refined cells
2540 
2541  // if there are any active cells with level less than n_global_levels()-1, then
2542  // there is obviously also one with level n_global_levels()-1, and
2543  // consequently there must be a hanging node somewhere.
2544  //
2545  // The problem is that we cannot just ask for the first active cell, but
2546  // instead need to filter over locally owned cells.
2547  bool have_coarser_cell = false;
2549  cell != this->end(this->n_global_levels()-2);
2550  ++cell)
2551  if (cell->is_locally_owned())
2552  {
2553  have_coarser_cell = true;
2554  break;
2555  }
2556 
2557  // return true if at least one process has a coarser cell
2558  return 0<Utilities::MPI::max(have_coarser_cell?1:0, this->mpi_communicator);
2559  }
2560 
2561 
2562 
2563  template <int dim, int spacedim>
2564  void
2566  {
2567  DynamicSparsityPattern cell_connectivity;
2568  GridTools::get_vertex_connectivity_of_cells (*this, cell_connectivity);
2569  coarse_cell_to_p4est_tree_permutation.resize (this->n_cells(0));
2571  reorder_hierarchical (cell_connectivity,
2572  coarse_cell_to_p4est_tree_permutation);
2573 
2574  p4est_tree_to_coarse_cell_permutation
2575  = Utilities::invert_permutation (coarse_cell_to_p4est_tree_permutation);
2576  }
2577 
2578 
2579 
2580  template <int dim, int spacedim>
2581  void
2582  Triangulation<dim,spacedim>::write_mesh_vtk (const char *file_basename) const
2583  {
2584  Assert (parallel_forest != 0,
2585  ExcMessage ("Can't produce output when no forest is created yet."));
2587  vtk_write_file (parallel_forest, 0, file_basename);
2588  }
2589 
2590 
2591  template <int dim, int spacedim>
2592  void
2594  save(const char *filename) const
2595  {
2597  ExcMessage ("not all SolutionTransfer's got deserialized after the last load()"));
2598  int real_data_size = 0;
2599  if (attached_data_size>0)
2600  real_data_size = attached_data_size+sizeof(CellStatus);
2601 
2602  Assert(this->n_cells()>0, ExcMessage("Can not save() an empty Triangulation."));
2603 
2604  if (this->my_subdomain==0)
2605  {
2606  std::string fname=std::string(filename)+".info";
2607  std::ofstream f(fname.c_str());
2608  f << "version nproc attached_bytes n_attached_objs n_coarse_cells" << std::endl
2609  << 2 << " "
2611  << real_data_size << " "
2612  << attached_data_pack_callbacks.size() << " "
2613  << this->n_cells(0)
2614  << std::endl;
2615  }
2616 
2617  if (attached_data_size>0)
2618  {
2620  ->attach_mesh_data();
2621  }
2622 
2623  ::internal::p4est::functions<dim>::save(filename, parallel_forest, attached_data_size>0);
2624 
2627 
2628  tria->n_attached_datas = 0;
2629  tria->attached_data_size = 0;
2630  tria->attached_data_pack_callbacks.clear();
2631 
2632  // and release the data
2633  void *userptr = parallel_forest->user_pointer;
2634  ::internal::p4est::functions<dim>::reset_data (parallel_forest, 0, NULL, NULL);
2635  parallel_forest->user_pointer = userptr;
2636  }
2637 
2638 
2639  template <int dim, int spacedim>
2640  void
2642  load (const char *filename,
2643  const bool autopartition)
2644  {
2645  Assert(this->n_cells()>0, ExcMessage("load() only works if the Triangulation already contains a coarse mesh!"));
2646  Assert(this->n_levels()==1, ExcMessage("Triangulation may only contain coarse cells when calling load()."));
2647 
2648 
2649  if (parallel_ghost != 0)
2650  {
2652  parallel_ghost = 0;
2653  }
2655  parallel_forest = 0;
2657  connectivity = 0;
2658 
2659  unsigned int version, numcpus, attached_size, attached_count, n_coarse_cells;
2660  {
2661  std::string fname=std::string(filename)+".info";
2662  std::ifstream f(fname.c_str());
2663  std::string firstline;
2664  getline(f, firstline); //skip first line
2665  f >> version >> numcpus >> attached_size >> attached_count >> n_coarse_cells;
2666  }
2667 
2668  Assert(version == 2, ExcMessage("Incompatible version found in .info file."));
2669  Assert(this->n_cells(0) == n_coarse_cells, ExcMessage("Number of coarse cells differ!"));
2670 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
2671 #else
2673  ExcMessage("parallel::distributed::Triangulation::load() only supports loading "
2674  "saved data with a greater or equal number of processes than were used to "
2675  "save() when using p4est 0.3.4.2."));
2676 #endif
2677 
2678  attached_data_size = 0;
2679  n_attached_datas = 0;
2680  n_attached_deserialize = attached_count;
2681 
2682 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3)
2684  filename, this->mpi_communicator,
2685  attached_size, attached_size>0,
2686  autopartition, 0,
2687  this,
2688  &connectivity);
2689 #else
2690  (void)autopartition;
2691  parallel_forest = ::internal::p4est::functions<dim>::load (
2692  filename, this->mpi_communicator,
2693  attached_size, attached_size>0,
2694  this,
2695  &connectivity);
2696 #endif
2697  if (numcpus != Utilities::MPI::n_mpi_processes (this->mpi_communicator))
2698  // We are changing the number of CPUs so we need to repartition.
2699  // Note that p4est actually distributes the cells between the changed
2700  // number of CPUs and so everything works without this call, but
2701  // this command changes the distribution for some reason, so we
2702  // will leave it in here.
2703  repartition();
2704 
2705  try
2706  {
2708  }
2709  catch (const typename Triangulation<dim>::DistortedCellList &)
2710  {
2711  // the underlying
2712  // triangulation should not
2713  // be checking for
2714  // distorted cells
2715  AssertThrow (false, ExcInternalError());
2716  }
2717 
2718  this->update_number_cache ();
2719  }
2720 
2721 
2722 
2723  template <int dim, int spacedim>
2724  unsigned int
2726  {
2727  Assert (parallel_forest != 0,
2728  ExcMessage ("Can't produce a check sum when no forest is created yet."));
2729  return ::internal::p4est::functions<dim>::checksum (parallel_forest);
2730  }
2731 
2732  template <int dim, int spacedim>
2733  void
2736  {
2738 
2739  if (this->n_levels() == 0)
2740  return;
2741 
2743  {
2744  // find level ghost owners
2746  cell = this->begin();
2747  cell != this->end();
2748  ++cell)
2749  if (cell->level_subdomain_id() != numbers::artificial_subdomain_id
2750  && cell->level_subdomain_id() != this->locally_owned_subdomain())
2751  this->number_cache.level_ghost_owners.insert(cell->level_subdomain_id());
2752 
2753  Assert(this->number_cache.level_ghost_owners.size() < Utilities::MPI::n_mpi_processes(this->mpi_communicator), ExcInternalError());
2754  }
2755  }
2756 
2757 
2758  template <int dim, int spacedim>
2759  typename ::internal::p4est::types<dim>::tree *
2761  init_tree(const int dealii_coarse_cell_index) const
2762  {
2763  const unsigned int tree_index
2764  = coarse_cell_to_p4est_tree_permutation[dealii_coarse_cell_index];
2765  typename ::internal::p4est::types<dim>::tree *tree
2766  = static_cast<typename ::internal::p4est::types<dim>::tree *>
2767  (sc_array_index (parallel_forest->trees,
2768  tree_index));
2769 
2770  return tree;
2771  }
2772 
2773 
2774 
2775  template <>
2776  void
2778  {
2779  const unsigned int dim = 2, spacedim = 2;
2780  Assert (this->n_cells(0) > 0, ExcInternalError());
2781  Assert (this->n_levels() == 1, ExcInternalError());
2782 
2783  // data structures that counts how many cells touch each vertex
2784  // (vertex_touch_count), and which cells touch a given vertex (together
2785  // with the local numbering of that vertex within the cells that touch
2786  // it)
2787  std::vector<unsigned int> vertex_touch_count;
2788  std::vector<
2789  std::list<
2790  std::pair<Triangulation<dim,spacedim>::active_cell_iterator,
2791  unsigned int> > >
2792  vertex_to_cell;
2793  get_vertex_to_cell_mappings (*this,
2794  vertex_touch_count,
2795  vertex_to_cell);
2796  const ::internal::p4est::types<2>::locidx
2797  num_vtt = std::accumulate (vertex_touch_count.begin(),
2798  vertex_touch_count.end(),
2799  0);
2800 
2801  // now create a connectivity object with the right sizes for all
2802  // arrays. set vertex information only in debug mode (saves a few bytes
2803  // in optimized mode)
2804  const bool set_vertex_info
2805 #ifdef DEBUG
2806  = true
2807 #else
2808  = false
2809 #endif
2810  ;
2811 
2812  connectivity
2814  connectivity_new ((set_vertex_info == true ? this->n_vertices() : 0),
2815  this->n_cells(0),
2816  this->n_vertices(),
2817  num_vtt);
2818 
2819  set_vertex_and_cell_info (*this,
2820  vertex_touch_count,
2821  vertex_to_cell,
2822  coarse_cell_to_p4est_tree_permutation,
2823  set_vertex_info,
2824  connectivity);
2825 
2826  Assert (p4est_connectivity_is_valid (connectivity) == 1,
2827  ExcInternalError());
2828 
2829  // now create a forest out of the connectivity data structure
2830  parallel_forest
2833  connectivity,
2834  /* minimum initial number of quadrants per tree */ 0,
2835  /* minimum level of upfront refinement */ 0,
2836  /* use uniform upfront refinement */ 1,
2837  /* user_data_size = */ 0,
2838  /* user_data_constructor = */ NULL,
2839  /* user_pointer */ this);
2840  }
2841 
2842 
2843 
2844  // TODO: This is a verbatim copy of the 2,2 case. However, we can't just
2845  // specialize the dim template argument, but let spacedim open
2846  template <>
2847  void
2849  {
2850  const unsigned int dim = 2, spacedim = 3;
2851  Assert (this->n_cells(0) > 0, ExcInternalError());
2852  Assert (this->n_levels() == 1, ExcInternalError());
2853 
2854  // data structures that counts how many cells touch each vertex
2855  // (vertex_touch_count), and which cells touch a given vertex (together
2856  // with the local numbering of that vertex within the cells that touch
2857  // it)
2858  std::vector<unsigned int> vertex_touch_count;
2859  std::vector<
2860  std::list<
2861  std::pair<Triangulation<dim,spacedim>::active_cell_iterator,
2862  unsigned int> > >
2863  vertex_to_cell;
2864  get_vertex_to_cell_mappings (*this,
2865  vertex_touch_count,
2866  vertex_to_cell);
2867  const ::internal::p4est::types<2>::locidx
2868  num_vtt = std::accumulate (vertex_touch_count.begin(),
2869  vertex_touch_count.end(),
2870  0);
2871 
2872  // now create a connectivity object with the right sizes for all
2873  // arrays. set vertex information only in debug mode (saves a few bytes
2874  // in optimized mode)
2875  const bool set_vertex_info
2876 #ifdef DEBUG
2877  = true
2878 #else
2879  = false
2880 #endif
2881  ;
2882 
2883  connectivity
2885  connectivity_new ((set_vertex_info == true ? this->n_vertices() : 0),
2886  this->n_cells(0),
2887  this->n_vertices(),
2888  num_vtt);
2889 
2890  set_vertex_and_cell_info (*this,
2891  vertex_touch_count,
2892  vertex_to_cell,
2893  coarse_cell_to_p4est_tree_permutation,
2894  set_vertex_info,
2895  connectivity);
2896 
2897  Assert (p4est_connectivity_is_valid (connectivity) == 1,
2898  ExcInternalError());
2899 
2900  // now create a forest out of the connectivity data structure
2901  parallel_forest
2904  connectivity,
2905  /* minimum initial number of quadrants per tree */ 0,
2906  /* minimum level of upfront refinement */ 0,
2907  /* use uniform upfront refinement */ 1,
2908  /* user_data_size = */ 0,
2909  /* user_data_constructor = */ NULL,
2910  /* user_pointer */ this);
2911  }
2912 
2913 
2914 
2915  template <>
2916  void
2918  {
2919  const int dim = 3, spacedim = 3;
2920  Assert (this->n_cells(0) > 0, ExcInternalError());
2921  Assert (this->n_levels() == 1, ExcInternalError());
2922 
2923  // data structures that counts how many cells touch each vertex
2924  // (vertex_touch_count), and which cells touch a given vertex (together
2925  // with the local numbering of that vertex within the cells that touch
2926  // it)
2927  std::vector<unsigned int> vertex_touch_count;
2928  std::vector<
2929  std::list<
2930  std::pair<Triangulation<3>::active_cell_iterator,
2931  unsigned int> > >
2932  vertex_to_cell;
2933  get_vertex_to_cell_mappings (*this,
2934  vertex_touch_count,
2935  vertex_to_cell);
2936  const ::internal::p4est::types<2>::locidx
2937  num_vtt = std::accumulate (vertex_touch_count.begin(),
2938  vertex_touch_count.end(),
2939  0);
2940 
2941  std::vector<unsigned int> edge_touch_count;
2942  std::vector<
2943  std::list<
2944  std::pair<Triangulation<3>::active_cell_iterator,
2945  unsigned int> > >
2946  edge_to_cell;
2947  get_edge_to_cell_mappings (*this,
2948  edge_touch_count,
2949  edge_to_cell);
2950  const ::internal::p4est::types<2>::locidx
2951  num_ett = std::accumulate (edge_touch_count.begin(),
2952  edge_touch_count.end(),
2953  0);
2954 
2955  // now create a connectivity object with the right sizes for all arrays
2956  const bool set_vertex_info
2957 #ifdef DEBUG
2958  = true
2959 #else
2960  = false
2961 #endif
2962  ;
2963 
2964  connectivity
2966  connectivity_new ((set_vertex_info == true ? this->n_vertices() : 0),
2967  this->n_cells(0),
2968  this->n_active_lines(),
2969  num_ett,
2970  this->n_vertices(),
2971  num_vtt);
2972 
2973  set_vertex_and_cell_info (*this,
2974  vertex_touch_count,
2975  vertex_to_cell,
2976  coarse_cell_to_p4est_tree_permutation,
2977  set_vertex_info,
2978  connectivity);
2979 
2980  // next to tree-to-edge
2981  // data. note that in p4est lines
2982  // are ordered as follows
2983  // *---3---* *---3---*
2984  // /| | / /|
2985  // 6 | 11 6 7 11
2986  // / 10 | / / |
2987  // * | | *---2---* |
2988  // | *---1---* | | *
2989  // | / / | 9 /
2990  // 8 4 5 8 | 5
2991  // |/ / | |/
2992  // *---0---* *---0---*
2993  // whereas in deal.II they are like this:
2994  // *---7---* *---7---*
2995  // /| | / /|
2996  // 4 | 11 4 5 11
2997  // / 10 | / / |
2998  // * | | *---6---* |
2999  // | *---3---* | | *
3000  // | / / | 9 /
3001  // 8 0 1 8 | 1
3002  // |/ / | |/
3003  // *---2---* *---2---*
3004 
3005  const unsigned int deal_to_p4est_line_index[12]
3006  = { 4, 5, 0, 1, 6, 7, 2, 3, 8, 9, 10, 11 } ;
3007 
3008  for (Triangulation<dim,spacedim>::active_cell_iterator
3009  cell = this->begin_active();
3010  cell != this->end(); ++cell)
3011  {
3012  const unsigned int
3013  index = coarse_cell_to_p4est_tree_permutation[cell->index()];
3014  for (unsigned int e=0; e<GeometryInfo<3>::lines_per_cell; ++e)
3015  connectivity->tree_to_edge[index*GeometryInfo<3>::lines_per_cell+
3016  deal_to_p4est_line_index[e]]
3017  = cell->line(e)->index();
3018  }
3019 
3020  // now also set edge-to-tree
3021  // information
3022  connectivity->ett_offset[0] = 0;
3023  std::partial_sum (edge_touch_count.begin(),
3024  edge_touch_count.end(),
3025  &connectivity->ett_offset[1]);
3026 
3027  Assert (connectivity->ett_offset[this->n_active_lines()] ==
3028  num_ett,
3029  ExcInternalError());
3030 
3031  for (unsigned int v=0; v<this->n_active_lines(); ++v)
3032  {
3033  Assert (edge_to_cell[v].size() == edge_touch_count[v],
3034  ExcInternalError());
3035 
3036  std::list<std::pair
3037  <Triangulation<dim,spacedim>::active_cell_iterator,
3038  unsigned int> >::const_iterator
3039  p = edge_to_cell[v].begin();
3040  for (unsigned int c=0; c<edge_touch_count[v]; ++c, ++p)
3041  {
3042  connectivity->edge_to_tree[connectivity->ett_offset[v]+c]
3043  = coarse_cell_to_p4est_tree_permutation[p->first->index()];
3044  connectivity->edge_to_edge[connectivity->ett_offset[v]+c]
3045  = deal_to_p4est_line_index[p->second];
3046  }
3047  }
3048 
3049  Assert (p8est_connectivity_is_valid (connectivity) == 1,
3050  ExcInternalError());
3051 
3052  // now create a forest out of the connectivity data structure
3053  parallel_forest
3056  connectivity,
3057  /* minimum initial number of quadrants per tree */ 0,
3058  /* minimum level of upfront refinement */ 0,
3059  /* use uniform upfront refinement */ 1,
3060  /* user_data_size = */ 0,
3061  /* user_data_constructor = */ NULL,
3062  /* user_pointer */ this);
3063  }
3064 
3065 
3066 
3067  namespace
3068  {
3069  // this function combines vertices that have different locations (and
3070  // thus, different vertex_index) but represent the same topological
3071  // entity over periodic boundaries. The vector
3072  // topological_vertex_numbering contains a linear map from 0 to
3073  // n_vertices at input and at output relates periodic vertices with only
3074  // one vertex index. The output is used to always identify the same
3075  // vertex according to the periodicity, e.g. when finding the maximum
3076  // cell level around a vertex.
3077  //
3078  // Example: On a 3D cell with vertices numbered from 0 to 7 and periodic
3079  // boundary conditions in x direction, the vector
3080  // topological_vertex_numbering will contain the numbers
3081  // {0,0,2,2,4,4,6,6} (because the vertex pairs {0,1}, {2,3}, {4,5},
3082  // {6,7} belong together, respectively). If periodicity is set in x and
3083  // z direction, the output is {0,0,2,2,0,0,2,2}, and if periodicity is
3084  // in all directions, the output is simply {0,0,0,0,0,0,0,0}.
3085  template <typename ITERATOR>
3086  void
3087  identify_periodic_vertices_recursively(const GridTools::PeriodicFacePair<ITERATOR> &periodic,
3088  std::vector<unsigned int> &topological_vertex_numbering)
3089  {
3090  const unsigned int dim = ITERATOR::AccessorType::dimension;
3091 
3092  // for hanging nodes we will consider all necessary coupling already
3093  // on the parent level, so we just need to consider neighbors of the
3094  // same level
3095  if (periodic.cell[0]->has_children() &&
3096  periodic.cell[1]->has_children())
3097  {
3098  // copy orientations etc. from parent to child
3099  GridTools::PeriodicFacePair<ITERATOR> periodic_child = periodic;
3100 
3101  // find appropriate pairs of child elements
3102  for (unsigned int cf=0; cf<periodic.cell[0]->face(periodic.face_idx[0])->n_children(); ++cf)
3103  {
3104  const unsigned int child_index_0 =
3105  GeometryInfo<dim>::child_cell_on_face(periodic.cell[0]->refinement_case(),
3106  periodic.face_idx[0], cf,
3107  periodic.orientation[0],
3108  periodic.orientation[1],
3109  periodic.orientation[2]);
3110  periodic_child.cell[0] = periodic.cell[0]->child(child_index_0);
3111  periodic_child.face_idx[0] = periodic.face_idx[0];
3112 
3113  // the second face is in standard orientation in terms of the
3114  // periodic face pair
3115  const unsigned int child_index_1 =
3116  GeometryInfo<dim>::child_cell_on_face(periodic.cell[1]->refinement_case(),
3117  periodic.face_idx[1], cf);
3118 
3119  periodic_child.cell[1] = periodic.cell[1]->child(child_index_1);
3120  periodic_child.face_idx[1] = periodic.face_idx[1];
3121 
3122  // recursive call into children
3123  identify_periodic_vertices_recursively (periodic_child,
3124  topological_vertex_numbering);
3125  }
3126  }
3127 
3128  for (unsigned int v=0; v<GeometryInfo<dim-1>::vertices_per_cell; ++v)
3129  {
3130  // take possible non-standard orientation of face on cell[0] into
3131  // account
3132  const unsigned int vface0 =
3134  periodic.orientation[1],
3135  periodic.orientation[2]);
3136  const unsigned int vi0 = topological_vertex_numbering[periodic.cell[0]->face(periodic.face_idx[0])->vertex_index(vface0)];
3137  const unsigned int vi1 = topological_vertex_numbering[periodic.cell[1]->face(periodic.face_idx[1])->vertex_index(v)];
3138  const unsigned int min_index = std::min(vi0, vi1);
3139  topological_vertex_numbering[periodic.cell[0]->face(periodic.face_idx[0])->vertex_index(vface0)]
3140  = topological_vertex_numbering[periodic.cell[1]->face(periodic.face_idx[1])->vertex_index(v)]
3141  = min_index;
3142  }
3143  }
3144 
3145 
3146 
3147  // ensures the 2:1 mesh balance for periodic boundary conditions in the
3148  // artificial cell layer (the active cells are taken care of by p4est)
3149  template <int dim, int spacedim>
3150  bool enforce_mesh_balance_over_periodic_boundaries
3152  const std::vector<GridTools::PeriodicFacePair<typename ::Triangulation<dim,spacedim>::cell_iterator> > periodic_face_pairs_level_0)
3153  {
3154  if (periodic_face_pairs_level_0.empty())
3155  return false;
3156 
3157  std::vector<bool> flags_before[2];
3158  tria.save_coarsen_flags (flags_before[0]);
3159  tria.save_refine_flags (flags_before[1]);
3160 
3161  std::vector<unsigned int> topological_vertex_numbering(tria.n_vertices());
3162  for (unsigned int i=0; i<topological_vertex_numbering.size(); ++i)
3163  topological_vertex_numbering[i] = i;
3164  for (unsigned int i=0; i<periodic_face_pairs_level_0.size(); ++i)
3165  {
3166  identify_periodic_vertices_recursively(periodic_face_pairs_level_0[i],
3167  topological_vertex_numbering);
3168  }
3169 
3170  // this code is replicated from grid/tria.cc but using an indirection
3171  // for periodic boundary conditions
3172  bool continue_iterating = true;
3173  std::vector<int> vertex_level(tria.n_vertices());
3174  while (continue_iterating)
3175  {
3176  // store highest level one of the cells adjacent to a vertex
3177  // belongs to
3178  std::fill (vertex_level.begin(), vertex_level.end(), 0);
3179  typename Triangulation<dim,spacedim>::active_cell_iterator
3180  cell = tria.begin_active(), endc = tria.end();
3181  for (; cell!=endc; ++cell)
3182  {
3183  if (cell->refine_flag_set())
3184  for (unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_cell;
3185  ++vertex)
3186  vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]]
3187  = std::max (vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]],
3188  cell->level()+1);
3189  else if (!cell->coarsen_flag_set())
3190  for (unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_cell;
3191  ++vertex)
3192  vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]]
3193  = std::max (vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]],
3194  cell->level());
3195  else
3196  {
3197  // if coarsen flag is set then tentatively assume
3198  // that the cell will be coarsened. this isn't
3199  // always true (the coarsen flag could be removed
3200  // again) and so we may make an error here. we try
3201  // to correct this by iterating over the entire
3202  // process until we are converged
3203  Assert (cell->coarsen_flag_set(), ExcInternalError());
3204  for (unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_cell;
3205  ++vertex)
3206  vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]]
3207  = std::max (vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]],
3208  cell->level()-1);
3209  }
3210  }
3211 
3212  continue_iterating = false;
3213 
3214  // loop over all cells in reverse order. do so because we
3215  // can then update the vertex levels on the adjacent
3216  // vertices and maybe already flag additional cells in this
3217  // loop
3218  //
3219  // note that not only may we have to add additional
3220  // refinement flags, but we will also have to remove
3221  // coarsening flags on cells adjacent to vertices that will
3222  // see refinement
3223  for (cell=tria.last_active(); cell != endc; --cell)
3224  if (cell->refine_flag_set() == false)
3225  {
3226  for (unsigned int vertex=0;
3227  vertex<GeometryInfo<dim>::vertices_per_cell; ++vertex)
3228  if (vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]] >=
3229  cell->level()+1)
3230  {
3231  // remove coarsen flag...
3232  cell->clear_coarsen_flag();
3233 
3234  // ...and if necessary also refine the current
3235  // cell, at the same time updating the level
3236  // information about vertices
3237  if (vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]] >
3238  cell->level()+1)
3239  {
3240  cell->set_refine_flag();
3241  continue_iterating = true;
3242 
3243  for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell;
3244  ++v)
3245  vertex_level[topological_vertex_numbering[cell->vertex_index(v)]]
3246  = std::max (vertex_level[topological_vertex_numbering[cell->vertex_index(v)]],
3247  cell->level()+1);
3248  }
3249 
3250  // continue and see whether we may, for example,
3251  // go into the inner 'if' above based on a
3252  // different vertex
3253  }
3254  }
3255 
3256  // clear coarsen flag if not all children were marked
3257  for (typename Triangulation<dim,spacedim>::cell_iterator cell = tria.begin();
3258  cell!=tria.end(); ++cell)
3259  {
3260  // nothing to do if we are already on the finest level
3261  if (cell->active())
3262  continue;
3263 
3264  const unsigned int n_children=cell->n_children();
3265  unsigned int flagged_children=0;
3266  for (unsigned int child=0; child<n_children; ++child)
3267  if (cell->child(child)->active() &&
3268  cell->child(child)->coarsen_flag_set())
3269  ++flagged_children;
3270 
3271  // if not all children were flagged for coarsening, remove
3272  // coarsen flags
3273  if (flagged_children < n_children)
3274  for (unsigned int child=0; child<n_children; ++child)
3275  if (cell->child(child)->active())
3276  cell->child(child)->clear_coarsen_flag();
3277  }
3278  }
3279  std::vector<bool> flags_after[2];
3280  tria.save_coarsen_flags (flags_after[0]);
3281  tria.save_refine_flags (flags_after[1]);
3282  return ((flags_before[0] != flags_after[0]) ||
3283  (flags_before[1] != flags_after[1]));
3284  }
3285  }
3286 
3287 
3288 
3289  template <int dim, int spacedim>
3290  bool
3292  {
3293  std::vector<bool> flags_before[2];
3294  this->save_coarsen_flags (flags_before[0]);
3295  this->save_refine_flags (flags_before[1]);
3296 
3297  bool mesh_changed = false;
3298  do
3299  {
3301 
3302  // enforce 2:1 mesh balance over periodic boundaries
3303  if (this->smooth_grid &
3305  mesh_changed = enforce_mesh_balance_over_periodic_boundaries(*this,
3307  }
3308  while (mesh_changed);
3309 
3310  // check if any of the refinement flags were changed during this
3311  // function and return that value
3312  std::vector<bool> flags_after[2];
3313  this->save_coarsen_flags (flags_after[0]);
3314  this->save_refine_flags (flags_after[1]);
3315  return ((flags_before[0] != flags_after[0]) ||
3316  (flags_before[1] != flags_after[1]));
3317  }
3318 
3319 
3320 
3321  template <int dim, int spacedim>
3322  void
3324  {
3325  // disable mesh smoothing for recreating the deal.II triangulation,
3326  // otherwise we might not be able to reproduce the p4est mesh
3327  // exactly. We restore the original smoothing at the end of this
3328  // function. Note that the smoothing flag is used in the normal
3329  // refinement process.
3331  save_smooth = this->smooth_grid;
3332 
3333  // We will refine manually to match the p4est further down, which
3334  // obeys a level difference of 2 at each vertex (see the balance call
3335  // to p4est). We can disable this here so we store fewer artificial
3336  // cells (in some cases). For geometric multigrid it turns out that
3337  // we will miss level cells at shared vertices if we ignore this.
3338  // See tests/mpi/mg_06.
3341  else
3343 
3344  bool mesh_changed = false;
3345 
3346  // remove all deal.II refinements. Note that we could skip this and
3347  // start from our current state, because the algorithm later coarsens as
3348  // necessary. This has the advantage of being faster when large parts
3349  // of the local partition changes (likely) and gives a deterministic
3350  // ordering of the cells (useful for snapshot/resume).
3351  // TODO: is there a more efficient way to do this?
3353  while (this->begin_active()->level() > 0)
3354  {
3356  cell = this->begin_active();
3357  cell != this->end();
3358  ++cell)
3359  {
3360  cell->set_coarsen_flag();
3361  }
3362 
3364  const bool saved_refinement_in_progress = refinement_in_progress;
3365  refinement_in_progress = true;
3366 
3367  try
3368  {
3370  }
3371  catch (const typename Triangulation<dim, spacedim>::DistortedCellList &)
3372  {
3373  // the underlying triangulation should not be checking for
3374  // distorted cells
3375  AssertThrow (false, ExcInternalError());
3376  }
3377 
3378  refinement_in_progress = saved_refinement_in_progress;
3379  }
3380 
3381 
3382  // query p4est for the ghost cells
3383  if (parallel_ghost != 0)
3384  {
3386  parallel_ghost = 0;
3387  }
3389  (dim == 2
3390  ?
3391  typename ::internal::p4est::types<dim>::
3392  balance_type(P4EST_CONNECT_CORNER)
3393  :
3394  typename ::internal::p4est::types<dim>::
3395  balance_type(P8EST_CONNECT_CORNER)));
3396 
3397  Assert (parallel_ghost, ExcInternalError());
3398 
3399 
3400  // set all cells to artificial. we will later set it to the correct
3401  // subdomain in match_tree_recursively
3403  cell = this->begin(0);
3404  cell != this->end(0);
3405  ++cell)
3406  cell->recursively_set_subdomain_id(numbers::artificial_subdomain_id);
3407 
3408  do
3409  {
3411  cell = this->begin(0);
3412  cell != this->end(0);
3413  ++cell)
3414  {
3415  // if this processor stores no part of the forest that comes out
3416  // of this coarse grid cell, then we need to delete all children
3417  // of this cell (the coarse grid cell remains)
3418  if (tree_exists_locally<dim,spacedim>(parallel_forest,
3419  coarse_cell_to_p4est_tree_permutation[cell->index()])
3420  == false)
3421  {
3422  delete_all_children<dim,spacedim> (cell);
3423  if (!cell->has_children())
3424  cell->set_subdomain_id (numbers::artificial_subdomain_id);
3425  }
3426 
3427  else
3428  {
3429  // this processor stores at least a part of the tree that
3430  // comes out of this cell.
3431 
3432  typename ::internal::p4est::types<dim>::quadrant p4est_coarse_cell;
3433  typename ::internal::p4est::types<dim>::tree *tree =
3434  init_tree(cell->index());
3435 
3436  ::internal::p4est::init_coarse_quadrant<dim>(p4est_coarse_cell);
3437 
3438  match_tree_recursively<dim,spacedim> (*tree, cell,
3439  p4est_coarse_cell,
3440  *parallel_forest,
3441  this->my_subdomain);
3442  }
3443  }
3444 
3445  // check mesh for ghostcells, refine as necessary. iterate over
3446  // every ghostquadrant, find corresponding deal coarsecell and
3447  // recurse.
3448  typename ::internal::p4est::types<dim>::quadrant *quadr;
3449  unsigned int ghost_owner=0;
3450  typename ::internal::p4est::types<dim>::topidx ghost_tree=0;
3451 
3452  for (unsigned int g_idx=0; g_idx<parallel_ghost->ghosts.elem_count; ++g_idx)
3453  {
3454  while (g_idx >= (unsigned int)parallel_ghost->proc_offsets[ghost_owner+1])
3455  ++ghost_owner;
3456  while (g_idx >= (unsigned int)parallel_ghost->tree_offsets[ghost_tree+1])
3457  ++ghost_tree;
3458 
3459  quadr = static_cast<typename ::internal::p4est::types<dim>::quadrant *>
3460  ( sc_array_index(&parallel_ghost->ghosts, g_idx) );
3461 
3462  unsigned int coarse_cell_index =
3463  p4est_tree_to_coarse_cell_permutation[ghost_tree];
3464 
3465  match_quadrant<dim,spacedim> (this, coarse_cell_index, *quadr, ghost_owner);
3466  }
3467 
3468  // fix all the flags to make sure we have a consistent mesh
3470 
3471  // see if any flags are still set
3472  mesh_changed = false;
3473  for (typename Triangulation<dim,spacedim>::active_cell_iterator
3474  cell = this->begin_active();
3475  cell != this->end();
3476  ++cell)
3477  if (cell->refine_flag_set() || cell->coarsen_flag_set())
3478  {
3479  mesh_changed = true;
3480  break;
3481  }
3482 
3483  // actually do the refinement but prevent the refinement hook below
3484  // from taking over
3485  const bool saved_refinement_in_progress = refinement_in_progress;
3486  refinement_in_progress = true;
3487 
3488  try
3489  {
3491  }
3492  catch (const typename Triangulation<dim,spacedim>::DistortedCellList &)
3493  {
3494  // the underlying triangulation should not be checking for
3495  // distorted cells
3496  AssertThrow (false, ExcInternalError());
3497  }
3498 
3499  refinement_in_progress = saved_refinement_in_progress;
3500  }
3501  while (mesh_changed);
3502 
3503 #ifdef DEBUG
3504  // check if correct number of ghosts is created
3505  unsigned int num_ghosts = 0;
3506 
3507  for (typename Triangulation<dim,spacedim>::active_cell_iterator
3508  cell = this->begin_active();
3509  cell != this->end();
3510  ++cell)
3511  {
3512  if (cell->subdomain_id() != this->my_subdomain
3513  &&
3514  cell->subdomain_id() != numbers::artificial_subdomain_id)
3515  ++num_ghosts;
3516  }
3517 
3518  Assert( num_ghosts == parallel_ghost->ghosts.elem_count, ExcInternalError());
3519 #endif
3520 
3521 
3522 
3523  // fill level_subdomain_ids for geometric multigrid
3524  // the level ownership of a cell is defined as the owner if the cell is active or as the owner of child(0)
3525  // we need this information for all our ancestors and the same-level neighbors of our own cells (=level ghosts)
3526  if (settings & construct_multigrid_hierarchy)
3527  {
3528  // step 1: We set our own ids all the way down and all the others to
3529  // -1. Note that we do not fill other cells we could figure out the
3530  // same way, because we might accidentally set an id for a cell that
3531  // is not a ghost cell.
3532  for (unsigned int lvl=this->n_levels(); lvl>0; )
3533  {
3534  --lvl;
3535  for (typename Triangulation<dim,spacedim>::cell_iterator cell = this->begin(lvl); cell!=this->end(lvl); ++cell)
3536  {
3537  if ((!cell->has_children() && cell->subdomain_id()==this->locally_owned_subdomain())
3538  || (cell->has_children() && cell->child(0)->level_subdomain_id()==this->locally_owned_subdomain()))
3539  cell->set_level_subdomain_id(this->locally_owned_subdomain());
3540  else
3541  {
3542  //not our cell
3543  cell->set_level_subdomain_id(numbers::artificial_subdomain_id);
3544  }
3545  }
3546  }
3547 
3548  //step 2: make sure all the neighbors to our level_cells exist. Need
3549  //to look up in p4est...
3550  std::vector<std::vector<bool> > marked_vertices(this->n_levels());
3551  for (unsigned int lvl=0; lvl < this->n_levels(); ++lvl)
3552  marked_vertices[lvl] = mark_locally_active_vertices_on_level(lvl);
3553 
3554  for (typename Triangulation<dim,spacedim>::cell_iterator cell = this->begin(0); cell!=this->end(0); ++cell)
3555  {
3556  typename ::internal::p4est::types<dim>::quadrant p4est_coarse_cell;
3557  const unsigned int tree_index
3558  = coarse_cell_to_p4est_tree_permutation[cell->index()];
3559  typename ::internal::p4est::types<dim>::tree *tree =
3560  init_tree(cell->index());
3561 
3562  ::internal::p4est::init_coarse_quadrant<dim>(p4est_coarse_cell);
3563 
3564  determine_level_subdomain_id_recursively<dim,spacedim> (*tree, tree_index, cell,
3565  p4est_coarse_cell,
3566  *parallel_forest,
3567  this->my_subdomain,
3568  marked_vertices);
3569  }
3570 
3571  //step 3: make sure we have the parent of our level cells
3572  for (unsigned int lvl=this->n_levels(); lvl>0;)
3573  {
3574  --lvl;
3575  for (typename Triangulation<dim,spacedim>::cell_iterator cell = this->begin(lvl); cell!=this->end(lvl); ++cell)
3576  {
3577  if (cell->has_children())
3578  for (unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
3579  {
3580  if (cell->child(c)->level_subdomain_id()==this->locally_owned_subdomain())
3581  {
3582  //at least one of the children belongs to us, so
3583  //make sure we set the level subdomain id
3585  mark = cell->child(0)->level_subdomain_id();
3586  Assert(mark != numbers::artificial_subdomain_id, ExcInternalError()); //we should know the child(0)
3587  cell->set_level_subdomain_id(mark);
3588  break;
3589  }
3590  }
3591  }
3592  }
3593 
3594  //step 4: Special case: on each level we need all the face neighbors
3595  // of our own level cells these are normally on the same level,
3596  // unless the neighbor is active and coarser. It can end up on a
3597  // different processor. Luckily, the level_subdomain_id can be
3598  // figured out without communication, because the cell is active
3599  // (and so level_subdomain_id=subdomain_id). Finally, also consider
3600  // the opposite case: if we are the coarser neighbor for another
3601  // processor, also mark them.
3602  for (typename Triangulation<dim,spacedim>::cell_iterator cell = this->begin(); cell!=this->end(); ++cell)
3603  {
3604  bool cell_level_mine = cell->level_subdomain_id() == this->locally_owned_subdomain();
3605 
3606  for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
3607  {
3608  if (cell->face(f)->at_boundary() || cell->neighbor(f)->level() >= cell->level())
3609  continue;
3610 
3611  bool neighbor_level_mine = cell->neighbor(f)->level_subdomain_id() == this->locally_owned_subdomain();
3612 
3613  if (cell_level_mine && !neighbor_level_mine)
3614  {
3615  // set the neighbor level_id up
3616  Assert(cell->neighbor(f)->active(), ExcInternalError());
3617  Assert(cell->neighbor(f)->subdomain_id() != numbers::artificial_subdomain_id, ExcInternalError());
3618  Assert(cell->neighbor(f)->level_subdomain_id() == numbers::artificial_subdomain_id
3619  || cell->neighbor(f)->level_subdomain_id() == cell->neighbor(f)->subdomain_id(), ExcInternalError());
3620  cell->neighbor(f)->set_level_subdomain_id(cell->neighbor(f)->subdomain_id());
3621  }
3622  else if (!cell_level_mine && neighbor_level_mine)
3623  {
3624  // set the current cell up because it is a neighbor for us
3625  Assert(cell->active(), ExcInternalError());
3626  Assert(cell->subdomain_id() != numbers::artificial_subdomain_id, ExcInternalError());
3627  Assert(cell->level_subdomain_id() == numbers::artificial_subdomain_id
3628  || cell->level_subdomain_id() == cell->subdomain_id(), ExcInternalError());
3629  cell->set_level_subdomain_id(cell->subdomain_id());
3630  }
3631  }
3632 
3633  }
3634 
3635  }
3636 
3637 
3638 
3639  // check that our local copy has exactly as many cells as the p4est
3640  // original (at least if we are on only one processor); for parallel
3641  // computations, we want to check that we have at least as many as p4est
3642  // stores locally (in the future we should check that we have exactly as
3643  // many non-artificial cells as parallel_forest->local_num_quadrants)
3644  {
3645  const unsigned int total_local_cells = this->n_active_cells();
3646  (void)total_local_cells;
3647 
3649  Assert (static_cast<unsigned int>(parallel_forest->local_num_quadrants) ==
3650  total_local_cells,
3651  ExcInternalError())
3652  else
3653  Assert (static_cast<unsigned int>(parallel_forest->local_num_quadrants) <=
3654  total_local_cells,
3655  ExcInternalError());
3656 
3657  // count the number of owned, active cells and compare with p4est.
3658  unsigned int n_owned = 0;
3659  for (typename Triangulation<dim,spacedim>::active_cell_iterator
3660  cell = this->begin_active();
3661  cell != this->end(); ++cell)
3662  {
3663  if (cell->subdomain_id() == this->my_subdomain)
3664  ++n_owned;
3665  }
3666 
3667  Assert(static_cast<unsigned int>(parallel_forest->local_num_quadrants) ==
3668  n_owned, ExcInternalError());
3669 
3670  }
3671 
3672  this->smooth_grid = save_smooth;
3673  }
3674 
3675 
3676 
3677  template <int dim, int spacedim>
3678  void
3680  {
3681  // first make sure that recursive calls are handled correctly
3682  if (refinement_in_progress == true)
3683  {
3685  return;
3686  }
3687 
3688  // do not allow anisotropic refinement
3689 #ifdef DEBUG
3690  for (typename Triangulation<dim,spacedim>::active_cell_iterator
3691  cell = this->begin_active();
3692  cell != this->end(); ++cell)
3693  if (cell->is_locally_owned() && cell->refine_flag_set())
3694  Assert (cell->refine_flag_set() ==
3696  ExcMessage ("This class does not support anisotropic refinement"));
3697 #endif
3698 
3699 
3700  // safety check: p4est has an upper limit on the level of a cell
3702  {
3703  for (typename Triangulation<dim,spacedim>::active_cell_iterator
3705  cell != this->end(::internal::p4est::functions<dim>::max_level-1); ++cell)
3706  {
3707  AssertThrow(!(cell->refine_flag_set()),
3708  ExcMessage("Fatal Error: maximum refinement level of p4est reached."));
3709  }
3710  }
3711 
3712  // now do the work we're supposed to do when we are in charge
3713  refinement_in_progress = true;
3715 
3716  // make sure all flags are cleared on cells we don't own, since nothing
3717  // good can come of that if they are still around
3718  for (typename Triangulation<dim,spacedim>::active_cell_iterator
3719  cell = this->begin_active();
3720  cell != this->end(); ++cell)
3721  if (cell->is_ghost() || cell->is_artificial())
3722  {
3723  cell->clear_refine_flag ();
3724  cell->clear_coarsen_flag ();
3725  }
3726 
3727 
3728  // count how many cells will be refined and coarsened, and allocate that
3729  // much memory
3730  RefineAndCoarsenList<dim,spacedim>
3731  refine_and_coarsen_list (*this,
3732  p4est_tree_to_coarse_cell_permutation,
3733  this->my_subdomain);
3734 
3735  // copy refine and coarsen flags into p4est and execute the refinement
3736  // and coarsening. this uses the refine_and_coarsen_list just built,
3737  // which is communicated to the callback functions through
3738  // p4est's user_pointer object
3739  Assert (parallel_forest->user_pointer == this,
3740  ExcInternalError());
3741  parallel_forest->user_pointer = &refine_and_coarsen_list;
3742 
3743  if (parallel_ghost != 0)
3744  {
3746  parallel_ghost = 0;
3747  }
3749  refine (parallel_forest, /* refine_recursive */ false,
3750  &RefineAndCoarsenList<dim,spacedim>::refine_callback,
3751  /*init_callback=*/NULL);
3753  coarsen (parallel_forest, /* coarsen_recursive */ false,
3754  &RefineAndCoarsenList<dim,spacedim>::coarsen_callback,
3755  /*init_callback=*/NULL);
3756 
3757  // make sure all cells in the lists have been consumed
3758  Assert (refine_and_coarsen_list.pointers_are_at_end(),
3759  ExcInternalError());
3760 
3761  // reset the pointer
3762  parallel_forest->user_pointer = this;
3763 
3764  // enforce 2:1 hanging node condition
3766  balance (parallel_forest,
3767  /* face and corner balance */
3768  (dim == 2
3769  ?
3770  typename ::internal::p4est::types<dim>::
3771  balance_type(P4EST_CONNECT_FULL)
3772  :
3773  typename ::internal::p4est::types<dim>::
3774  balance_type(P8EST_CONNECT_FULL)),
3775  /*init_callback=*/NULL);
3776 
3777  // before repartitioning the mesh let others attach mesh related info
3778  // (such as SolutionTransfer data) to the p4est
3779  attach_mesh_data();
3780 
3782  {
3783  // partition the new mesh between all processors. If cell weights have
3784  // not been given balance the number of cells.
3785  if (this->signals.cell_weight.num_slots() == 0)
3787  partition (parallel_forest,
3788  /* prepare coarsening */ 1,
3789  /* weight_callback */ NULL);
3790  else
3791  {
3792  // get cell weights for a weighted repartitioning.
3793  const std::vector<unsigned int> cell_weights = get_cell_weights();
3794 
3795  PartitionWeights<dim,spacedim> partition_weights (cell_weights);
3796 
3797  // attach (temporarily) a pointer to the cell weights through p4est's
3798  // user_pointer object
3799  Assert (parallel_forest->user_pointer == this,
3800  ExcInternalError());
3801  parallel_forest->user_pointer = &partition_weights;
3802 
3804  partition (parallel_forest,
3805  /* prepare coarsening */ 1,
3806  /* weight_callback */ &PartitionWeights<dim,spacedim>::cell_weight);
3807 
3808  // reset the user pointer to its previous state
3809  parallel_forest->user_pointer = this;
3810  }
3811  }
3812 
3813  // finally copy back from local part of tree to deal.II
3814  // triangulation. before doing so, make sure there are no refine or
3815  // coarsen flags pending
3816  for (typename Triangulation<dim,spacedim>::active_cell_iterator
3817  cell = this->begin_active();
3818  cell != this->end(); ++cell)
3819  {
3820  cell->clear_refine_flag();
3821  cell->clear_coarsen_flag();
3822  }
3823 
3824  try
3825  {
3827  }
3828  catch (const typename Triangulation<dim>::DistortedCellList &)
3829  {
3830  // the underlying triangulation should not be checking for distorted
3831  // cells
3832  AssertThrow (false, ExcInternalError());
3833  }
3834 
3835 
3836  refinement_in_progress = false;
3837 
3838  this->update_number_cache ();
3839  }
3840 
3841  template <int dim, int spacedim>
3842  void
3844  {
3845 
3846 #ifdef DEBUG
3847  for (typename Triangulation<dim,spacedim>::active_cell_iterator
3848  cell = this->begin_active();
3849  cell != this->end(); ++cell)
3850  if (cell->is_locally_owned())
3851  Assert (
3852  !cell->refine_flag_set() && !cell->coarsen_flag_set(),
3853  ExcMessage ("Error: There shouldn't be any cells flagged for coarsening/refinement when calling repartition()."));
3854 #endif
3855 
3856  refinement_in_progress = true;
3857 
3858  // before repartitioning the mesh let others attach mesh related info
3859  // (such as SolutionTransfer data) to the p4est
3860  attach_mesh_data();
3861 
3862  if (this->signals.cell_weight.num_slots() == 0)
3863  {
3864  // no cell weights given -- call p4est's 'partition' without a
3865  // callback for cell weights
3867  partition (parallel_forest,
3868  /* prepare coarsening */ 1,
3869  /* weight_callback */ NULL);
3870  }
3871  else
3872  {
3873  // get cell weights for a weighted repartitioning.
3874  const std::vector<unsigned int> cell_weights = get_cell_weights();
3875 
3876  PartitionWeights<dim,spacedim> partition_weights (cell_weights);
3877 
3878  // attach (temporarily) a pointer to the cell weights through p4est's
3879  // user_pointer object
3880  Assert (parallel_forest->user_pointer == this,
3881  ExcInternalError());
3882  parallel_forest->user_pointer = &partition_weights;
3883 
3885  partition (parallel_forest,
3886  /* prepare coarsening */ 1,
3887  /* weight_callback */ &PartitionWeights<dim,spacedim>::cell_weight);
3888 
3889  // reset the user pointer to its previous state
3890  parallel_forest->user_pointer = this;
3891  }
3892 
3893  try
3894  {
3896  }
3897  catch (const typename Triangulation<dim>::DistortedCellList &)
3898  {
3899  // the underlying triangulation should not be checking for distorted
3900  // cells
3901  AssertThrow (false, ExcInternalError());
3902  }
3903 
3904  refinement_in_progress = false;
3905 
3906  // update how many cells, edges, etc, we store locally
3907  this->update_number_cache ();
3908  }
3909 
3910 
3911 
3912  template <int dim, int spacedim>
3913  void
3915  communicate_locally_moved_vertices (const std::vector<bool> &vertex_locally_moved)
3916  {
3917  Assert (vertex_locally_moved.size() == this->n_vertices(),
3918  ExcDimensionMismatch(vertex_locally_moved.size(),
3919  this->n_vertices()));
3920 #ifdef DEBUG
3921  {
3922  const std::vector<bool> locally_owned_vertices
3924  for (unsigned int i=0; i<locally_owned_vertices.size(); ++i)
3925  Assert ((vertex_locally_moved[i] == false)
3926  ||
3927  (locally_owned_vertices[i] == true),
3928  ExcMessage ("The vertex_locally_moved argument must not "
3929  "contain vertices that are not locally owned"));
3930  }
3931 #endif
3932 
3933  std::map<unsigned int, std::set<::types::subdomain_id> >
3934  vertices_with_ghost_neighbors;
3935 
3936  // First find out which process should receive which vertices.
3937  // these are specifically the ones that sit on ghost cells and,
3938  // among these, the ones that we own locally
3939  for (typename Triangulation<dim,spacedim>::active_cell_iterator
3940  cell=this->begin_active(); cell!=this->end();
3941  ++cell)
3942  if (cell->is_ghost())
3943  for (unsigned int vertex_no=0;
3944  vertex_no<GeometryInfo<dim>::vertices_per_cell; ++vertex_no)
3945  {
3946  const unsigned int process_local_vertex_no = cell->vertex_index(vertex_no);
3947  vertices_with_ghost_neighbors[process_local_vertex_no].insert
3948  (cell->subdomain_id());
3949  }
3950 
3951  // now collect cells and their vertices
3952  // for the interested neighbors
3953  typedef
3954  std::map<::types::subdomain_id, CommunicateLocallyMovedVertices::CellInfo<dim,spacedim> > cellmap_t;
3955  cellmap_t needs_to_get_cells;
3956 
3958  cell = this->begin(0);
3959  cell != this->end(0);
3960  ++cell)
3961  {
3962  typename ::internal::p4est::types<dim>::quadrant p4est_coarse_cell;
3963  ::internal::p4est::init_coarse_quadrant<dim>(p4est_coarse_cell);
3964 
3965  CommunicateLocallyMovedVertices::fill_vertices_recursively<dim,spacedim>
3966  (*this,
3967  this->get_coarse_cell_to_p4est_tree_permutation()[cell->index()],
3968  cell,
3969  p4est_coarse_cell,
3970  vertices_with_ghost_neighbors,
3971  vertex_locally_moved,
3972  needs_to_get_cells);
3973  }
3974 
3975  // sending
3976  std::vector<std::vector<char> > sendbuffers (needs_to_get_cells.size());
3977  std::vector<std::vector<char> >::iterator buffer = sendbuffers.begin();
3978  std::vector<MPI_Request> requests (needs_to_get_cells.size());
3979  std::vector<unsigned int> destinations;
3980 
3981  unsigned int idx=0;
3982 
3983  for (typename cellmap_t::iterator it=needs_to_get_cells.begin();
3984  it!=needs_to_get_cells.end();
3985  ++it, ++buffer, ++idx)
3986  {
3987  const unsigned int num_cells = it->second.tree_index.size();
3988  (void)num_cells;
3989  destinations.push_back(it->first);
3990 
3991  Assert(num_cells==it->second.quadrants.size(), ExcInternalError());
3992  Assert(num_cells>0, ExcInternalError());
3993 
3994  // pack all the data into
3995  // the buffer for this
3996  // recipient and send
3997  // it. keep data around
3998  // till we can make sure
3999  // that the packet has been
4000  // received
4001  it->second.pack_data (*buffer);
4002  MPI_Isend(&(*buffer)[0], buffer->size(),
4003  MPI_BYTE, it->first,
4004  123, this->get_communicator(), &requests[idx]);
4005  }
4006 
4007  Assert(destinations.size()==needs_to_get_cells.size(), ExcInternalError());
4008 
4009  // collect the neighbors
4010  // that are going to send stuff to us
4011  const std::vector<unsigned int> senders
4013  (this->get_communicator(), destinations);
4014 
4015  // receive ghostcelldata
4016  std::vector<char> receive;
4017  CommunicateLocallyMovedVertices::CellInfo<dim,spacedim> cellinfo;
4018  for (unsigned int i=0; i<senders.size(); ++i)
4019  {
4020  MPI_Status status;
4021  int len;
4022  MPI_Probe(MPI_ANY_SOURCE, 123, this->get_communicator(), &status);
4023  MPI_Get_count(&status, MPI_BYTE, &len);
4024  receive.resize(len);
4025 
4026  char *ptr = &receive[0];
4027  MPI_Recv(ptr, len, MPI_BYTE, status.MPI_SOURCE, status.MPI_TAG,
4028  this->get_communicator(), &status);
4029 
4030  cellinfo.unpack_data(receive);
4031  const unsigned int cells = cellinfo.tree_index.size();
4032  for (unsigned int c=0; c<cells; ++c)
4033  {
4034  typename ::parallel::distributed::Triangulation<dim,spacedim>::cell_iterator
4035  cell (this,
4036  0,
4037  this->get_p4est_tree_to_coarse_cell_permutation()[cellinfo.tree_index[c]]);
4038 
4039  typename ::internal::p4est::types<dim>::quadrant p4est_coarse_cell;
4040  ::internal::p4est::init_coarse_quadrant<dim>(p4est_coarse_cell);
4041 
4042  CommunicateLocallyMovedVertices::set_vertices_recursively<dim,spacedim> (*this,
4043  p4est_coarse_cell,
4044  cell,
4045  cellinfo.quadrants[c],
4046  cellinfo.first_vertices[c],
4047  cellinfo.first_vertex_indices[c]);
4048  }
4049  }
4050 
4051  // complete all sends, so that we can
4052  // safely destroy the buffers.
4053  if (requests.size() > 0)
4054  MPI_Waitall(requests.size(), &requests[0], MPI_STATUSES_IGNORE);
4055 
4056  //check all msgs got sent and received
4057  Assert(Utilities::MPI::sum(needs_to_get_cells.size(), this->get_communicator())
4058  == Utilities::MPI::sum(senders.size(), this->get_communicator()),
4059  ExcInternalError());
4060  }
4061 
4062  template <int dim, int spacedim>
4063  unsigned int
4065  register_data_attach (const std::size_t size,
4066  const std_cxx11::function<void(const cell_iterator &,
4067  const CellStatus,
4068  void *)> &pack_callback)
4069  {
4070  Assert(size>0, ExcMessage("register_data_attach(), size==0"));
4072  ExcMessage("register_data_attach(), not all data has been unpacked last time?"));
4073 
4074  unsigned int offset = attached_data_size+sizeof(CellStatus);
4075  ++n_attached_datas;
4076  attached_data_size+=size;
4077  attached_data_pack_callbacks.push_back(
4078  std::pair<unsigned int, pack_callback_t> (offset, pack_callback)
4079  );
4080  return offset;
4081  }
4082 
4083 
4084 
4085  template <int dim, int spacedim>
4086  void
4088  notify_ready_to_unpack (const unsigned int offset,
4089  const std_cxx11::function<void (const cell_iterator &,
4090  const CellStatus,
4091  const void *)> &unpack_callback)
4092  {
4093  Assert (offset >= sizeof(CellStatus),
4094  ExcMessage ("invalid offset in notify_ready_to_unpack()"));
4095  Assert (offset < sizeof(CellStatus)+attached_data_size,
4096  ExcMessage ("invalid offset in notify_ready_to_unpack()"));
4097  Assert (n_attached_datas > 0, ExcMessage ("notify_ready_to_unpack() called too often"));
4098 
4099  // Recurse over p4est and hand the caller the data back
4101  cell = this->begin (0);
4102  cell != this->end (0);
4103  ++cell)
4104  {
4105  //skip coarse cells, that are not ours
4106  if (tree_exists_locally<dim, spacedim> (parallel_forest,
4107  coarse_cell_to_p4est_tree_permutation[cell->index() ])
4108  == false)
4109  continue;
4110 
4111  typename ::internal::p4est::types<dim>::quadrant p4est_coarse_cell;
4112  typename ::internal::p4est::types<dim>::tree *tree =
4113  init_tree (cell->index());
4114 
4115  ::internal::p4est::init_coarse_quadrant<dim> (p4est_coarse_cell);
4116 
4117  // parent_cell is not correct here, but is only used in a refined
4118  // cell
4119  post_mesh_data_recursively<dim, spacedim> (*tree,
4120  cell,
4121  cell,
4122  p4est_coarse_cell,
4123  offset,
4124  unpack_callback);
4125  }
4126 
4127  --n_attached_datas;
4128  if (n_attached_deserialize > 0)
4129  {
4131  attached_data_pack_callbacks.pop_front();
4132  }
4133 
4134  // important: only remove data if we are not in the deserialization
4135  // process. There, each SolutionTransfer registers and unpacks before
4136  // the next one does this, so n_attached_datas is only 1 here. This
4137  // would destroy the saved data before the second SolutionTransfer can
4138  // get it. This created a bug that is documented in
4139  // tests/mpi/p4est_save_03 with more than one SolutionTransfer.
4140  if (!n_attached_datas && n_attached_deserialize == 0)
4141  {
4142  // everybody got his data, time for cleanup!
4143  attached_data_size = 0;
4145 
4146  // and release the data
4147  void *userptr = parallel_forest->user_pointer;
4148  ::internal::p4est::functions<dim>::reset_data (parallel_forest, 0, NULL, NULL);
4149  parallel_forest->user_pointer = userptr;
4150  }
4151  }
4152 
4153 
4154  template <int dim, int spacedim>
4155  const std::vector<types::global_dof_index> &
4157  {
4158  return p4est_tree_to_coarse_cell_permutation;
4159  }
4160 
4161 
4162 
4163  template <int dim, int spacedim>
4164  const std::vector<types::global_dof_index> &
4166  {
4168  }
4169 
4170 
4171 
4172  namespace
4173  {
4178  template <int dim, int spacedim>
4179  struct FindGhosts
4180  {
4181  typename ::parallel::distributed::Triangulation<dim,spacedim> *triangulation;
4182  sc_array_t *subids;
4183  std::map<unsigned int, std::set<::types::subdomain_id> >
4184  *vertices_with_ghost_neighbors;
4185  };
4186 
4192  template <int dim, int spacedim>
4193  void
4194  find_ghosts_corner
4195  (typename ::internal::p4est::iter<dim>::corner_info *info,
4196  void *user_data)
4197  {
4198  int i, j;
4199  int nsides = info->sides.elem_count;
4200  typename ::internal::p4est::iter<dim>::corner_side *sides =
4201  (typename ::internal::p4est::iter<dim>::corner_side *)
4202  (info->sides.array);
4203  FindGhosts<dim,spacedim> *fg = static_cast<FindGhosts<dim,spacedim> *>(user_data);
4204  sc_array_t *subids = fg->subids;
4205  typename ::parallel::distributed::Triangulation<dim,spacedim> *triangulation = fg->triangulation;
4206  int nsubs;
4207  ::types::subdomain_id *subdomain_ids;
4208  std::map<unsigned int, std::set<::types::subdomain_id> >
4209  *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
4210 
4211  subids->elem_count = 0;
4212  for (i = 0; i < nsides; i++)
4213  {
4214  if (sides[i].is_ghost)
4215  {
4216  typename ::parallel::distributed::Triangulation<dim,spacedim>::cell_iterator cell = cell_from_quad (triangulation, sides[i].treeid, *(sides[i].quad));
4217  Assert (cell->is_ghost(), ExcMessage ("ghost quad did not find ghost cell"));
4218  ::types::subdomain_id *subid =
4219  static_cast<::types::subdomain_id *>(sc_array_push (subids));
4220  *subid = cell->subdomain_id();
4221  }
4222  }
4223 
4224  if (!subids->elem_count)
4225  {
4226  return;
4227  }
4228 
4229  nsubs = (int) subids->elem_count;
4230  subdomain_ids = (::types::subdomain_id *) (subids->array);
4231 
4232  for (i = 0; i < nsides; i++)
4233  {
4234  if (!sides[i].is_ghost)
4235  {
4236  typename ::parallel::distributed::Triangulation<dim,spacedim>::cell_iterator cell = cell_from_quad (triangulation, sides[i].treeid, *(sides[i].quad));
4237 
4238  Assert (!cell->is_ghost(), ExcMessage ("local quad found ghost cell"));
4239 
4240  for (j = 0; j < nsubs; j++)
4241  {
4242  (*vertices_with_ghost_neighbors)[cell->vertex_index(sides[i].corner)]
4243  .insert (subdomain_ids[j]);
4244  }
4245  }
4246  }
4247 
4248  subids->elem_count = 0;
4249  }
4250 
4254  template <int dim, int spacedim>
4255  void
4256  find_ghosts_edge
4257  (typename ::internal::p4est::iter<dim>::edge_info *info,
4258  void *user_data)
4259  {
4260  int i, j, k;
4261  int nsides = info->sides.elem_count;
4262  typename ::internal::p4est::iter<dim>::edge_side *sides =
4263  (typename ::internal::p4est::iter<dim>::edge_side *)
4264  (info->sides.array);
4265  FindGhosts<dim,spacedim> *fg = static_cast<FindGhosts<dim,spacedim> *>(user_data);
4266  sc_array_t *subids = fg->subids;
4267  typename ::parallel::distributed::Triangulation<dim,spacedim> *triangulation = fg->triangulation;
4268  int nsubs;
4269  ::types::subdomain_id *subdomain_ids;
4270  std::map<unsigned int, std::set<::types::subdomain_id> >
4271  *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
4272 
4273  subids->elem_count = 0;
4274  for (i = 0; i < nsides; i++)
4275  {
4276  if (sides[i].is_hanging)
4277  {
4278  for (j = 0; j < 2; j++)
4279  {
4280  if (sides[i].is.hanging.is_ghost[j])
4281  {
4282  typename ::parallel::distributed::Triangulation<dim,spacedim>::cell_iterator cell = cell_from_quad (triangulation, sides[i].treeid, *(sides[i].is.hanging.quad[j]));
4283  ::types::subdomain_id *subid =
4284  static_cast<::types::subdomain_id *>(sc_array_push (subids));
4285  *subid = cell->subdomain_id();
4286  }
4287  }
4288  }
4289  }
4290 
4291  if (!subids->elem_count)
4292  {
4293  return;
4294  }
4295 
4296  nsubs = (int) subids->elem_count;
4297  subdomain_ids = (::types::subdomain_id *) (subids->array);
4298 
4299  for (i = 0; i < nsides; i++)
4300  {
4301  if (sides[i].is_hanging)
4302  {
4303  for (j = 0; j < 2; j++)
4304  {
4305  if (!sides[i].is.hanging.is_ghost[j])
4306  {
4307  typename ::parallel::distributed::Triangulation<dim,spacedim>::cell_iterator cell = cell_from_quad (triangulation, sides[i].treeid, *(sides[i].is.hanging.quad[j]));
4308 
4309  for (k = 0; k < nsubs; k++)
4310  {
4311  (*vertices_with_ghost_neighbors)[cell->vertex_index(p8est_edge_corners[sides[i].edge][1^j])]
4312  .insert (subdomain_ids[k]);
4313  }
4314  }
4315  }
4316  }
4317  }
4318 
4319  subids->elem_count = 0;
4320  }
4321 
4325  template <int dim, int spacedim>
4326  void
4327  find_ghosts_face
4328  (typename ::internal::p4est::iter<dim>::face_info *info,
4329  void *user_data)
4330  {
4331  int i, j, k;
4332  int nsides = info->sides.elem_count;
4333  typename ::internal::p4est::iter<dim>::face_side *sides =
4334  (typename ::internal::p4est::iter<dim>::face_side *)
4335  (info->sides.array);
4336  FindGhosts<dim,spacedim> *fg = static_cast<FindGhosts<dim,spacedim> *>(user_data);
4337  sc_array_t *subids = fg->subids;
4338  typename ::parallel::distributed::Triangulation<dim,spacedim> *triangulation = fg->triangulation;
4339  int nsubs;
4340  ::types::subdomain_id *subdomain_ids;
4341  std::map<unsigned int, std::set<::types::subdomain_id> >
4342  *vertices_with_ghost_neighbors = fg->vertices_with_ghost_neighbors;
4343  int limit = (dim == 2) ? 2 : 4;
4344 
4345  subids->elem_count = 0;
4346  for (i = 0; i < nsides; i++)
4347  {
4348  if (sides[i].is_hanging)
4349  {
4350  for (j = 0; j < limit; j++)
4351  {
4352  if (sides[i].is.hanging.is_ghost[j])
4353  {
4354  typename ::parallel::distributed::Triangulation<dim,spacedim>::cell_iterator cell = cell_from_quad (triangulation, sides[i].treeid, *(sides[i].is.hanging.quad[j]));
4355  ::types::subdomain_id *subid =
4356  static_cast<::types::subdomain_id *>(sc_array_push (subids));
4357  *subid = cell->subdomain_id();
4358  }
4359  }
4360  }
4361  }
4362 
4363  if (!subids->elem_count)
4364  {
4365  return;
4366  }
4367 
4368  nsubs = (int) subids->elem_count;
4369  subdomain_ids = (::types::subdomain_id *) (subids->array);
4370 
4371  for (i = 0; i < nsides; i++)
4372  {
4373  if (sides[i].is_hanging)
4374  {
4375  for (j = 0; j < limit; j++)
4376  {
4377  if (!sides[i].is.hanging.is_ghost[j])
4378  {
4379  typename ::parallel::distributed::Triangulation<dim,spacedim>::cell_iterator cell = cell_from_quad (triangulation, sides[i].treeid, *(sides[i].is.hanging.quad[j]));
4380 
4381  for (k = 0; k < nsubs; k++)
4382  {
4383  if (dim == 2)
4384  {
4385  (*vertices_with_ghost_neighbors)[cell->vertex_index(p4est_face_corners[sides[i].face][(limit - 1)^j])]
4386  .insert (subdomain_ids[k]);
4387  }
4388  else
4389  {
4390  (*vertices_with_ghost_neighbors)[cell->vertex_index(p8est_face_corners[sides[i].face][(limit - 1)^j])]
4391  .insert (subdomain_ids[k]);
4392  }
4393  }
4394  }
4395  }
4396  }
4397  }
4398 
4399  subids->elem_count = 0;
4400  }
4401  }
4402 
4403 
4404 
4405  namespace
4406  {
4412  template <typename ITERATOR>
4413  void
4414  mark_periodic_vertices_recursively(const GridTools::PeriodicFacePair<ITERATOR> &periodic,
4415  const int target_level,
4416  std::vector<bool> &active_vertices_on_level)
4417  {
4418  if (periodic.cell[0]->level() > target_level)
4419  return;
4420 
4421  const unsigned int dim = ITERATOR::AccessorType::dimension;
4422  // for hanging nodes there is nothing to do since we are interested in
4423  // the connections on the same level...
4424  if (periodic.cell[0]->level() < target_level &&
4425  periodic.cell[0]->has_children() &&
4426  periodic.cell[1]->has_children())
4427  {
4428  // copy orientations etc. from parent to child
4429  GridTools::PeriodicFacePair<ITERATOR> periodic_child = periodic;
4430 
4431  // find appropriate pairs of child elements
4432  for (unsigned int cf=0; cf<periodic.cell[0]->face(periodic.face_idx[0])->n_children(); ++cf)
4433  {
4434  const unsigned int child_index_0 =
4435  GeometryInfo<dim>::child_cell_on_face(periodic.cell[0]->refinement_case(),
4436  periodic.face_idx[0], cf,
4437  periodic.orientation[0],
4438  periodic.orientation[1],
4439  periodic.orientation[2]);
4440  periodic_child.cell[0] = periodic.cell[0]->child(child_index_0);
4441  periodic_child.face_idx[0] = periodic.face_idx[0];
4442 
4443  // the second face is in standard orientation in terms of the
4444  // periodic face pair
4445  const unsigned int child_index_1 =
4446  GeometryInfo<dim>::child_cell_on_face(periodic.cell[1]->refinement_case(),
4447  periodic.face_idx[1], cf);
4448 
4449  periodic_child.cell[1] = periodic.cell[1]->child(child_index_1);
4450  periodic_child.face_idx[1] = periodic.face_idx[1];
4451 
4452  // recursive call into children
4453  mark_periodic_vertices_recursively (periodic_child, target_level,
4454  active_vertices_on_level);
4455  }
4456  return;
4457  }
4458 
4459  if (periodic.cell[0]->level() != target_level)
4460  return;
4461 
4462  for (unsigned int v=0; v<GeometryInfo<dim-1>::vertices_per_cell; ++v)
4463  {
4464  // take possible non-standard orientation of face on cell[0] into
4465  // account
4466  const unsigned int vface0 =
4468  periodic.orientation[1],
4469  periodic.orientation[2]);
4470  if (active_vertices_on_level[periodic.cell[0]->face(periodic.face_idx[0])->vertex_index(vface0)] ||
4471  active_vertices_on_level[periodic.cell[1]->face(periodic.face_idx[1])->vertex_index(v)])
4472  active_vertices_on_level[periodic.cell[0]->face(periodic.face_idx[0])->vertex_index(vface0)]
4473  = active_vertices_on_level[periodic.cell[1]->face(periodic.face_idx[1])->vertex_index(v)]
4474  = true;
4475  }
4476  }
4477 
4478 
4479 
4484  template <typename ITERATOR>
4485  void
4486  set_periodic_ghost_neighbors_recursively(const GridTools::PeriodicFacePair<ITERATOR> &periodic,
4487  const int target_level,
4488  std::map<unsigned int, std::set<::types::subdomain_id> > &vertices_with_ghost_neighbors)
4489  {
4490  if (periodic.cell[0]->level() > target_level)
4491  return;
4492 
4493  // for hanging nodes there is nothing to do since we are interested in
4494  // the connections on the same level...
4495  if (periodic.cell[0]->level() < target_level &&
4496  periodic.cell[0]->has_children() &&
4497  periodic.cell[1]->has_children())
4498  {
4499  // copy orientations etc. from parent to child
4500  GridTools::PeriodicFacePair<ITERATOR> periodic_child = periodic;
4501 
4502  // find appropriate pairs of child elements
4503  for (unsigned int cf=0; cf<periodic.cell[0]->face(periodic.face_idx[0])->n_children(); ++cf)
4504  {
4505  unsigned int c=0;
4506  for (; c<periodic.cell[0]->n_children(); ++c)
4507  {
4508  if (periodic.cell[0]->child(c)->face(periodic.face_idx[0]) ==
4509  periodic.cell[0]->face(periodic.face_idx[0])->child(cf))
4510  break;
4511  }
4512  Assert(c < periodic.cell[0]->n_children(),
4513  ExcMessage("Face child not found"));
4514  periodic_child.cell[0] = periodic.cell[0]->child(c);
4515  periodic_child.face_idx[0] = periodic.face_idx[0];
4516 
4517  c=0;
4518  for (; c<periodic.cell[1]->n_children(); ++c)
4519  {
4520  if (periodic.cell[1]->child(c)->face(periodic.face_idx[1]) ==
4521  periodic.cell[1]->face(periodic.face_idx[1])->child(cf))
4522  break;
4523  }
4524  Assert(c < periodic.cell[1]->n_children(),
4525  ExcMessage("Face child not found"));
4526  periodic_child.cell[1] = periodic.cell[1]->child(c);
4527  periodic_child.face_idx[1] = periodic.face_idx[1];
4528 
4529  // recursive call into children
4530  set_periodic_ghost_neighbors_recursively (periodic_child, target_level,
4531  vertices_with_ghost_neighbors);
4532  }
4533  return;
4534  }
4535 
4536  if (periodic.cell[0]->level() != target_level)
4537  return;
4538 
4539  // TODO: fix non-standard orientation
4540  Assert(periodic.orientation[0] == true &&
4541  periodic.orientation[1] == false &&
4542  periodic.orientation[2] == false,
4543  ExcNotImplemented());
4544 
4545  for (unsigned int v=0; v<GeometryInfo<ITERATOR::AccessorType::dimension-1>::vertices_per_cell; ++v)
4546  {
4547  const unsigned int
4548  idx0 = periodic.cell[0]->face(periodic.face_idx[0])->vertex_index(v),
4549  idx1 = periodic.cell[1]->face(periodic.face_idx[1])->vertex_index(v);
4550  if (vertices_with_ghost_neighbors.find(idx0) !=
4551  vertices_with_ghost_neighbors.end())
4552  vertices_with_ghost_neighbors[idx1].
4553  insert(vertices_with_ghost_neighbors[idx0].begin(),
4554  vertices_with_ghost_neighbors[idx0].end());
4555  if (vertices_with_ghost_neighbors.find(idx1) !=
4556  vertices_with_ghost_neighbors.end())
4557  vertices_with_ghost_neighbors[idx0].
4558  insert(vertices_with_ghost_neighbors[idx1].begin(),
4559  vertices_with_ghost_neighbors[idx1].end());
4560  }
4561  }
4562  }
4563 
4564 
4565 
4570  template <int dim, int spacedim>
4571  void
4574  (std::map<unsigned int, std::set<::types::subdomain_id> >
4575  &vertices_with_ghost_neighbors)
4576  {
4577  Assert (dim>1, ExcNotImplemented());
4578 
4579  FindGhosts<dim,spacedim> fg;
4580  fg.subids = sc_array_new (sizeof (::types::subdomain_id));
4581  fg.triangulation = this;
4582  fg.vertices_with_ghost_neighbors = &vertices_with_ghost_neighbors;
4583 
4584  // switch between functions. to make the compiler happy, we need to cast
4585  // the first two arguments to the type p[48]est_iterate wants to see. this
4586  // cast is the identity cast in each of the two branches, so it is safe.
4587  switch (dim)
4588  {
4589  case 2:
4590  p4est_iterate (reinterpret_cast<::internal::p4est::types<2>::forest *>(this->parallel_forest),
4591  reinterpret_cast<::internal::p4est::types<2>::ghost *>(this->parallel_ghost),
4592  static_cast<void *>(&fg),
4593  NULL, find_ghosts_face<2,spacedim>, find_ghosts_corner<2,spacedim>);
4594  break;
4595 
4596  case 3:
4597  p8est_iterate (reinterpret_cast<::internal::p4est::types<3>::forest *>(this->parallel_forest),
4598  reinterpret_cast<::internal::p4est::types<3>::ghost *>(this->parallel_ghost),
4599  static_cast<void *>(&fg),
4600  NULL, find_ghosts_face<3,spacedim>, find_ghosts_edge<3,spacedim>, find_ghosts_corner<3,spacedim>);
4601  break;
4602 
4603  default:
4604  Assert (false, ExcNotImplemented());
4605  }
4606 
4607  sc_array_destroy (fg.subids);
4608  }
4609 
4610 
4611 
4616  template <int dim, int spacedim>
4617  void
4620  (const unsigned int level,
4621  std::map<unsigned int, std::set<::types::subdomain_id> >
4622  &vertices_with_ghost_neighbors)
4623  {
4624  const std::vector<bool> locally_active_vertices =
4626  for (cell_iterator cell = this->begin(level); cell != this->end(level); ++cell)
4627  if (cell->level_subdomain_id() != ::numbers::artificial_subdomain_id
4628  && cell->level_subdomain_id() != this->locally_owned_subdomain())
4629  for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
4630  if (locally_active_vertices[cell->vertex_index(v)])
4631  vertices_with_ghost_neighbors[cell->vertex_index(v)]
4632  .insert (cell->level_subdomain_id());
4633 
4634  for (unsigned int i=0; i<periodic_face_pairs_level_0.size(); ++i)
4635  set_periodic_ghost_neighbors_recursively(periodic_face_pairs_level_0[i],
4636  level, vertices_with_ghost_neighbors);
4637  }
4638 
4639 
4640 
4641  template<int dim, int spacedim>
4642  std::vector<bool>
4644  ::mark_locally_active_vertices_on_level (const unsigned int level) const
4645  {
4646  Assert (dim>1, ExcNotImplemented());
4647 
4648  std::vector<bool> marked_vertices(this->n_vertices(), false);
4649  for (cell_iterator cell = this->begin(level); cell != this->end(level); ++cell)
4650  if (cell->level_subdomain_id() == this->locally_owned_subdomain())
4651  for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
4652  marked_vertices[cell->vertex_index(v)] = true;
4653 
4654  for (unsigned int i=0; i<periodic_face_pairs_level_0.size(); ++i)
4655  mark_periodic_vertices_recursively(periodic_face_pairs_level_0[i],
4656  level, marked_vertices);
4657 
4658  return marked_vertices;
4659  }
4660 
4661 
4662 
4663  template<int dim, int spacedim>
4664  void
4667  periodicity_vector)
4668  {
4669 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,1)
4671  ExcMessage ("The triangulation is empty!"));
4672  Assert (this->n_levels() == 1,
4673  ExcMessage ("The triangulation is refined!"));
4674 
4675  typedef std::vector<GridTools::PeriodicFacePair<cell_iterator> >
4676  FaceVector;
4677  typename FaceVector::const_iterator it, periodic_end;
4678  it = periodicity_vector.begin();
4679  periodic_end = periodicity_vector.end();
4680 
4681  for (; it<periodic_end; ++it)
4682  {
4683  const cell_iterator first_cell = it->cell[0];
4684  const cell_iterator second_cell = it->cell[1];
4685  const unsigned int face_left = it->face_idx[0];
4686  const unsigned int face_right = it->face_idx[1];
4687 
4688  //respective cells of the matching faces in p4est
4689  const unsigned int tree_left
4690  = coarse_cell_to_p4est_tree_permutation[std::distance(this->begin(),
4691  first_cell)];
4692  const unsigned int tree_right
4693  = coarse_cell_to_p4est_tree_permutation[std::distance(this->begin(),
4694  second_cell)];
4695 
4696  // p4est wants to know which corner the first corner on
4697  // the face with the lower id is mapped to on the face with
4698  // with the higher id. For d==2 there are only two possibilities
4699  // that are determined by it->orientation[1].
4700  // For d==3 we have to use GridTools::OrientationLookupTable.
4701  // The result is given below.
4702 
4703  unsigned int p4est_orientation = 0;
4704  if (dim==2)
4705  p4est_orientation = it->orientation[1];
4706  else
4707  {
4708  const unsigned int face_idx_list[] = {face_left, face_right};
4709  const cell_iterator cell_list[] = {first_cell, second_cell};
4710  unsigned int lower_idx, higher_idx;
4711  if (face_left<=face_right)
4712  {
4713  higher_idx = 1;
4714  lower_idx = 0;
4715  }
4716  else
4717  {
4718  higher_idx = 0;
4719  lower_idx = 1;
4720  }
4721 
4722  // get the cell index of the first index on the face with the lower id
4723  unsigned int first_p4est_idx_on_cell = p8est_face_corners[face_idx_list[lower_idx]][0];
4724  unsigned int first_dealii_idx_on_face = numbers::invalid_unsigned_int;
4725  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_face; ++i)
4726  {
4727  const unsigned int first_dealii_idx_on_cell
4729  (face_idx_list[lower_idx], i,
4730  cell_list[lower_idx]->face_orientation(face_idx_list[lower_idx]),
4731  cell_list[lower_idx]->face_flip(face_idx_list[lower_idx]),
4732  cell_list[lower_idx]->face_rotation(face_idx_list[lower_idx]));
4733  if (first_p4est_idx_on_cell == first_dealii_idx_on_cell)
4734  {
4735  first_dealii_idx_on_face = i;
4736  break;
4737  }
4738  }
4739  Assert( first_dealii_idx_on_face != numbers::invalid_unsigned_int, ExcInternalError());
4740  // Now map dealii_idx_on_face according to the orientation
4741  const unsigned int left_to_right [8][4] = {{0,2,1,3},{0,1,2,3},{3,1,2,0},{3,2,1,0},
4742  {2,3,0,1},{1,3,0,2},{1,0,3,2},{2,0,3,1}
4743  };
4744  const unsigned int right_to_left [8][4] = {{0,2,1,3},{0,1,2,3},{3,1,2,0},{3,2,1,0},
4745  {2,3,0,1},{2,0,3,1},{1,0,3,2},{1,3,0,2}
4746  };
4747  const unsigned int second_dealii_idx_on_face
4748  = lower_idx==0?left_to_right[it->orientation.to_ulong()][first_dealii_idx_on_face]:
4749  right_to_left[it->orientation.to_ulong()][first_dealii_idx_on_face];
4750  const unsigned int second_dealii_idx_on_cell
4752  (face_idx_list[higher_idx], second_dealii_idx_on_face,
4753  cell_list[higher_idx]->face_orientation(face_idx_list[higher_idx]),
4754  cell_list[higher_idx]->face_flip(face_idx_list[higher_idx]),
4755  cell_list[higher_idx]->face_rotation(face_idx_list[higher_idx]));
4756  //map back to p4est
4757  const unsigned int second_p4est_idx_on_face
4758  = p8est_corner_face_corners[second_dealii_idx_on_cell][face_idx_list[higher_idx]];
4759  p4est_orientation = second_p4est_idx_on_face;
4760  }
4761 
4763  connectivity_join_faces (connectivity,
4764  tree_left,
4765  tree_right,
4766  face_left,
4767  face_right,
4768  p4est_orientation);
4769  }
4770 
4771 
4773  (connectivity) == 1, ExcInternalError());
4774 
4775  // now create a forest out of the connectivity data structure
4777  parallel_forest
4780  connectivity,
4781  /* minimum initial number of quadrants per tree */ 0,
4782  /* minimum level of upfront refinement */ 0,
4783  /* use uniform upfront refinement */ 1,
4784  /* user_data_size = */ 0,
4785  /* user_data_constructor = */ NULL,
4786  /* user_pointer */ this);
4787 
4788 
4789  try
4790  {
4792  }
4793  catch (const typename Triangulation<dim>::DistortedCellList &)
4794  {
4795  // the underlying triangulation should not be checking for distorted
4796  // cells
4797  AssertThrow (false, ExcInternalError());
4798  }
4799 
4801  periodicity_vector.begin(),
4802  periodicity_vector.end());
4803 
4804 #else
4805  Assert(false, ExcMessage ("Need p4est version >= 0.3.4.1!"));
4806 #endif
4807  }
4808 
4809 
4810 
4811  template <int dim, int spacedim>
4812  std::size_t
4814  {
4815  std::size_t mem=
4819  + MemoryConsumption::memory_consumption(parallel_forest)
4823 // + MemoryConsumption::memory_consumption(attached_data_pack_callbacks) //TODO[TH]: how?
4824  + MemoryConsumption::memory_consumption(coarse_cell_to_p4est_tree_permutation)
4825  + MemoryConsumption::memory_consumption(p4est_tree_to_coarse_cell_permutation)
4827 
4828  return mem;
4829  }
4830 
4831 
4832 
4833  template <int dim, int spacedim>
4834  std::size_t
4836  {
4837  return ::internal::p4est::functions<dim>::forest_memory_used(parallel_forest)
4839  }
4840 
4841 
4842 
4843  template <int dim, int spacedim>
4844  void
4846  copy_triangulation (const ::Triangulation<dim, spacedim> &old_tria)
4847  {
4848  try
4849  {
4851  copy_triangulation (old_tria);
4852  }
4853  catch (const typename ::Triangulation<dim,spacedim>::DistortedCellList &)
4854  {
4855  // the underlying triangulation should not be checking for distorted
4856  // cells
4857  AssertThrow (false, ExcInternalError());
4858  }
4859 
4860  // note that now we have some content in the p4est objects and call the
4861  // functions that do the actual work (which are dimension dependent, so
4862  // separate)
4864 
4865  Assert (old_tria.n_levels() == 1,
4866  ExcMessage ("Parallel distributed triangulations can only be copied, "
4867  "if they are not refined!"));
4868 
4869  if (const ::parallel::distributed::Triangulation<dim,spacedim> *
4870  old_tria_x = dynamic_cast<const ::parallel::distributed::Triangulation<dim,spacedim> *>(&old_tria))
4871  {
4872  Assert (!old_tria_x->refinement_in_progress,
4873  ExcMessage ("Parallel distributed triangulations can only "
4874  "be copied, if no refinement is in progress!"));
4875 
4876  // duplicate MPI communicator, stored in the base class
4878 
4879  coarse_cell_to_p4est_tree_permutation = old_tria_x->coarse_cell_to_p4est_tree_permutation;
4880  p4est_tree_to_coarse_cell_permutation = old_tria_x->p4est_tree_to_coarse_cell_permutation;
4881  attached_data_size = old_tria_x->attached_data_size;
4882  n_attached_datas = old_tria_x->n_attached_datas;
4883 
4884  settings = old_tria_x->settings;
4885  }
4886  else
4887  {
4889  }
4890 
4892 
4893  try
4894  {
4896  }
4897  catch (const typename Triangulation<dim>::DistortedCellList &)
4898  {
4899  // the underlying triangulation should not be checking for distorted
4900  // cells
4901  AssertThrow (false, ExcInternalError());
4902  }
4903 
4904  this->update_number_cache ();
4905  }
4906 
4907 
4908  template <int dim, int spacedim>
4909  void
4912  {
4913  // determine size of memory in bytes to attach to each cell. This needs
4914  // to be constant because of p4est.
4915  if (attached_data_size==0)
4916  {
4917  Assert(n_attached_datas==0, ExcInternalError());
4918 
4919  //nothing to do
4920  return;
4921  }
4922 
4923  // realloc user_data in p4est
4924  void *userptr = parallel_forest->user_pointer;
4927  NULL, NULL);
4928  parallel_forest->user_pointer = userptr;
4929 
4930 
4931  // Recurse over p4est and Triangulation
4932  // to find refined/coarsened/kept
4933  // cells. Then query and attach the data.
4935  cell = this->begin(0);
4936  cell != this->end(0);
4937  ++cell)
4938  {
4939  //skip coarse cells, that are not ours
4940  if (tree_exists_locally<dim,spacedim>(parallel_forest,
4941  coarse_cell_to_p4est_tree_permutation[cell->index()])
4942  == false)
4943  continue;
4944 
4945  typename ::internal::p4est::types<dim>::quadrant p4est_coarse_cell;
4946  typename ::internal::p4est::types<dim>::tree *tree =
4947  init_tree(cell->index());
4948 
4949  ::internal::p4est::init_coarse_quadrant<dim>(p4est_coarse_cell);
4950 
4951  attach_mesh_data_recursively<dim,spacedim>(*tree,
4952  cell,
4953  p4est_coarse_cell,
4955  }
4956  }
4957 
4958  template <int dim, int spacedim>
4959  std::vector<unsigned int>
4962  {
4963  // Allocate the space for the weights. In fact we do not know yet, how
4964  // many cells we own after the refinement (only p4est knows that
4965  // at this point). We simply reserve n_active_cells space and if many
4966  // more cells are refined than coarsened than additional reallocation
4967  // will be done inside get_cell_weights_recursively.
4968  std::vector<unsigned int> weights;
4969  weights.reserve(this->n_active_cells());
4970 
4971  // Recurse over p4est and Triangulation
4972  // to find refined/coarsened/kept
4973  // cells. Then append cell_weight.
4974  // Note that we need to follow the p4est ordering
4975  // instead of the deal.II ordering to get the cell_weights
4976  // in the same order p4est will encounter them during repartitioning.
4977  for (unsigned int c=0; c<this->n_cells(0); ++c)
4978  {
4979  // skip coarse cells, that are not ours
4980  if (tree_exists_locally<dim,spacedim>(parallel_forest,c) == false)
4981  continue;
4982 
4983  const unsigned int coarse_cell_index =
4984  p4est_tree_to_coarse_cell_permutation[c];
4985 
4987  dealii_coarse_cell (this, 0, coarse_cell_index);
4988 
4989  typename ::internal::p4est::types<dim>::quadrant p4est_coarse_cell;
4991  quadrant_set_morton (&p4est_coarse_cell,
4992  /*level=*/0,
4993  /*index=*/0);
4994  p4est_coarse_cell.p.which_tree = c;
4995 
4996  const typename ::internal::p4est::types<dim>::tree *tree =
4997  init_tree(coarse_cell_index);
4998 
4999  get_cell_weights_recursively<dim,spacedim>(*tree,
5000  dealii_coarse_cell,
5001  p4est_coarse_cell,
5002  this->signals,
5003  weights);
5004  }
5005 
5006  return weights;
5007  }
5008 
5009  template <int dim, int spacedim>
5010  typename ::Triangulation<dim,spacedim>::cell_iterator
5011  cell_from_quad
5012  (typename ::parallel::distributed::Triangulation<dim,spacedim> *triangulation,
5013  typename ::internal::p4est::types<dim>::topidx treeidx,
5014  typename ::internal::p4est::types<dim>::quadrant &quad)
5015  {
5016  int i, l = quad.level;
5017  int child_id;
5018  types::global_dof_index dealii_index =
5019  triangulation->get_p4est_tree_to_coarse_cell_permutation()[treeidx];
5020 
5021  for (i = 0; i < l; i++)
5022  {
5023  typename ::Triangulation<dim,spacedim>::cell_iterator cell (triangulation, i, dealii_index);
5025  Assert (cell->has_children (), ExcMessage ("p4est quadrant does not correspond to a cell!"));
5026  dealii_index = cell->child_index(child_id);
5027  }
5028 
5029  typename ::Triangulation<dim,spacedim>::cell_iterator out_cell (triangulation, l, dealii_index);
5030 
5031  return out_cell;
5032  }
5033 
5034 
5035 
5036  template <int spacedim>
5038  :
5039  ::parallel::Triangulation<1,spacedim>(MPI_COMM_WORLD,
5040  typename ::Triangulation<1,spacedim>::MeshSmoothing(),
5041  false)
5042  {
5043  Assert (false, ExcNotImplemented());
5044  }
5045 
5046 
5047  template <int spacedim>
5049  {
5050  Assert (false, ExcNotImplemented());
5051  }
5052 
5053 
5054 
5055  template <int spacedim>
5056  void
5058  (const std::vector<bool> &/*vertex_locally_moved*/)
5059  {
5060  Assert (false, ExcNotImplemented());
5061  }
5062 
5063 
5064  template <int spacedim>
5065  const std::vector<types::global_dof_index> &
5067  {
5068  static std::vector<types::global_dof_index> a;
5069  return a;
5070  }
5071 
5072  template <int spacedim>
5073  void
5076  (std::map<unsigned int, std::set<::types::subdomain_id> >
5077  &/*vertices_with_ghost_neighbors*/)
5078  {
5079  Assert (false, ExcNotImplemented());
5080  }
5081 
5082  template <int spacedim>
5083  void
5086  (const unsigned int /*level*/,
5087  std::map<unsigned int, std::set<::types::subdomain_id> >
5088  &/*vertices_with_ghost_neighbors*/)
5089  {
5090  Assert (false, ExcNotImplemented());
5091  }
5092 
5093  template <int spacedim>
5094  std::vector<bool>
5096  mark_locally_active_vertices_on_level (const unsigned int) const
5097  {
5098  Assert (false, ExcNotImplemented());
5099  return std::vector<bool>();
5100  }
5101  }
5102 }
5103 
5104 
5105 #else // DEAL_II_WITH_P4EST
5106 
5107 namespace parallel
5108 {
5109  namespace distributed
5110  {
5111  template <int dim, int spacedim>
5113  :
5115  {
5116  Assert (false, ExcNotImplemented());
5117  }
5118  }
5119 }
5120 
5121 #endif // DEAL_II_WITH_P4EST
5122 
5123 
5124 
5125 
5126 /*-------------- Explicit Instantiations -------------------------------*/
5127 #include "tria.inst"
5128 
5129 
5130 DEAL_II_NAMESPACE_CLOSE
unsigned int n_active_cells() const
Definition: tria.cc:10973
virtual void copy_triangulation(const ::Triangulation< dim, spacedim > &old_tria)
Definition: tria.cc:4846
std::vector< unsigned int > get_cell_weights()
Definition: tria.cc:4961
static unsigned int standard_to_real_face_vertex(const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
unsigned int n_vertices() const
static const unsigned int invalid_unsigned_int
Definition: types.h:164
std::vector< bool > mark_locally_active_vertices_on_level(const unsigned int level) const
Definition: tria.cc:4644
void fill_vertices_with_ghost_neighbors(std::map< unsigned int, std::set<::types::subdomain_id > > &vertices_with_ghost_neighbors)
Definition: tria.cc:4574
types::subdomain_id my_subdomain
Definition: tria_base.h:168
unsigned int register_data_attach(const std::size_t size, const std_cxx11::function< void(const cell_iterator &, const CellStatus, void *)> &pack_callback)
Definition: tria.cc:4065
unsigned int n_cells() const
Definition: tria.cc:10966
virtual void copy_triangulation(const ::Triangulation< dim, spacedim > &old_tria)
Definition: tria_base.cc:60
static unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
std::vector< Point< spacedim > > vertices
Definition: tria.h:3296
::ExceptionBase & ExcMessage(std::string arg1)
virtual std::size_t memory_consumption() const
Definition: tria.cc:4813
boost::signals2::signal< unsigned int(const cell_iterator &, const CellStatus), CellWeightSum< unsigned int > > cell_weight
Definition: tria.h:2072
virtual void copy_triangulation(const Triangulation< dim, spacedim > &old_tria)
Definition: tria.cc:9043
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: tria.cc:10397
void notify_ready_to_unpack(const unsigned int offset, const std_cxx11::function< void(const cell_iterator &, const CellStatus, const void *)> &unpack_callback)
Definition: tria.cc:4088
void copy_new_triangulation_to_p4est(::internal::int2type< 2 >)
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
void fill_level_vertices_with_ghost_neighbors(const unsigned int level, std::map< unsigned int, std::set<::types::subdomain_id > > &vertices_with_ghost_neighbors)
Definition: tria.cc:4620
virtual std::size_t memory_consumption_p4est() const
Definition: tria.cc:4835
void load(const char *filename, const bool autopartition=true)
Definition: tria.cc:2642
callback_list_t attached_data_pack_callbacks
Definition: tria.h:863
cell_iterator begin(const unsigned int level=0) const
Definition: tria.cc:10377
std::bitset< 3 > orientation
Definition: grid_tools.h:1175
const std::vector< types::global_dof_index > & get_p4est_tree_to_coarse_cell_permutation() const
Definition: tria.cc:4156
void add_periodicity(const std::vector< GridTools::PeriodicFacePair< cell_iterator > > &)
Definition: tria.cc:4666
typename ::internal::p4est::types< dim >::ghost * parallel_ghost
Definition: tria.h:818
unsigned int n_levels() const
typename ::internal::p4est::types< dim >::forest * parallel_forest
Definition: tria.h:813
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
Definition: tria.cc:2152
virtual bool prepare_coarsening_and_refinement()
Definition: tria.cc:3291
void get_vertex_connectivity_of_cells(const Triangulation< dim, spacedim > &triangulation, DynamicSparsityPattern &connectivity)
Definition: grid_tools.cc:1933
cell_iterator end() const
Definition: tria.cc:10465
virtual void execute_coarsening_and_refinement()
Definition: tria.cc:11602
unsigned int n_active_lines() const
Definition: tria.cc:11188
virtual bool prepare_coarsening_and_refinement()
Definition: tria.cc:12269
Triangulation(const MeshSmoothing smooth_grid=none, const bool check_for_distorted_cells=false)
Definition: tria.cc:8751
T sum(const T &t, const MPI_Comm &mpi_communicator)
Definition: mpi.h:611
unsigned int global_dof_index
Definition: types.h:88
static unsigned int child_cell_on_face(const RefinementCase< dim > &ref_case, const unsigned int face, const unsigned int subface, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false, const RefinementCase< dim-1 > &face_refinement_case=RefinementCase< dim-1 >::isotropic_refinement)
#define Assert(cond, exc)
Definition: exceptions.h:294
virtual void update_number_cache()
Definition: tria_base.cc:143
Signals signals
Definition: tria.h:2078
const std::vector< types::global_dof_index > & get_coarse_cell_to_p4est_tree_permutation() const
Definition: tria.cc:4165
virtual MPI_Comm get_communicator() const
Definition: tria_base.cc:135
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
Definition: tria.cc:9130
void save_coarsen_flags(std::ostream &out) const
Definition: tria.cc:9435
const std::vector< Point< spacedim > > & get_vertices() const
virtual bool has_hanging_nodes() const
Definition: tria.cc:2536
virtual std::size_t memory_consumption() const
Definition: tria_base.cc:76
void save_refine_flags(std::ostream &out) const
Definition: tria.cc:9370
unsigned int face_idx[2]
Definition: grid_tools.h:1168
void save(const char *filename) const
Definition: tria.cc:2594
typename ::internal::p4est::types< dim >::tree * init_tree(const int dealii_coarse_cell_index) const
Definition: tria.cc:2761
virtual void execute_coarsening_and_refinement()
Definition: tria.cc:3679
std::vector< unsigned int > invert_permutation(const std::vector< unsigned int > &permutation)
Definition: utilities.cc:483
unsigned int subdomain_id
Definition: types.h:42
virtual unsigned int n_global_levels() const
Definition: tria_base.cc:114
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:99
::Triangulation< dim, spacedim >::active_cell_iterator active_cell_iterator
Definition: tria.h:368
void reorder_hierarchical(const DynamicSparsityPattern &sparsity, std::vector< DynamicSparsityPattern::size_type > &new_indices)
const types::subdomain_id artificial_subdomain_id
Definition: types.h:255
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
void communicate_locally_moved_vertices(const std::vector< bool > &vertex_locally_moved)
Definition: tria.cc:3915
const std::vector< bool > & get_used_vertices() const
Definition: tria.cc:11511
unsigned int get_checksum() const
Definition: tria.cc:2725
MeshSmoothing smooth_grid
Definition: tria.h:3053
void setup_coarse_cell_to_p4est_tree_permutation()
Definition: tria.cc:2565
::Triangulation< dim, spacedim >::cell_iterator cell_iterator
Definition: tria.h:348
std::vector< bool > get_locally_owned_vertices(const Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:2075
std::vector< GridTools::PeriodicFacePair< cell_iterator > > periodic_face_pairs_level_0
Definition: tria.h:887
active_cell_iterator last_active() const
Definition: tria.cc:10444
Triangulation(MPI_Comm mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing smooth_grid=(::Triangulation< dim, spacedim >::none), const Settings settings=default_setting)
Definition: tria.cc:2111
unsigned int n_attached_deserialize
Definition: tria.h:849
types::subdomain_id locally_owned_subdomain() const
Definition: tria_base.cc:218
virtual void update_number_cache()
Definition: tria.cc:2735
std::vector< unsigned int > compute_point_to_point_communication_pattern(const MPI_Comm &mpi_comm, const std::vector< unsigned int > &destinations)
Definition: mpi.cc:126
T max(const T &t, const MPI_Comm &mpi_communicator)
Definition: mpi.h:693
virtual void clear()
Definition: tria.cc:8822
std::vector< types::global_dof_index > coarse_cell_to_p4est_tree_permutation
Definition: tria.h:879
void write_mesh_vtk(const char *file_basename) const
Definition: tria.cc:2582