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 #include <deal.II/distributed/p4est_wrappers.h> 35 DEAL_II_NAMESPACE_OPEN
38 #ifdef DEAL_II_WITH_P4EST 42 template <
int dim,
int spacedim>
45 std::vector<unsigned int> &vertex_touch_count,
46 std::vector<std::list<
50 vertex_touch_count.resize (triangulation.n_vertices());
51 vertex_to_cell.resize (triangulation.n_vertices());
54 cell = triangulation.begin_active();
55 cell != triangulation.end(); ++cell)
56 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
58 ++vertex_touch_count[cell->vertex_index(v)];
59 vertex_to_cell[cell->vertex_index(v)].emplace_back (cell, v);
65 template <
int dim,
int spacedim>
68 std::vector<unsigned int> &edge_touch_count,
69 std::vector<std::list<
75 edge_touch_count.resize (triangulation.n_active_lines());
76 edge_to_cell.resize (triangulation.n_active_lines());
79 cell = triangulation.begin_active();
80 cell != triangulation.end(); ++cell)
81 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++
l)
83 ++edge_touch_count[cell->line(l)->index()];
84 edge_to_cell[cell->line(l)->index()].emplace_back (cell, l);
94 template <
int dim,
int spacedim>
97 const std::vector<unsigned int> &vertex_touch_count,
98 const std::vector<std::list<
101 const std::vector<types::global_dof_index> &coarse_cell_to_p4est_tree_permutation,
102 const bool set_vertex_info,
111 Assert (triangulation.get_used_vertices().size() ==
112 triangulation.get_vertices().size(),
114 Assert (std::find (triangulation.get_used_vertices().begin(),
115 triangulation.get_used_vertices().end(),
117 == triangulation.get_used_vertices().end(),
119 if (set_vertex_info ==
true)
120 for (
unsigned int v=0; v<triangulation.n_vertices(); ++v)
122 connectivity->vertices[3*v ] = triangulation.get_vertices()[v][0];
123 connectivity->vertices[3*v+1] = triangulation.get_vertices()[v][1];
124 connectivity->vertices[3*v+2] = (spacedim == 2 ?
127 triangulation.get_vertices()[v][2]);
136 cell = triangulation.begin_active(),
137 endc = triangulation.end();
138 for (; cell != endc; ++cell)
141 index = coarse_cell_to_p4est_tree_permutation[cell->index()];
143 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
145 if (set_vertex_info ==
true)
152 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
153 if (cell->face(f)->at_boundary() ==
false)
155 = coarse_cell_to_p4est_tree_permutation[cell->neighbor(f)->index()];
158 = coarse_cell_to_p4est_tree_permutation[cell->index()];
162 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
163 if (cell->face(f)->at_boundary() ==
false)
170 = cell->neighbor_of_neighbor (f);
189 connectivity->tree_to_face[index*6 + f]
190 = cell->neighbor_of_neighbor (f);
192 unsigned int face_idx_list[2] =
193 {f, cell->neighbor_of_neighbor (f)};
195 cell_list[2] = {cell, cell->neighbor(f)};
196 unsigned int smaller_idx = 0;
198 if (f>cell->neighbor_of_neighbor (f))
201 unsigned int larger_idx = (smaller_idx+1) % 2;
210 cell_list[smaller_idx]->vertex_index(
212 face_idx_list[smaller_idx],
214 cell_list[smaller_idx]->face_orientation(face_idx_list[smaller_idx]),
215 cell_list[smaller_idx]->face_flip(face_idx_list[smaller_idx]),
216 cell_list[smaller_idx]->face_rotation(face_idx_list[smaller_idx]))
221 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_face; ++i)
225 cell_list[larger_idx]->vertex_index(
227 face_idx_list[larger_idx],
238 connectivity->tree_to_face[index*6 + f] += 6*v;
251 connectivity->ctt_offset[0] = 0;
252 std::partial_sum (vertex_touch_count.begin(),
253 vertex_touch_count.end(),
254 &connectivity->ctt_offset[1]);
257 num_vtt = std::accumulate (vertex_touch_count.begin(),
258 vertex_touch_count.end(),
261 Assert (connectivity->ctt_offset[triangulation.n_vertices()] ==
265 for (
unsigned int v=0; v<triangulation.n_vertices(); ++v)
267 Assert (vertex_to_cell[v].size() == vertex_touch_count[v],
270 typename std::list<std::pair
272 unsigned int> >::const_iterator
273 p = vertex_to_cell[v].begin();
274 for (
unsigned int c=0; c<vertex_touch_count[v]; ++c, ++p)
276 connectivity->corner_to_tree[connectivity->ctt_offset[v]+c]
277 = coarse_cell_to_p4est_tree_permutation[p->first->index()];
278 connectivity->corner_to_corner[connectivity->ctt_offset[v]+c]
286 template <
int dim,
int spacedim>
291 Assert (coarse_grid_cell < parallel_forest->connectivity->num_trees,
293 return ((coarse_grid_cell >= parallel_forest->first_local_tree)
295 (coarse_grid_cell <= parallel_forest->last_local_tree));
299 template <
int dim,
int spacedim>
303 if (cell->has_children())
304 for (
unsigned int c=0; c<cell->n_children(); ++c)
305 delete_all_children_and_self<dim,spacedim> (cell->child(c));
307 cell->set_coarsen_flag ();
312 template <
int dim,
int spacedim>
316 if (cell->has_children())
317 for (
unsigned int c=0; c<cell->n_children(); ++c)
318 delete_all_children_and_self<dim,spacedim> (cell->child(c));
322 template <
int dim,
int spacedim>
330 const std::vector<std::vector<bool> > &marked_vertices)
337 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
339 if (marked_vertices[dealii_cell->level()][dealii_cell->vertex_index(v)])
358 if (!used && dealii_cell->active() && dealii_cell->is_artificial()==
false 359 && dealii_cell->level()+1<
static_cast<int>(marked_vertices.size()))
361 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
363 if (marked_vertices[dealii_cell->level()+1][dealii_cell->vertex_index(v)])
372 if (!used && dealii_cell->active() && dealii_cell->is_artificial()==
false 373 && dealii_cell->level()>0)
375 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
377 if (marked_vertices[dealii_cell->level()-1][dealii_cell->vertex_index(v)])
391 Assert((owner!=-2) && (owner!=-1),
ExcMessage(
"p4est should know the owner."));
392 dealii_cell->set_level_subdomain_id(owner);
397 if (dealii_cell->has_children ())
401 for (
unsigned int c=0;
402 c<GeometryInfo<dim>::max_children_per_cell; ++c)
406 P4EST_QUADRANT_INIT(&p4est_child[c]);
409 P8EST_QUADRANT_INIT(&p4est_child[c]);
420 for (
unsigned int c=0;
421 c<GeometryInfo<dim>::max_children_per_cell; ++c)
423 determine_level_subdomain_id_recursively <dim,spacedim> (tree,tree_index,
424 dealii_cell->child(c),
434 template <
int dim,
int spacedim>
443 if (sc_array_bsearch(const_cast<sc_array_t *>(&tree.quadrants),
449 delete_all_children<dim,spacedim> (dealii_cell);
450 if (!dealii_cell->has_children())
451 dealii_cell->set_subdomain_id(my_subdomain);
460 if (dealii_cell->has_children () ==
false)
461 dealii_cell->set_refine_flag ();
466 for (
unsigned int c=0;
467 c<GeometryInfo<dim>::max_children_per_cell; ++c)
471 P4EST_QUADRANT_INIT(&p4est_child[c]);
474 P8EST_QUADRANT_INIT(&p4est_child[c]);
485 for (
unsigned int c=0;
486 c<GeometryInfo<dim>::max_children_per_cell; ++c)
496 delete_all_children<dim,spacedim> (dealii_cell->child(c));
497 dealii_cell->child(c)
504 match_tree_recursively<dim,spacedim> (tree,
505 dealii_cell->child(c),
515 template <
int dim,
int spacedim>
517 match_quadrant (const ::Triangulation<dim,spacedim> *tria,
518 unsigned int dealii_index,
522 const int l = ghost_quadrant.level;
524 for (
int i = 0; i <
l; ++i)
527 if (cell->has_children () ==
false)
529 cell->clear_coarsen_flag();
530 cell->set_refine_flag ();
535 dealii_index = cell->child_index(child_id);
539 if (cell->has_children())
540 delete_all_children<dim,spacedim> (cell);
543 cell->clear_coarsen_flag();
544 cell->set_subdomain_id(ghost_owner);
550 template <
int dim,
int spacedim>
555 const std::map<
unsigned int,
typename std::function<
557 typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus,
561 int idx = sc_array_bsearch(const_cast<sc_array_t *>(&tree.quadrants),
571 bool p4est_has_children = (idx == -1);
573 if (p4est_has_children && dealii_cell->has_children())
578 for (
unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
582 P4EST_QUADRANT_INIT(&p4est_child[c]);
585 P8EST_QUADRANT_INIT(&p4est_child[c]);
594 for (
unsigned int c=0;
595 c<GeometryInfo<dim>::max_children_per_cell; ++c)
597 attach_mesh_data_recursively<dim,spacedim> (tree,
598 dealii_cell->child(c),
603 else if (!p4est_has_children && !dealii_cell->has_children())
608 sc_array_index (const_cast<sc_array_t *>(&tree.quadrants), idx)
617 for (
const auto &it : pack_callbacks)
619 void *ptr =
static_cast<char *
>(q->p.user_data) + it.first;
620 (it.second)(dealii_cell,
625 else if (p4est_has_children)
634 for (
unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
638 P4EST_QUADRANT_INIT(&p4est_child[c]);
641 P8EST_QUADRANT_INIT(&p4est_child[c]);
649 int child0_idx = sc_array_bsearch(const_cast<sc_array_t *>(&tree.quadrants),
652 Assert(child0_idx != -1,
ExcMessage(
"the first child should exist as an active quadrant!"));
656 sc_array_index (const_cast<sc_array_t *>(&tree.quadrants), child0_idx)
665 for (
const auto &it : pack_callbacks)
669 void *ptr =
static_cast<char *
>(q->p.user_data) + it.first;
671 (it.second)(dealii_cell,
677 for (
unsigned int i=1; i<GeometryInfo<dim>::max_children_per_cell; ++i)
679 int child_idx = sc_array_bsearch(const_cast<sc_array_t *>(&tree.quadrants),
683 sc_array_index (const_cast<sc_array_t *>(&tree.quadrants), child_idx)
696 sc_array_index (const_cast<sc_array_t *>(&tree.quadrants), idx)
705 for (
unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
709 for (
const auto &it : pack_callbacks)
711 void *ptr =
static_cast<char *
>(q->p.user_data) + it.first;
712 (it.second)(dealii_cell,
719 template <
int dim,
int spacedim>
725 std::vector<unsigned int> &weight)
727 const int idx = sc_array_bsearch(const_cast<sc_array_t *>(&tree.quadrants),
737 const bool p4est_has_children = (idx == -1);
739 if (p4est_has_children && dealii_cell->has_children())
744 for (
unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
748 P4EST_QUADRANT_INIT(&p4est_child[c]);
751 P8EST_QUADRANT_INIT(&p4est_child[c]);
760 for (
unsigned int c=0;
761 c<GeometryInfo<dim>::max_children_per_cell; ++c)
763 get_cell_weights_recursively<dim,spacedim> (tree,
764 dealii_cell->child(c),
770 else if (!p4est_has_children && !dealii_cell->has_children())
773 weight.push_back(1000);
777 else if (p4est_has_children)
780 unsigned int parent_weight(1000);
784 for (
unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
787 weight.push_back(parent_weight);
793 weight.push_back(1000);
800 template <
int dim,
int spacedim>
806 const unsigned int handle,
807 const typename std::function<
811 int idx = sc_array_bsearch(const_cast<sc_array_t *>(&tree.quadrants),
822 const bool p4est_has_children = (idx == -1);
823 if (p4est_has_children)
830 for (
unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
834 P4EST_QUADRANT_INIT(&p4est_child[c]);
837 P8EST_QUADRANT_INIT(&p4est_child[c]);
846 for (
unsigned int c=0;
847 c<GeometryInfo<dim>::max_children_per_cell; ++c)
849 post_mesh_data_recursively<dim,spacedim> (tree,
850 dealii_cell->child(c),
863 sc_array_index (const_cast<sc_array_t *>(&tree.quadrants), idx)
866 void *ptr =
static_cast<char *
>(q->p.user_data) + handle;
867 typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus
868 status = *
static_cast< 869 typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus *
875 unpack_callback(dealii_cell, status, ptr);
880 unpack_callback(parent_cell, status, ptr);
885 unpack_callback(dealii_cell, status, ptr);
905 template <
int dim,
int spacedim>
906 class RefineAndCoarsenList
910 const std::vector<types::global_dof_index> &p4est_tree_to_coarse_cell_permutation,
937 bool pointers_are_at_end ()
const;
940 std::vector<typename internal::p4est::types<dim>::quadrant> refine_list;
941 typename std::vector<typename internal::p4est::types<dim>::quadrant>::const_iterator current_refine_pointer;
943 std::vector<typename internal::p4est::types<dim>::quadrant> coarsen_list;
944 typename std::vector<typename internal::p4est::types<dim>::quadrant>::const_iterator current_coarsen_pointer;
953 template <
int dim,
int spacedim>
955 RefineAndCoarsenList<dim,spacedim>::
956 pointers_are_at_end ()
const 958 return ((current_refine_pointer == refine_list.end())
960 (current_coarsen_pointer == coarsen_list.end()));
965 template <
int dim,
int spacedim>
966 RefineAndCoarsenList<dim,spacedim>::
968 const std::vector<types::global_dof_index> &p4est_tree_to_coarse_cell_permutation,
972 unsigned int n_refine_flags = 0,
975 cell = triangulation.begin_active();
976 cell != triangulation.end(); ++cell)
979 if (cell->subdomain_id() != my_subdomain)
982 if (cell->refine_flag_set())
984 else if (cell->coarsen_flag_set())
988 refine_list.reserve (n_refine_flags);
989 coarsen_list.reserve (n_coarsen_flags);
999 for (
unsigned int c=0; c<triangulation.n_cells(0); ++c)
1001 unsigned int coarse_cell_index =
1002 p4est_tree_to_coarse_cell_permutation[c];
1005 cell (&triangulation, 0, coarse_cell_index);
1012 p4est_cell.p.which_tree = c;
1013 build_lists (cell, p4est_cell, my_subdomain);
1017 Assert(refine_list.size() == n_refine_flags,
1019 Assert(coarsen_list.size() == n_coarsen_flags,
1023 for (
unsigned int i=1; i<refine_list.size(); ++i)
1024 Assert (refine_list[i].p.which_tree >=
1025 refine_list[i-1].p.which_tree,
1027 for (
unsigned int i=1; i<coarsen_list.size(); ++i)
1028 Assert (coarsen_list[i].p.which_tree >=
1029 coarsen_list[i-1].p.which_tree,
1032 current_refine_pointer = refine_list.begin();
1033 current_coarsen_pointer = coarsen_list.begin();
1038 template <
int dim,
int spacedim>
1040 RefineAndCoarsenList<dim,spacedim>::
1045 if (!cell->has_children())
1047 if (cell->subdomain_id() == my_subdomain)
1049 if (cell->refine_flag_set())
1050 refine_list.push_back (p4est_cell);
1051 else if (cell->coarsen_flag_set())
1052 coarsen_list.push_back (p4est_cell);
1059 for (
unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
1063 P4EST_QUADRANT_INIT(&p4est_child[c]);
1066 P8EST_QUADRANT_INIT(&p4est_child[c]);
1074 for (
unsigned int c=0;
1075 c<GeometryInfo<dim>::max_children_per_cell; ++c)
1077 p4est_child[c].p.which_tree = p4est_cell.p.which_tree;
1078 build_lists (cell->child(c),
1086 template <
int dim,
int spacedim>
1088 RefineAndCoarsenList<dim,spacedim>::
1093 RefineAndCoarsenList<dim,spacedim> *this_object
1094 =
reinterpret_cast<RefineAndCoarsenList<dim,spacedim>*
>(forest->user_pointer);
1098 if (this_object->current_refine_pointer == this_object->refine_list.end())
1101 Assert (coarse_cell_index <=
1102 this_object->current_refine_pointer->p.which_tree,
1107 if (coarse_cell_index <
1108 this_object->current_refine_pointer->p.which_tree)
1112 Assert (coarse_cell_index <=
1113 this_object->current_refine_pointer->p.which_tree,
1119 quadrant_compare (quadrant,
1120 &*this_object->current_refine_pointer) <= 0,
1125 quadrant_is_equal (quadrant, &*this_object->current_refine_pointer))
1127 ++this_object->current_refine_pointer;
1137 template <
int dim,
int spacedim>
1139 RefineAndCoarsenList<dim,spacedim>::
1144 RefineAndCoarsenList<dim,spacedim> *this_object
1145 =
reinterpret_cast<RefineAndCoarsenList<dim,spacedim>*
>(forest->user_pointer);
1149 if (this_object->current_coarsen_pointer ==
1150 this_object->coarsen_list.end())
1153 Assert (coarse_cell_index <=
1154 this_object->current_coarsen_pointer->p.which_tree,
1159 if (coarse_cell_index <
1160 this_object->current_coarsen_pointer->p.which_tree)
1164 Assert (coarse_cell_index <=
1165 this_object->current_coarsen_pointer->p.which_tree,
1171 quadrant_compare (children[0],
1172 &*this_object->current_coarsen_pointer) <= 0,
1178 quadrant_is_equal (children[0],
1179 &*this_object->current_coarsen_pointer))
1182 ++this_object->current_coarsen_pointer;
1186 for (
unsigned int c=1; c<GeometryInfo<dim>::max_children_per_cell; ++c)
1189 quadrant_is_equal (children[c],
1190 &*this_object->current_coarsen_pointer),
1192 ++this_object->current_coarsen_pointer;
1210 template <
int dim,
int spacedim>
1211 class PartitionWeights
1219 explicit PartitionWeights (
const std::vector<unsigned int> &cell_weights);
1235 std::vector<unsigned int> cell_weights_list;
1236 std::vector<unsigned int>::const_iterator current_pointer;
1240 template <
int dim,
int spacedim>
1241 PartitionWeights<dim,spacedim>::
1242 PartitionWeights (
const std::vector<unsigned int> &cell_weights)
1244 cell_weights_list(cell_weights)
1248 current_pointer = cell_weights_list.begin();
1252 template <
int dim,
int spacedim>
1254 PartitionWeights<dim,spacedim>::
1263 PartitionWeights<dim,spacedim> *this_object
1264 =
reinterpret_cast<PartitionWeights<dim,spacedim>*
>(forest->user_pointer);
1266 Assert (this_object->current_pointer >= this_object->cell_weights_list.begin(),
1268 Assert (this_object->current_pointer < this_object->cell_weights_list.end(),
1272 return *this_object->current_pointer++;
1279 namespace distributed
1285 template <
int dim,
int spacedim>
1288 const typename ::Triangulation<dim,spacedim>::MeshSmoothing smooth_grid,
1296 (settings_ & construct_multigrid_hierarchy) ?
1298 (smooth_grid |
Triangulation<dim,spacedim>::limit_level_difference_at_vertices) :
1301 settings(settings_),
1302 triangulation_has_content (false),
1303 connectivity (nullptr),
1304 parallel_forest (nullptr),
1305 cell_attached_data ({0, 0, 0, {}})
1307 parallel_ghost =
nullptr;
1312 template <
int dim,
int spacedim>
1334 template <
int dim,
int spacedim>
1346 catch (
const typename ::Triangulation<dim,spacedim>::DistortedCellList &)
1356 triangulation_has_content =
true;
1358 setup_coarse_cell_to_p4est_tree_permutation ();
1360 copy_new_triangulation_to_p4est (std::integral_constant<int, dim>());
1364 copy_local_forest_to_triangulation ();
1373 this->update_number_cache ();
1374 this->update_periodic_face_map();
1380 namespace CommunicateLocallyMovedVertices
1389 template <
int dim,
int spacedim>
1393 std::vector<unsigned int> tree_index;
1395 std::vector<typename ::internal::p4est::types<dim>::quadrant> quadrants;
1398 std::vector<unsigned int> vertex_indices;
1401 std::vector<::Point<spacedim> > vertices;
1405 std::vector<unsigned int * > first_vertex_indices;
1406 std::vector<::Point<spacedim>* > first_vertices;
1408 unsigned int bytes_for_buffer ()
const 1410 return sizeof(
unsigned int) +
1411 tree_index.size() *
sizeof(
unsigned int) +
1412 quadrants.size() *
sizeof(typename ::internal::p4est
1413 ::types<dim>::quadrant) +
1414 vertex_indices.size() *
sizeof(
unsigned int) +
1418 void pack_data (std::vector<char> &buffer)
const 1420 buffer.resize(bytes_for_buffer());
1422 char *ptr = buffer.data();
1424 const unsigned int num_cells = tree_index.size();
1425 std::memcpy(ptr, &num_cells,
sizeof(
unsigned int));
1426 ptr +=
sizeof(
unsigned int);
1430 num_cells*
sizeof(
unsigned int));
1431 ptr += num_cells*
sizeof(
unsigned int);
1435 num_cells *
sizeof(typename ::internal::p4est::
1436 types<dim>::quadrant));
1437 ptr += num_cells*
sizeof(typename ::internal::p4est::types<dim>::
1441 vertex_indices.data(),
1442 vertex_indices.size() *
sizeof(
unsigned int));
1443 ptr += vertex_indices.size() *
sizeof(
unsigned int);
1449 Assert (ptr == buffer.data()+buffer.size(),
1454 void unpack_data (
const std::vector<char> &buffer)
1456 const char *ptr = buffer.data();
1458 memcpy(&cells, ptr,
sizeof(
unsigned int));
1459 ptr +=
sizeof(
unsigned int);
1461 tree_index.resize(cells);
1462 memcpy(tree_index.data(),ptr,
sizeof(
unsigned int)*cells);
1463 ptr +=
sizeof(
unsigned int)*cells;
1465 quadrants.resize(cells);
1466 memcpy(quadrants.data(),ptr,
1467 sizeof(typename ::internal::p4est::types<dim>::quadrant)*cells);
1468 ptr +=
sizeof(typename ::internal::p4est::types<dim>::quadrant)*cells;
1470 vertex_indices.clear();
1471 first_vertex_indices.resize(cells);
1472 std::vector<unsigned int> n_vertices_on_cell(cells);
1473 std::vector<unsigned int> first_indices (cells);
1474 for (
unsigned int c=0; c<cells; ++c)
1479 const unsigned int *
const vertex_index
1480 =
reinterpret_cast<const unsigned int *
>(ptr);
1481 first_indices[c] = vertex_indices.size();
1482 vertex_indices.push_back(*vertex_index);
1483 n_vertices_on_cell[c] = *vertex_index;
1484 ptr +=
sizeof(
unsigned int);
1486 vertex_indices.resize(vertex_indices.size() + n_vertices_on_cell[c]);
1487 memcpy(&vertex_indices[vertex_indices.size() - n_vertices_on_cell[c]],
1488 ptr, n_vertices_on_cell[c]*
sizeof(
unsigned int));
1489 ptr += n_vertices_on_cell[c]*
sizeof(
unsigned int);
1491 for (
unsigned int c=0; c<cells; ++c)
1492 first_vertex_indices[c] = &vertex_indices[first_indices[c]];
1495 first_vertices.resize(cells);
1496 for (
unsigned int c=0; c<cells; ++c)
1498 first_indices[c] = vertices.size();
1499 vertices.resize(vertices.size() + n_vertices_on_cell[c]);
1500 memcpy(&vertices[vertices.size() - n_vertices_on_cell[c]],
1504 for (
unsigned int c=0; c<cells; ++c)
1505 first_vertices[c] = &vertices[first_indices[c]];
1507 Assert (ptr == buffer.data() + buffer.size(),
1526 template <
int dim,
int spacedim>
1529 const unsigned int tree_index,
1531 const typename ::internal::p4est::types<dim>::quadrant &p4est_cell,
1532 const std::map<
unsigned int, std::set<::types::subdomain_id> > &vertices_with_ghost_neighbors,
1533 const std::vector<bool> &vertex_locally_moved,
1538 if (dealii_cell->has_children())
1540 typename ::internal::p4est::types<dim>::quadrant
1542 ::internal::p4est::init_quadrant_children<dim>(p4est_cell, p4est_child);
1545 for (
unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
1546 fill_vertices_recursively<dim,spacedim>(tria,
1548 dealii_cell->child(c),
1550 vertices_with_ghost_neighbors,
1551 vertex_locally_moved,
1563 if (dealii_cell->is_locally_owned())
1565 std::set<::types::subdomain_id> send_to;
1566 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
1568 const std::map<unsigned int, std::set<::types::subdomain_id> >::const_iterator
1569 neighbor_subdomains_of_vertex
1570 = vertices_with_ghost_neighbors.find (dealii_cell->vertex_index(v));
1572 if (neighbor_subdomains_of_vertex
1573 != vertices_with_ghost_neighbors.end())
1575 Assert(neighbor_subdomains_of_vertex->second.size()!=0,
1577 send_to.insert(neighbor_subdomains_of_vertex->second.begin(),
1578 neighbor_subdomains_of_vertex->second.end());
1582 if (send_to.size() > 0)
1584 std::vector<unsigned int> vertex_indices;
1585 std::vector<::Point<spacedim> > local_vertices;
1586 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
1587 if (vertex_locally_moved[dealii_cell->vertex_index(v)])
1589 vertex_indices.push_back(v);
1590 local_vertices.push_back(dealii_cell->vertex(v));
1593 if (vertex_indices.size()>0)
1594 for (std::set<::types::subdomain_id>::iterator it=send_to.begin();
1595 it!=send_to.end(); ++it)
1597 const ::types::subdomain_id subdomain = *it;
1601 const typename std::map<::types::subdomain_id, CellInfo<dim, spacedim> >::iterator
1603 = needs_to_get_cell.insert (std::make_pair(subdomain,
1604 CellInfo<dim,spacedim>()))
1607 p->second.tree_index.push_back(tree_index);
1608 p->second.quadrants.push_back(p4est_cell);
1610 p->second.vertex_indices.push_back(vertex_indices.size());
1611 p->second.vertex_indices.insert(p->second.vertex_indices.end(),
1612 vertex_indices.begin(),
1613 vertex_indices.end());
1615 p->second.vertices.insert(p->second.vertices.end(),
1616 local_vertices.begin(),
1617 local_vertices.end());
1634 template <
int dim,
int spacedim>
1636 set_vertices_recursively (
1638 const typename ::internal::p4est::types<dim>::quadrant &p4est_cell,
1640 const typename ::internal::p4est::types<dim>::quadrant &quadrant,
1641 const ::Point<spacedim> *
const vertices,
1642 const unsigned int *
const vertex_indices)
1644 if (::internal::p4est::quadrant_is_equal<dim>(p4est_cell, quadrant))
1650 const unsigned int n_vertices = vertex_indices[0];
1653 for (
unsigned int i=0; i<n_vertices; ++i)
1654 dealii_cell->vertex(vertex_indices[i+1]) = vertices[i];
1659 if (! dealii_cell->has_children())
1662 if (! ::internal::p4est::quadrant_is_ancestor<dim> (p4est_cell, quadrant))
1665 typename ::internal::p4est::types<dim>::quadrant
1667 ::internal::p4est::init_quadrant_children<dim>(p4est_cell, p4est_child);
1669 for (
unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
1670 set_vertices_recursively<dim,spacedim> (tria, p4est_child[c],
1671 dealii_cell->child(c),
1680 template <
int dim,
int spacedim>
1684 triangulation_has_content =
false;
1686 if (parallel_ghost !=
nullptr)
1689 parallel_ghost =
nullptr;
1692 if (parallel_forest !=
nullptr)
1695 parallel_forest =
nullptr;
1698 if (connectivity !=
nullptr)
1701 connectivity =
nullptr;
1704 coarse_cell_to_p4est_tree_permutation.resize (0);
1705 p4est_tree_to_coarse_cell_permutation.resize (0);
1709 this->update_number_cache ();
1714 template <
int dim,
int spacedim>
1718 if (this->n_global_levels()<=1)
1727 bool have_coarser_cell =
false;
1729 cell != this->end(this->n_global_levels()-2);
1731 if (cell->is_locally_owned())
1733 have_coarser_cell =
true;
1743 template <
int dim,
int spacedim>
1748 ::GridTools::get_vertex_connectivity_of_cells (*
this, cell_connectivity);
1749 coarse_cell_to_p4est_tree_permutation.resize (this->n_cells(0));
1752 coarse_cell_to_p4est_tree_permutation);
1754 p4est_tree_to_coarse_cell_permutation
1760 template <
int dim,
int spacedim>
1764 Assert (parallel_forest !=
nullptr,
1765 ExcMessage (
"Can't produce output when no forest is created yet."));
1772 template <
int dim,
int spacedim>
1777 Assert(cell_attached_data.n_attached_deserialize==0,
1778 ExcMessage (
"not all SolutionTransfer's got deserialized after the last load()"));
1779 Assert(this->n_cells()>0,
ExcMessage(
"Can not save() an empty Triangulation."));
1784 if (this->my_subdomain==0)
1786 const unsigned int buffer_size_per_cell
1787 = (cell_attached_data.cumulative_fixed_data_size > 0
1789 cell_attached_data.cumulative_fixed_data_size+
sizeof(
CellStatus)
1793 std::string fname=std::string(filename)+
".info";
1794 std::ofstream f(fname.c_str());
1795 f <<
"version nproc attached_bytes n_attached_objs n_coarse_cells" << std::endl
1798 << buffer_size_per_cell <<
" " 1799 << cell_attached_data.pack_callbacks.size() <<
" " 1804 if (cell_attached_data.cumulative_fixed_data_size>0)
1807 ->attach_mesh_data();
1811 cell_attached_data.cumulative_fixed_data_size>0);
1819 tria->cell_attached_data.n_attached_datas = 0;
1820 tria->cell_attached_data.cumulative_fixed_data_size = 0;
1821 tria->cell_attached_data.pack_callbacks.clear();
1825 void *userptr = parallel_forest->user_pointer;
1827 parallel_forest->user_pointer = userptr;
1832 template <
int dim,
int spacedim>
1836 const bool autopartition)
1838 Assert(this->n_cells()>0,
ExcMessage(
"load() only works if the Triangulation already contains a coarse mesh!"));
1839 Assert(this->n_levels()==1,
ExcMessage(
"Triangulation may only contain coarse cells when calling load()."));
1842 if (parallel_ghost !=
nullptr)
1845 parallel_ghost =
nullptr;
1848 parallel_forest =
nullptr;
1850 connectivity =
nullptr;
1852 unsigned int version, numcpus, attached_size, attached_count, n_coarse_cells;
1854 std::string fname=std::string(filename)+
".info";
1855 std::ifstream f(fname.c_str());
1857 std::string firstline;
1858 getline(f, firstline);
1859 f >> version >> numcpus >> attached_size >> attached_count >> n_coarse_cells;
1862 Assert(version == 2,
ExcMessage(
"Incompatible version found in .info file."));
1863 Assert(this->n_cells(0) == n_coarse_cells,
ExcMessage(
"Number of coarse cells differ!"));
1864 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3) 1867 ExcMessage(
"parallel::distributed::Triangulation::load() only supports loading " 1868 "saved data with a greater or equal number of processes than were used to " 1869 "save() when using p4est 0.3.4.2."));
1874 cell_attached_data.cumulative_fixed_data_size = 0;
1875 cell_attached_data.n_attached_datas = 0;
1876 cell_attached_data.n_attached_deserialize = attached_count;
1878 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,3) 1880 filename, this->mpi_communicator,
1881 attached_size, attached_size>0,
1886 (void)autopartition;
1888 filename, this->mpi_communicator,
1889 attached_size, attached_size>0,
1903 copy_local_forest_to_triangulation ();
1914 this->update_number_cache ();
1919 this->update_periodic_face_map();
1924 template <
int dim,
int spacedim>
1928 Assert (parallel_forest !=
nullptr,
1929 ExcMessage (
"Can't produce a check sum when no forest is created yet."));
1930 return ::internal::p4est::functions<dim>::checksum (parallel_forest);
1935 template <
int dim,
int spacedim>
1942 if (settings & construct_multigrid_hierarchy)
1948 template <
int dim,
int spacedim>
1949 typename ::internal::p4est::types<dim>::tree *
1953 const unsigned int tree_index
1954 = coarse_cell_to_p4est_tree_permutation[dealii_coarse_cell_index];
1955 typename ::internal::p4est::types<dim>::tree *tree
1956 =
static_cast<typename ::internal::p4est::types<dim>::tree *
> 1957 (sc_array_index (parallel_forest->trees,
1969 const unsigned int dim = 2, spacedim = 2;
1977 std::vector<unsigned int> vertex_touch_count;
1980 std::pair<Triangulation<dim,spacedim>::active_cell_iterator,
1983 get_vertex_to_cell_mappings (*
this,
1986 const ::internal::p4est::types<2>::locidx
1987 num_vtt = std::accumulate (vertex_touch_count.begin(),
1988 vertex_touch_count.end(),
1994 const bool set_vertex_info
2009 set_vertex_and_cell_info (*
this,
2012 coarse_cell_to_p4est_tree_permutation,
2016 Assert (p4est_connectivity_is_valid (connectivity) == 1,
2040 const unsigned int dim = 2, spacedim = 3;
2048 std::vector<unsigned int> vertex_touch_count;
2051 std::pair<Triangulation<dim,spacedim>::active_cell_iterator,
2054 get_vertex_to_cell_mappings (*
this,
2057 const ::internal::p4est::types<2>::locidx
2058 num_vtt = std::accumulate (vertex_touch_count.begin(),
2059 vertex_touch_count.end(),
2065 const bool set_vertex_info
2080 set_vertex_and_cell_info (*
this,
2083 coarse_cell_to_p4est_tree_permutation,
2087 Assert (p4est_connectivity_is_valid (connectivity) == 1,
2109 const int dim = 3, spacedim = 3;
2117 std::vector<unsigned int> vertex_touch_count;
2120 std::pair<Triangulation<3>::active_cell_iterator,
2123 get_vertex_to_cell_mappings (*
this,
2126 const ::internal::p4est::types<2>::locidx
2127 num_vtt = std::accumulate (vertex_touch_count.begin(),
2128 vertex_touch_count.end(),
2131 std::vector<unsigned int> edge_touch_count;
2134 std::pair<Triangulation<3>::active_cell_iterator,
2137 get_edge_to_cell_mappings (*
this,
2140 const ::internal::p4est::types<2>::locidx
2141 num_ett = std::accumulate (edge_touch_count.begin(),
2142 edge_touch_count.end(),
2146 const bool set_vertex_info
2158 this->n_active_lines(),
2163 set_vertex_and_cell_info (*
this,
2166 coarse_cell_to_p4est_tree_permutation,
2195 const unsigned int deal_to_p4est_line_index[12]
2196 = { 4, 5, 0, 1, 6, 7, 2, 3, 8, 9, 10, 11 } ;
2199 cell = this->begin_active();
2200 cell != this->end(); ++cell)
2203 index = coarse_cell_to_p4est_tree_permutation[cell->index()];
2204 for (
unsigned int e=0; e<GeometryInfo<3>::lines_per_cell; ++
e)
2206 deal_to_p4est_line_index[e]]
2207 = cell->line(e)->index();
2212 connectivity->ett_offset[0] = 0;
2213 std::partial_sum (edge_touch_count.begin(),
2214 edge_touch_count.end(),
2215 &connectivity->ett_offset[1]);
2217 Assert (connectivity->ett_offset[this->n_active_lines()] ==
2221 for (
unsigned int v=0; v<this->n_active_lines(); ++v)
2223 Assert (edge_to_cell[v].size() == edge_touch_count[v],
2228 unsigned int> >::const_iterator
2229 p = edge_to_cell[v].begin();
2230 for (
unsigned int c=0; c<edge_touch_count[v]; ++c, ++p)
2232 connectivity->edge_to_tree[connectivity->ett_offset[v]+c]
2233 = coarse_cell_to_p4est_tree_permutation[p->first->index()];
2234 connectivity->edge_to_edge[connectivity->ett_offset[v]+c]
2235 = deal_to_p4est_line_index[p->second];
2239 Assert (p8est_connectivity_is_valid (connectivity) == 1,
2261 template <
int dim,
int spacedim>
2262 bool enforce_mesh_balance_over_periodic_boundaries
2268 std::vector<bool> flags_before[2];
2272 std::vector<unsigned int> topological_vertex_numbering(tria.
n_vertices());
2273 for (
unsigned int i=0; i<topological_vertex_numbering.size(); ++i)
2274 topological_vertex_numbering[i] = i;
2291 typename std::map<std::pair<cell_iterator, unsigned int>,
2292 std::pair<std::pair<cell_iterator,unsigned int>, std::bitset<3> > >::const_iterator it;
2295 const cell_iterator &cell_1 = it->first.first;
2296 const unsigned int face_no_1 = it->first.second;
2297 const cell_iterator &cell_2 = it->second.first.first;
2298 const unsigned int face_no_2 = it->second.first.second;
2299 const std::bitset<3> face_orientation = it->second.second;
2301 if (cell_1->level() == cell_2->level())
2303 for (
unsigned int v=0; v<
GeometryInfo<dim-1>::vertices_per_cell; ++v)
2307 const unsigned int vface0 =
2309 face_orientation[1],
2310 face_orientation[2]);
2311 const unsigned int vi0 = topological_vertex_numbering[cell_1->face(face_no_1)->vertex_index(vface0)];
2312 const unsigned int vi1 = topological_vertex_numbering[cell_2->face(face_no_2)->vertex_index(v)];
2313 const unsigned int min_index = std::min(vi0, vi1);
2314 topological_vertex_numbering[cell_1->face(face_no_1)->vertex_index(vface0)]
2315 = topological_vertex_numbering[cell_2->face(face_no_2)->vertex_index(v)]
2322 for (
unsigned int i=0; i<topological_vertex_numbering.size(); ++i)
2324 const unsigned int j = topological_vertex_numbering[i];
2326 Assert(topological_vertex_numbering[j] == j,
2333 bool continue_iterating =
true;
2334 std::vector<int> vertex_level(tria.
n_vertices());
2335 while (continue_iterating)
2339 std::fill (vertex_level.begin(), vertex_level.end(), 0);
2342 for (; cell!=endc; ++cell)
2344 if (cell->refine_flag_set())
2345 for (
unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_cell;
2347 vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]]
2348 = std::max (vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]],
2350 else if (!cell->coarsen_flag_set())
2351 for (
unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_cell;
2353 vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]]
2354 = std::max (vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]],
2365 for (
unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_cell;
2367 vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]]
2368 = std::max (vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]],
2373 continue_iterating =
false;
2384 for (cell=tria.
last_active(); cell != endc; --cell)
2385 if (cell->refine_flag_set() ==
false)
2387 for (
unsigned int vertex=0;
2388 vertex<GeometryInfo<dim>::vertices_per_cell; ++vertex)
2389 if (vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]] >=
2393 cell->clear_coarsen_flag();
2398 if (vertex_level[topological_vertex_numbering[cell->vertex_index(vertex)]] >
2401 cell->set_refine_flag();
2402 continue_iterating =
true;
2404 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell;
2406 vertex_level[topological_vertex_numbering[cell->vertex_index(v)]]
2407 = std::max (vertex_level[topological_vertex_numbering[cell->vertex_index(v)]],
2419 cell!=tria.
end(); ++cell)
2425 const unsigned int n_children=cell->n_children();
2426 unsigned int flagged_children=0;
2427 for (
unsigned int child=0; child<n_children; ++child)
2428 if (cell->child(child)->active() &&
2429 cell->child(child)->coarsen_flag_set())
2434 if (flagged_children < n_children)
2435 for (
unsigned int child=0; child<n_children; ++child)
2436 if (cell->child(child)->active())
2437 cell->child(child)->clear_coarsen_flag();
2440 std::vector<bool> flags_after[2];
2443 return ((flags_before[0] != flags_after[0]) ||
2444 (flags_before[1] != flags_after[1]));
2450 template <
int dim,
int spacedim>
2454 std::vector<bool> flags_before[2];
2455 this->save_coarsen_flags (flags_before[0]);
2456 this->save_refine_flags (flags_before[1]);
2458 bool mesh_changed =
false;
2462 this->update_periodic_face_map();
2464 if (this->smooth_grid &
2466 mesh_changed = enforce_mesh_balance_over_periodic_boundaries(*
this);
2468 while (mesh_changed);
2472 std::vector<bool> flags_after[2];
2473 this->save_coarsen_flags (flags_after[0]);
2474 this->save_refine_flags (flags_after[1]);
2475 return ((flags_before[0] != flags_after[0]) ||
2476 (flags_before[1] != flags_after[1]));
2481 template <
int dim,
int spacedim>
2491 save_smooth = this->smooth_grid;
2501 if (settings & construct_multigrid_hierarchy)
2506 bool mesh_changed =
false;
2514 if (settings & mesh_reconstruction_after_repartitioning)
2515 while (this->begin_active()->level() > 0)
2518 cell = this->begin_active();
2519 cell != this->end();
2522 cell->set_coarsen_flag();
2525 this->prepare_coarsening_and_refinement();
2540 if (parallel_ghost !=
nullptr)
2543 parallel_ghost =
nullptr;
2548 typename ::internal::p4est::types<dim>::
2549 balance_type(P4EST_CONNECT_CORNER)
2551 typename ::internal::p4est::types<dim>::
2552 balance_type(P8EST_CONNECT_CORNER)));
2560 cell = this->begin(0);
2561 cell != this->end(0);
2568 cell = this->begin(0);
2569 cell != this->end(0);
2575 if (tree_exists_locally<dim,spacedim>(parallel_forest,
2576 coarse_cell_to_p4est_tree_permutation[cell->index()])
2579 delete_all_children<dim,spacedim> (cell);
2580 if (!cell->has_children())
2589 typename ::internal::p4est::types<dim>::quadrant p4est_coarse_cell;
2590 typename ::internal::p4est::types<dim>::tree *tree =
2591 init_tree(cell->index());
2593 ::internal::p4est::init_coarse_quadrant<dim>(p4est_coarse_cell);
2595 match_tree_recursively<dim,spacedim> (*tree, cell,
2598 this->my_subdomain);
2605 typename ::internal::p4est::types<dim>::quadrant *quadr;
2607 typename ::internal::p4est::types<dim>::topidx ghost_tree=0;
2609 for (
unsigned int g_idx=0; g_idx<parallel_ghost->ghosts.elem_count; ++g_idx)
2611 while (g_idx >= (
unsigned int)parallel_ghost->proc_offsets[ghost_owner+1])
2613 while (g_idx >= (
unsigned int)parallel_ghost->tree_offsets[ghost_tree+1])
2616 quadr =
static_cast<typename ::internal::p4est::types<dim>::quadrant *
> 2617 ( sc_array_index(¶llel_ghost->ghosts, g_idx) );
2619 unsigned int coarse_cell_index =
2620 p4est_tree_to_coarse_cell_permutation[ghost_tree];
2622 match_quadrant<dim,spacedim> (
this, coarse_cell_index, *quadr, ghost_owner);
2626 this->prepare_coarsening_and_refinement ();
2629 mesh_changed =
false;
2631 cell = this->begin_active();
2632 cell != this->end();
2634 if (cell->refine_flag_set() || cell->coarsen_flag_set())
2636 mesh_changed =
true;
2653 while (mesh_changed);
2657 unsigned int num_ghosts = 0;
2660 cell = this->begin_active();
2661 cell != this->end();
2664 if (cell->subdomain_id() != this->my_subdomain
2678 if (settings & construct_multigrid_hierarchy)
2684 for (
unsigned int lvl=this->n_levels(); lvl>0; )
2688 for (cell = this->begin(lvl); cell!=endc; ++cell)
2690 if ((!cell->has_children() && cell->subdomain_id()==this->locally_owned_subdomain())
2691 || (cell->has_children() && cell->child(0)->level_subdomain_id()==this->locally_owned_subdomain()))
2692 cell->set_level_subdomain_id(this->locally_owned_subdomain());
2703 std::vector<std::vector<bool> > marked_vertices(this->n_levels());
2704 for (
unsigned int lvl=0; lvl < this->n_levels(); ++lvl)
2705 marked_vertices[lvl] = mark_locally_active_vertices_on_level(lvl);
2709 typename ::internal::p4est::types<dim>::quadrant p4est_coarse_cell;
2710 const unsigned int tree_index
2711 = coarse_cell_to_p4est_tree_permutation[cell->index()];
2712 typename ::internal::p4est::types<dim>::tree *tree =
2713 init_tree(cell->index());
2715 ::internal::p4est::init_coarse_quadrant<dim>(p4est_coarse_cell);
2717 determine_level_subdomain_id_recursively<dim,spacedim> (*tree, tree_index, cell,
2725 for (
unsigned int lvl=this->n_levels(); lvl>0;)
2729 for (cell = this->begin(lvl); cell!=endc; ++cell)
2731 if (cell->has_children())
2732 for (
unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
2734 if (cell->child(c)->level_subdomain_id()==this->locally_owned_subdomain())
2739 cell->child(0)->level_subdomain_id();
2741 cell->set_level_subdomain_id(mark);
2758 const unsigned int total_local_cells = this->n_active_cells();
2759 (void)total_local_cells;
2762 Assert (static_cast<unsigned int>(parallel_forest->local_num_quadrants) ==
2766 Assert (static_cast<unsigned int>(parallel_forest->local_num_quadrants) <=
2771 unsigned int n_owned = 0;
2773 cell = this->begin_active();
2774 cell != this->end(); ++cell)
2776 if (cell->subdomain_id() == this->my_subdomain)
2780 Assert(static_cast<unsigned int>(parallel_forest->local_num_quadrants) ==
2785 this->smooth_grid = save_smooth;
2790 template <
int dim,
int spacedim>
2797 cell = this->begin_active();
2798 cell != this->end(); ++cell)
2799 if (cell->is_locally_owned() && cell->refine_flag_set())
2800 Assert (cell->refine_flag_set() ==
2802 ExcMessage (
"This class does not support anisotropic refinement"));
2814 ExcMessage(
"Fatal Error: maximum refinement level of p4est reached."));
2822 this->prepare_coarsening_and_refinement ();
2827 cell = this->begin_active();
2828 cell != this->end(); ++cell)
2829 if (cell->is_ghost() || cell->is_artificial())
2831 cell->clear_refine_flag ();
2832 cell->clear_coarsen_flag ();
2838 RefineAndCoarsenList<dim,spacedim>
2839 refine_and_coarsen_list (*
this,
2840 p4est_tree_to_coarse_cell_permutation,
2841 this->my_subdomain);
2847 Assert (parallel_forest->user_pointer ==
this,
2849 parallel_forest->user_pointer = &refine_and_coarsen_list;
2851 if (parallel_ghost !=
nullptr)
2854 parallel_ghost =
nullptr;
2857 refine (parallel_forest,
false,
2858 &RefineAndCoarsenList<dim,spacedim>::refine_callback,
2861 coarsen (parallel_forest,
false,
2862 &RefineAndCoarsenList<dim,spacedim>::coarsen_callback,
2866 Assert (refine_and_coarsen_list.pointers_are_at_end(),
2870 parallel_forest->user_pointer =
this;
2878 typename ::internal::p4est::types<dim>::
2879 balance_type(P4EST_CONNECT_FULL)
2881 typename ::internal::p4est::types<dim>::
2882 balance_type(P8EST_CONNECT_FULL)),
2889 if (!(settings & no_automatic_repartitioning))
2895 partition (parallel_forest,
2901 const std::vector<unsigned int> cell_weights = get_cell_weights();
2903 PartitionWeights<dim,spacedim> partition_weights (cell_weights);
2907 Assert (parallel_forest->user_pointer ==
this,
2909 parallel_forest->user_pointer = &partition_weights;
2914 &PartitionWeights<dim,spacedim>::cell_weight);
2917 parallel_forest->user_pointer =
this;
2925 cell = this->begin_active();
2926 cell != this->end(); ++cell)
2928 cell->clear_refine_flag();
2929 cell->clear_coarsen_flag();
2934 copy_local_forest_to_triangulation ();
2959 if (settings & construct_multigrid_hierarchy)
2961 for (
unsigned int lvl=0; lvl<this->n_global_levels(); ++lvl)
2963 std::vector<bool> active_verts = this->mark_locally_active_vertices_on_level(lvl);
2965 const unsigned int maybe_coarser_lvl = (lvl>0) ? (lvl-1) : lvl;
2967 endc = this->end(lvl);
2968 for (; cell != endc; ++cell)
2969 if (cell->level() ==
static_cast<int>(lvl) || cell->active())
2971 const bool is_level_artificial =
2973 bool need_to_know =
false;
2974 for (
unsigned int vertex=0; vertex<GeometryInfo<dim>::vertices_per_cell;
2976 if (active_verts[cell->vertex_index(vertex)])
2978 need_to_know =
true;
2982 Assert(!need_to_know || !is_level_artificial,
2983 ExcMessage(
"Internal error: the owner of cell" 2984 + cell->id().to_string()
2985 +
" is unknown even though it is needed for geometric multigrid."));
2991 this->update_number_cache ();
2998 this->update_periodic_face_map();
3003 template <
int dim,
int spacedim>
3010 cell = this->begin_active();
3011 cell != this->end(); ++cell)
3012 if (cell->is_locally_owned())
3014 !cell->refine_flag_set() && !cell->coarsen_flag_set(),
3015 ExcMessage (
"Error: There shouldn't be any cells flagged for coarsening/refinement when calling repartition()."));
3034 const std::vector<unsigned int> cell_weights = get_cell_weights();
3036 PartitionWeights<dim,spacedim> partition_weights (cell_weights);
3040 Assert (parallel_forest->user_pointer ==
this,
3042 parallel_forest->user_pointer = &partition_weights;
3047 &PartitionWeights<dim,spacedim>::cell_weight);
3050 parallel_forest->user_pointer =
this;
3055 copy_local_forest_to_triangulation ();
3065 this->update_number_cache ();
3066 this->update_periodic_face_map();
3071 template <
int dim,
int spacedim>
3076 Assert (vertex_locally_moved.size() == this->n_vertices(),
3078 this->n_vertices()));
3081 const std::vector<bool> locally_owned_vertices
3082 = ::GridTools::get_locally_owned_vertices (*
this);
3083 for (
unsigned int i=0; i<locally_owned_vertices.size(); ++i)
3084 Assert ((vertex_locally_moved[i] ==
false)
3086 (locally_owned_vertices[i] ==
true),
3087 ExcMessage (
"The vertex_locally_moved argument must not " 3088 "contain vertices that are not locally owned"));
3098 const std::map<unsigned int, std::set<::types::subdomain_id> >
3099 vertices_with_ghost_neighbors = compute_vertices_with_ghost_neighbors ();
3104 std::map<::types::subdomain_id, CommunicateLocallyMovedVertices::CellInfo<dim,spacedim> > cellmap_t;
3105 cellmap_t needs_to_get_cells;
3108 cell = this->begin(0);
3109 cell != this->end(0);
3112 typename ::internal::p4est::types<dim>::quadrant p4est_coarse_cell;
3113 ::internal::p4est::init_coarse_quadrant<dim>(p4est_coarse_cell);
3115 CommunicateLocallyMovedVertices::fill_vertices_recursively<dim,spacedim>
3117 this->get_coarse_cell_to_p4est_tree_permutation()[cell->index()],
3120 vertices_with_ghost_neighbors,
3121 vertex_locally_moved,
3122 needs_to_get_cells);
3126 std::vector<std::vector<char> > sendbuffers (needs_to_get_cells.size());
3127 std::vector<std::vector<char> >::iterator buffer = sendbuffers.begin();
3128 std::vector<MPI_Request> requests (needs_to_get_cells.size());
3129 std::vector<unsigned int> destinations;
3133 for (
typename cellmap_t::iterator it=needs_to_get_cells.begin();
3134 it!=needs_to_get_cells.end();
3135 ++it, ++buffer, ++idx)
3137 const unsigned int num_cells = it->second.tree_index.size();
3139 destinations.push_back(it->first);
3151 it->second.pack_data (*buffer);
3152 const int ierr = MPI_Isend(&(*buffer)[0], buffer->size(),
3153 MPI_BYTE, it->first,
3154 123, this->get_communicator(), &requests[idx]);
3162 const std::vector<unsigned int> senders
3164 (this->get_communicator(), destinations);
3167 std::vector<char> receive;
3168 CommunicateLocallyMovedVertices::CellInfo<dim,spacedim> cellinfo;
3169 for (
unsigned int i=0; i<senders.size(); ++i)
3173 int ierr = MPI_Probe(MPI_ANY_SOURCE, 123, this->get_communicator(), &status);
3175 ierr = MPI_Get_count(&status, MPI_BYTE, &len);
3177 receive.resize(len);
3179 char *ptr = receive.data();
3180 ierr = MPI_Recv(ptr, len, MPI_BYTE, status.MPI_SOURCE, status.MPI_TAG,
3181 this->get_communicator(), &status);
3184 cellinfo.unpack_data(receive);
3185 const unsigned int cells = cellinfo.tree_index.size();
3186 for (
unsigned int c=0; c<cells; ++c)
3188 typename ::parallel::distributed::Triangulation<dim,spacedim>::cell_iterator
3191 this->get_p4est_tree_to_coarse_cell_permutation()[cellinfo.tree_index[c]]);
3193 typename ::internal::p4est::types<dim>::quadrant p4est_coarse_cell;
3194 ::internal::p4est::init_coarse_quadrant<dim>(p4est_coarse_cell);
3196 CommunicateLocallyMovedVertices::set_vertices_recursively<dim,spacedim> (*
this,
3199 cellinfo.quadrants[c],
3200 cellinfo.first_vertices[c],
3201 cellinfo.first_vertex_indices[c]);
3207 if (requests.size() > 0)
3209 const int ierr = MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
3221 template <
int dim,
int spacedim>
3227 void *)> &pack_callback)
3230 Assert(cell_attached_data.pack_callbacks.size()==cell_attached_data.n_attached_datas,
3231 ExcMessage(
"register_data_attach(), not all data has been unpacked last time?"));
3233 const unsigned int current_buffer_positition = cell_attached_data.cumulative_fixed_data_size+
sizeof(
CellStatus);
3234 ++cell_attached_data.n_attached_datas;
3235 cell_attached_data.cumulative_fixed_data_size += size;
3236 cell_attached_data.pack_callbacks[current_buffer_positition] = pack_callback;
3240 return current_buffer_positition;
3245 template <
int dim,
int spacedim>
3251 const void *)> &unpack_callback)
3255 Assert (handle < cell_attached_data.cumulative_fixed_data_size +
sizeof(
CellStatus),
3257 Assert (cell_attached_data.n_attached_datas > 0,
3258 ExcMessage (
"The notify_ready_to_unpack() has already been called " 3259 "once for each registered callback."));
3263 cell = this->begin (0);
3264 cell != this->end (0);
3268 if (tree_exists_locally<dim, spacedim> (parallel_forest,
3269 coarse_cell_to_p4est_tree_permutation[cell->index() ])
3273 typename ::internal::p4est::types<dim>::quadrant p4est_coarse_cell;
3274 typename ::internal::p4est::types<dim>::tree *tree =
3275 init_tree (cell->index());
3277 ::internal::p4est::init_coarse_quadrant<dim> (p4est_coarse_cell);
3281 post_mesh_data_recursively<dim, spacedim> (*tree,
3289 --cell_attached_data.n_attached_datas;
3290 if (cell_attached_data.n_attached_deserialize > 0)
3292 --cell_attached_data.n_attached_deserialize;
3297 const unsigned int n_elements_removed
3298 = cell_attached_data.pack_callbacks.erase (handle);
3299 (void)n_elements_removed;
3300 Assert (n_elements_removed == 1,
3310 if (cell_attached_data.n_attached_datas == 0
3312 cell_attached_data.n_attached_deserialize == 0)
3315 cell_attached_data.cumulative_fixed_data_size = 0;
3316 cell_attached_data.pack_callbacks.clear();
3319 void *userptr = parallel_forest->user_pointer;
3321 parallel_forest->user_pointer = userptr;
3327 template <
int dim,
int spacedim>
3328 const std::vector<types::global_dof_index> &
3331 return p4est_tree_to_coarse_cell_permutation;
3336 template <
int dim,
int spacedim>
3337 const std::vector<types::global_dof_index> &
3340 return coarse_cell_to_p4est_tree_permutation;
3345 template <
int dim,
int spacedim>
3346 std::map<unsigned int, std::set<::types::subdomain_id> >
3351 return ::internal::p4est::compute_vertices_with_ghost_neighbors<dim, spacedim> (*
this,
3352 this->parallel_forest,
3353 this->parallel_ghost);
3358 template <
int dim,
int spacedim>
3365 std::vector<bool> marked_vertices(this->n_vertices(),
false);
3367 endc = this->end(level);
3368 for ( ; cell != endc; ++cell)
3369 if (cell->level_subdomain_id() == this->locally_owned_subdomain())
3370 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
3371 marked_vertices[cell->vertex_index(v)] =
true;
3378 typename std::map<std::pair<cell_iterator, unsigned int>,
3379 std::pair<std::pair<cell_iterator,unsigned int>, std::bitset<3> > >::const_iterator it;
3391 for (
unsigned int repetition=0; repetition<dim; ++repetition)
3392 for (it = this->get_periodic_face_map().begin(); it!= this->get_periodic_face_map().end(); ++it)
3395 const unsigned int face_no_1 = it->first.second;
3397 const unsigned int face_no_2 = it->second.first.second;
3398 const std::bitset<3> &face_orientation = it->second.second;
3400 if (cell_1->level() == level &&
3401 cell_2->level() == level)
3403 for (
unsigned int v=0; v<
GeometryInfo<dim-1>::vertices_per_cell; ++v)
3406 const unsigned int vface0 =
3408 face_orientation[1],
3409 face_orientation[2]);
3410 if (marked_vertices[cell_1->face(face_no_1)->vertex_index(vface0)] ||
3411 marked_vertices[cell_2->face(face_no_2)->vertex_index(v)])
3412 marked_vertices[cell_1->face(face_no_1)->vertex_index(vface0)]
3413 = marked_vertices[cell_2->face(face_no_2)->vertex_index(v)]
3419 return marked_vertices;
3424 template <
int dim,
int spacedim>
3430 #if DEAL_II_P4EST_VERSION_GTE(0,3,4,1) 3431 Assert (triangulation_has_content ==
true,
3433 Assert (this->n_levels() == 1,
3434 ExcMessage (
"The triangulation is refined!"));
3436 typedef std::vector<::GridTools::PeriodicFacePair<cell_iterator> >
3438 typename FaceVector::const_iterator it, periodic_end;
3439 it = periodicity_vector.begin();
3440 periodic_end = periodicity_vector.end();
3442 for (; it<periodic_end; ++it)
3446 const unsigned int face_left = it->face_idx[0];
3447 const unsigned int face_right = it->face_idx[1];
3450 const unsigned int tree_left
3451 = coarse_cell_to_p4est_tree_permutation[std::distance(this->begin(),
3453 const unsigned int tree_right
3454 = coarse_cell_to_p4est_tree_permutation[std::distance(this->begin(),
3464 unsigned int p4est_orientation = 0;
3466 p4est_orientation = it->orientation[1];
3469 const unsigned int face_idx_list[] = {face_left, face_right};
3470 const cell_iterator cell_list[] = {first_cell, second_cell};
3471 unsigned int lower_idx, higher_idx;
3472 if (face_left<=face_right)
3484 unsigned int first_p4est_idx_on_cell = p8est_face_corners[face_idx_list[lower_idx]][0];
3486 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_face; ++i)
3488 const unsigned int first_dealii_idx_on_cell
3490 (face_idx_list[lower_idx], i,
3491 cell_list[lower_idx]->face_orientation(face_idx_list[lower_idx]),
3492 cell_list[lower_idx]->face_flip(face_idx_list[lower_idx]),
3493 cell_list[lower_idx]->face_rotation(face_idx_list[lower_idx]));
3494 if (first_p4est_idx_on_cell == first_dealii_idx_on_cell)
3496 first_dealii_idx_on_face = i;
3502 const unsigned int left_to_right [8][4] = {{0,2,1,3},{0,1,2,3},{3,1,2,0},{3,2,1,0},
3503 {2,3,0,1},{1,3,0,2},{1,0,3,2},{2,0,3,1}
3505 const unsigned int right_to_left [8][4] = {{0,2,1,3},{0,1,2,3},{3,1,2,0},{3,2,1,0},
3506 {2,3,0,1},{2,0,3,1},{1,0,3,2},{1,3,0,2}
3508 const unsigned int second_dealii_idx_on_face
3509 = lower_idx==0?left_to_right[it->orientation.to_ulong()][first_dealii_idx_on_face]:
3510 right_to_left[it->orientation.to_ulong()][first_dealii_idx_on_face];
3511 const unsigned int second_dealii_idx_on_cell
3513 (face_idx_list[higher_idx], second_dealii_idx_on_face,
3514 cell_list[higher_idx]->face_orientation(face_idx_list[higher_idx]),
3515 cell_list[higher_idx]->face_flip(face_idx_list[higher_idx]),
3516 cell_list[higher_idx]->face_rotation(face_idx_list[higher_idx]));
3518 const unsigned int second_p4est_idx_on_face
3519 = p8est_corner_face_corners[second_dealii_idx_on_cell][face_idx_list[higher_idx]];
3520 p4est_orientation = second_p4est_idx_on_face;
3552 copy_local_forest_to_triangulation ();
3565 this->update_number_cache();
3573 template <
int dim,
int spacedim>
3587 + memory_consumption_p4est();
3594 template <
int dim,
int spacedim>
3598 return ::internal::p4est::functions<dim>::forest_memory_used(parallel_forest)
3604 template <
int dim,
int spacedim>
3613 catch (
const typename ::Triangulation<dim,spacedim>::DistortedCellList &)
3623 triangulation_has_content =
true;
3625 Assert (other_tria.n_levels() == 1,
3626 ExcMessage (
"Parallel distributed triangulations can only be copied, " 3627 "if they are not refined!"));
3629 if (const ::parallel::distributed::Triangulation<dim,spacedim> *
3630 other_tria_x =
dynamic_cast<const ::parallel::distributed::Triangulation<dim,spacedim> *
>(&other_tria))
3632 coarse_cell_to_p4est_tree_permutation = other_tria_x->coarse_cell_to_p4est_tree_permutation;
3633 p4est_tree_to_coarse_cell_permutation = other_tria_x->p4est_tree_to_coarse_cell_permutation;
3634 cell_attached_data.cumulative_fixed_data_size = other_tria_x->cell_attached_data.cumulative_fixed_data_size;
3635 cell_attached_data.n_attached_datas = other_tria_x->cell_attached_data.n_attached_datas;
3637 settings = other_tria_x->settings;
3641 setup_coarse_cell_to_p4est_tree_permutation ();
3644 copy_new_triangulation_to_p4est (std::integral_constant<int, dim>());
3648 copy_local_forest_to_triangulation ();
3657 this->update_number_cache ();
3658 this->update_periodic_face_map();
3663 template <
int dim,
int spacedim>
3669 if (cell_attached_data.cumulative_fixed_data_size==0)
3676 void *userptr = parallel_forest->user_pointer;
3678 cell_attached_data.cumulative_fixed_data_size+
sizeof(
CellStatus),
3680 parallel_forest->user_pointer = userptr;
3687 cell = this->begin(0);
3688 cell != this->end(0);
3692 if (tree_exists_locally<dim,spacedim>(parallel_forest,
3693 coarse_cell_to_p4est_tree_permutation[cell->index()])
3697 typename ::internal::p4est::types<dim>::quadrant p4est_coarse_cell;
3698 typename ::internal::p4est::types<dim>::tree *tree =
3699 init_tree(cell->index());
3701 ::internal::p4est::init_coarse_quadrant<dim>(p4est_coarse_cell);
3703 attach_mesh_data_recursively<dim,spacedim>(*tree,
3706 cell_attached_data.pack_callbacks);
3712 template <
int dim,
int spacedim>
3713 std::vector<unsigned int>
3722 std::vector<unsigned int> weights;
3723 weights.reserve(this->n_active_cells());
3731 for (
unsigned int c=0; c<this->n_cells(0); ++c)
3734 if (tree_exists_locally<dim,spacedim>(parallel_forest,c) ==
false)
3737 const unsigned int coarse_cell_index =
3738 p4est_tree_to_coarse_cell_permutation[c];
3741 dealii_coarse_cell (
this, 0, coarse_cell_index);
3743 typename ::internal::p4est::types<dim>::quadrant p4est_coarse_cell;
3748 p4est_coarse_cell.p.which_tree = c;
3750 const typename ::internal::p4est::types<dim>::tree *tree =
3751 init_tree(coarse_cell_index);
3753 get_cell_weights_recursively<dim,spacedim>(*tree,
3765 template <
int spacedim>
3767 const typename ::Triangulation<1,spacedim>::MeshSmoothing smooth_grid,
3778 template <
int spacedim>
3786 template <
int spacedim>
3789 (
const std::vector<bool> &)
3796 template <
int spacedim>
3799 const std::function<
void (
const typename ::Triangulation<1,spacedim>::cell_iterator &,
3800 const typename ::Triangulation<1,spacedim>::CellStatus,
3809 template <
int spacedim>
3812 const std::function<
void (
const typename ::Triangulation<1,spacedim>::cell_iterator &,
3813 const typename ::Triangulation<1,spacedim>::CellStatus,
3821 template <
int spacedim>
3822 const std::vector<types::global_dof_index> &
3825 static std::vector<types::global_dof_index> a;
3831 template <
int spacedim>
3832 std::map<unsigned int, std::set<::types::subdomain_id> >
3837 return std::map<unsigned int, std::set<::types::subdomain_id> >();
3842 template <
int spacedim>
3843 std::map<unsigned int, std::set<::types::subdomain_id> >
3849 return std::map<unsigned int, std::set<::types::subdomain_id> >();
3854 template <
int spacedim>
3860 return std::vector<bool>();
3865 template <
int spacedim>
3875 template <
int spacedim>
3886 #endif // DEAL_II_WITH_P4EST 3892 #include "tria.inst" 3895 DEAL_II_NAMESPACE_CLOSE
virtual void copy_triangulation(const Triangulation< dim, spacedim > &other_tria)
unsigned int n_vertices() const
static const unsigned int invalid_unsigned_int
#define AssertNothrow(cond, exc)
virtual bool has_hanging_nodes() const
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)
static ::ExceptionBase & ExcIO()
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
virtual void copy_triangulation(const ::Triangulation< dim, spacedim > &old_tria)
boost::signals2::signal< unsigned int(const cell_iterator &, const CellStatus), CellWeightSum< unsigned int > > cell_weight
void fill_level_ghost_owners()
virtual std::map< unsigned int, std::set<::types::subdomain_id > > compute_vertices_with_ghost_neighbors() const
active_cell_iterator begin_active(const unsigned int level=0) const
#define AssertThrow(cond, exc)
void load(const char *filename, const bool autopartition=true)
cell_iterator begin(const unsigned int level=0) const
const std::vector< types::global_dof_index > & get_p4est_tree_to_coarse_cell_permutation() const
cell_iterator end() const
unsigned int register_data_attach(const std::size_t size, const std::function< void(const cell_iterator &, const CellStatus, void *)> &pack_callback)
virtual void execute_coarsening_and_refinement()
virtual bool prepare_coarsening_and_refinement()
static ::ExceptionBase & ExcMessage(std::string arg1)
boost::signals2::signal< void()> pre_distributed_save
Triangulation(const MeshSmoothing smooth_grid=none, const bool check_for_distorted_cells=false)
T sum(const T &t, const MPI_Comm &mpi_communicator)
#define Assert(cond, exc)
virtual void update_number_cache()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void load(Archive &ar, const unsigned int version)
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
void save_coarsen_flags(std::ostream &out) const
virtual std::size_t memory_consumption() const
void save_refine_flags(std::ostream &out) const
void notify_ready_to_unpack(const unsigned int handle, const std::function< void(const cell_iterator &, const CellStatus, const void *)> &unpack_callback)
void save(const char *filename) const
virtual std::size_t memory_consumption() const
boost::signals2::signal< void()> post_distributed_refinement
void save(Archive &ar, const unsigned int version) const
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)
std::vector< unsigned int > invert_permutation(const std::vector< unsigned int > &permutation)
unsigned int subdomain_id
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
::Triangulation< dim, spacedim >::active_cell_iterator active_cell_iterator
const types::subdomain_id artificial_subdomain_id
void communicate_locally_moved_vertices(const std::vector< bool > &vertex_locally_moved)
#define AssertThrowMPI(error_code)
virtual void add_periodicity(const std::vector< GridTools::PeriodicFacePair< cell_iterator > > &)
::Triangulation< dim, spacedim >::cell_iterator cell_iterator
static ::ExceptionBase & ExcNotImplemented()
active_cell_iterator last_active() const
Triangulation(MPI_Comm mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing smooth_grid=(::Triangulation< dim, spacedim >::none), const Settings settings=default_setting)
boost::signals2::signal< void()> post_distributed_load
std::vector< unsigned int > compute_point_to_point_communication_pattern(const MPI_Comm &mpi_comm, const std::vector< unsigned int > &destinations)
std::vector< bool > mark_locally_active_vertices_on_level(const int level) const
T max(const T &t, const MPI_Comm &mpi_communicator)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
const std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, std::bitset< 3 > > > & get_periodic_face_map() const
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
static ::ExceptionBase & ExcInternalError()
boost::signals2::signal< void()> pre_distributed_refinement