16 #ifndef dealii_grid_tools_h
17 # define dealii_grid_tools_h
45 # include <boost/archive/binary_iarchive.hpp>
46 # include <boost/archive/binary_oarchive.hpp>
47 # include <boost/geometry/index/rtree.hpp>
48 # include <boost/serialization/array.hpp>
49 # include <boost/serialization/vector.hpp>
51 # ifdef DEAL_II_WITH_ZLIB
52 # include <boost/iostreams/device/back_inserter.hpp>
53 # include <boost/iostreams/filter/gzip.hpp>
54 # include <boost/iostreams/filtering_stream.hpp>
55 # include <boost/iostreams/stream.hpp>
79 class MappingCollection;
87 template <
int dim,
int spacedim,
class MeshType>
92 using type =
typename MeshType::active_cell_iterator;
99 template <
int dim,
int spacedim>
107 template <
int dim,
int spacedim>
108 class ActiveCellIterator<dim, spacedim, ::
hp::DoFHandler<dim, spacedim>>
127 template <
int dim,
int spacedim>
141 template <
int dim,
int spacedim>
171 template <
int dim,
int spacedim>
187 template <
int dim,
int spacedim>
203 template <
int dim,
int spacedim>
229 template <
int dim,
typename T>
295 template <
int dim,
int spacedim>
318 template <
typename Iterator>
336 template <
int dim,
int spacedim>
338 tuple<std::vector<Point<spacedim>>, std::vector<CellData<dim>>,
SubCellData>
363 template <
int dim,
int spacedim>
386 template <
int dim,
int spacedim>
391 std::vector<unsigned int> & considered_vertices,
392 const double tol = 1
e-12);
453 template <
int dim,
typename Transformation,
int spacedim>
455 transform(
const Transformation & transformation,
463 template <
int dim,
int spacedim>
497 const unsigned int axis,
562 const bool solve_for_absolute_positions =
false);
569 template <
int dim,
int spacedim>
570 std::map<unsigned int, Point<spacedim>>
580 template <
int dim,
int spacedim>
582 scale(
const double scaling_factor,
595 template <
int dim,
int spacedim>
599 const bool keep_boundary =
true);
636 template <
int dim,
int spacedim>
639 const bool isotropic =
false,
640 const unsigned int max_iterations = 100);
668 template <
int dim,
int spacedim>
671 const double max_ratio = 1.6180339887,
672 const unsigned int max_iterations = 5);
765 template <
int dim,
int spacedim>
768 const double limit_angle_fraction = .75);
827 template <
int dim,
int spacedim>
830 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
831 std::vector<std::vector<Point<dim>>>,
832 std::vector<std::vector<unsigned int>>>
871 template <
int dim,
int spacedim>
874 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
875 std::vector<std::vector<Point<dim>>>,
876 std::vector<std::vector<unsigned int>>,
877 std::vector<unsigned int>>
952 template <
int dim,
int spacedim>
955 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
956 std::vector<std::vector<Point<dim>>>,
957 std::vector<std::vector<unsigned int>>,
958 std::vector<std::vector<Point<spacedim>>>,
959 std::vector<std::vector<unsigned int>>>
1006 template <
int dim,
int spacedim>
1007 std::map<unsigned int, Point<spacedim>>
1023 template <
int spacedim>
1054 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1058 const std::vector<bool> & marked_vertices = {});
1085 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1088 const MeshType<dim, spacedim> &mesh,
1090 const std::vector<bool> & marked_vertices = {});
1117 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1119 std::vector<typename MeshType<dim, spacedim>::active_cell_iterator>
1122 typename ::internal::
1123 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type>
1126 const unsigned int vertex_index);
1153 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1155 typename MeshType<dim, spacedim>::active_cell_iterator
1157 typename ::internal::ActiveCellIterator<dim,
1159 MeshType<dim, spacedim>>::type
1163 const std::vector<bool> &marked_vertices = {});
1252 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1254 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim>>
1256 std::pair<typename ::internal::
1257 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1261 const MeshType<dim, spacedim> &mesh,
1263 const std::vector<bool> &marked_vertices = {});
1276 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1278 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim>>
1280 std::pair<typename ::internal::
1281 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1286 const MeshType<dim, spacedim> &mesh,
1289 std::set<
typename MeshType<dim, spacedim>::active_cell_iterator>>
1292 const typename MeshType<dim, spacedim>::active_cell_iterator &cell_hint =
1293 typename MeshType<dim, spacedim>::active_cell_iterator(),
1294 const std::vector<bool> & marked_vertices = {},
1319 template <
int dim,
int spacedim>
1320 std::pair<typename hp::DoFHandler<dim, spacedim>::active_cell_iterator,
1333 template <
int dim,
int spacedim>
1334 std::pair<typename Triangulation<dim, spacedim>::active_cell_iterator,
1337 const Cache<dim, spacedim> &cache,
1341 const std::vector<bool> &marked_vertices = {});
1360 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1362 std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
1365 std::vector<std::pair<
1366 typename ::internal::
1367 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1372 const MeshType<dim, spacedim> &mesh,
1374 const double tolerance = 1
e-10,
1375 const std::vector<bool> & marked_vertices = {});
1398 template <
class MeshType>
1399 std::vector<typename MeshType::active_cell_iterator>
1426 template <
class MeshType>
1429 const typename MeshType::active_cell_iterator & cell,
1430 std::vector<typename MeshType::active_cell_iterator> &active_neighbors);
1484 template <
class MeshType>
1485 std::vector<typename MeshType::active_cell_iterator>
1487 const MeshType &mesh,
1488 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1499 template <
class MeshType>
1500 std::vector<typename MeshType::cell_iterator>
1502 const MeshType &mesh,
1503 const std::function<
bool(
const typename MeshType::cell_iterator &)>
1505 const unsigned int level);
1523 template <
class MeshType>
1524 std::vector<typename MeshType::active_cell_iterator>
1578 template <
class MeshType>
1579 std::vector<typename MeshType::active_cell_iterator>
1581 const MeshType &mesh,
1582 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1584 const double layer_thickness);
1611 template <
class MeshType>
1612 std::vector<typename MeshType::active_cell_iterator>
1614 const double layer_thickness);
1631 template <
class MeshType>
1634 const MeshType &mesh,
1635 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1690 template <
class MeshType>
1691 std::vector<BoundingBox<MeshType::space_dimension>>
1693 const MeshType &mesh,
1694 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1696 const unsigned int refinement_level = 0,
1697 const bool allow_merge =
false,
1729 template <
int spacedim>
1731 std::tuple<std::vector<std::vector<unsigned int>>,
1732 std::map<unsigned int, unsigned int>,
1733 std::map<unsigned int, std::vector<unsigned int>>>
1776 template <
int spacedim>
1778 std::tuple<std::map<unsigned int, std::vector<unsigned int>>,
1779 std::map<unsigned int, unsigned int>,
1780 std::map<unsigned int, std::vector<unsigned int>>>
1797 template <
int dim,
int spacedim>
1799 std::set<typename Triangulation<dim, spacedim>::active_cell_iterator>>
1816 template <
int dim,
int spacedim>
1817 std::vector<std::vector<Tensor<1, spacedim>>>
1834 template <
int dim,
int spacedim>
1853 template <
int dim,
int spacedim>
1854 std::map<unsigned int, types::global_vertex_index>
1871 template <
int dim,
int spacedim>
1872 std::pair<unsigned int, double>
1890 template <
int dim,
int spacedim>
1904 template <
int dim,
int spacedim>
1918 template <
int dim,
int spacedim>
1922 const unsigned int level,
1945 template <
int dim,
int spacedim>
1962 template <
int dim,
int spacedim>
1965 const std::vector<unsigned int> &cell_weights,
2015 template <
int dim,
int spacedim>
2033 template <
int dim,
int spacedim>
2036 const std::vector<unsigned int> &cell_weights,
2056 template <
int dim,
int spacedim>
2060 const bool group_siblings =
true);
2073 template <
int dim,
int spacedim>
2087 template <
int dim,
int spacedim>
2090 std::vector<types::subdomain_id> & subdomain);
2106 template <
int dim,
int spacedim>
2142 template <
int dim,
int spacedim>
2185 template <
typename MeshType>
2186 std::list<std::pair<
typename MeshType::cell_iterator,
2187 typename MeshType::cell_iterator>>
2199 template <
int dim,
int spacedim>
2213 template <
typename MeshType>
2238 template <
int dim,
int spacedim>
2296 template <
class MeshType>
2297 std::vector<typename MeshType::active_cell_iterator>
2324 template <
class Container>
2325 std::vector<typename Container::cell_iterator>
2327 const std::vector<typename Container::active_cell_iterator> &patch_cells);
2397 template <
class Container>
2400 const std::vector<typename Container::active_cell_iterator> &patch,
2402 &local_triangulation,
2405 Container::space_dimension>::active_cell_iterator,
2406 typename Container::active_cell_iterator> &patch_to_global_tria_map);
2444 template <
class DoFHandlerType>
2446 std::vector<typename DoFHandlerType::active_cell_iterator>>
2462 template <
typename CellIterator>
2564 template <
typename FaceIterator>
2566 std::bitset<3> & orientation,
2567 const FaceIterator & face1,
2568 const FaceIterator & face2,
2569 const int direction,
2578 template <
typename FaceIterator>
2581 const FaceIterator & face1,
2582 const FaceIterator & face2,
2583 const int direction,
2647 template <
typename MeshType>
2650 const MeshType & mesh,
2653 const int direction,
2685 template <
typename MeshType>
2688 const MeshType & mesh,
2690 const int direction,
2693 const ::Tensor<1, MeshType::space_dimension> &offset =
2725 template <
int dim,
int spacedim>
2728 const bool reset_boundary_ids =
false);
2753 template <
int dim,
int spacedim>
2756 const std::vector<types::boundary_id> &src_boundary_ids,
2757 const std::vector<types::manifold_id> &dst_manifold_ids,
2759 const std::vector<types::boundary_id> &reset_boundary_ids = {});
2792 template <
int dim,
int spacedim>
2795 const bool compute_face_ids =
false);
2823 template <
int dim,
int spacedim>
2828 const std::set<types::manifold_id> &)> &disambiguation_function =
2829 [](
const std::set<types::manifold_id> &manifold_ids) {
2830 if (manifold_ids.size() == 1)
2831 return *manifold_ids.begin();
2835 bool overwrite_only_flat_manifold_ids =
true);
2920 template <
typename DataType,
typename MeshType>
2923 const MeshType & mesh,
2924 const std::function<std_cxx17::optional<DataType>(
2925 const typename MeshType::active_cell_iterator &)> &
pack,
2926 const std::function<
void(
const typename MeshType::active_cell_iterator &,
2927 const DataType &)> &
unpack);
2940 template <
int spacedim>
2941 std::vector<std::vector<BoundingBox<spacedim>>>
2980 template <
int spacedim>
3005 template <
int dim,
int spacedim>
3009 std::map<
unsigned int, std::vector<unsigned int>> &coinciding_vertex_groups,
3010 std::map<unsigned int, unsigned int> &vertex_to_coinciding_vertex_group);
3018 template <
int dim,
int spacedim>
3019 std::map<unsigned int, std::set<::types::subdomain_id>>
3035 template <
int dim,
typename T>
3055 template <
class Archive>
3057 save(Archive &ar,
const unsigned int version)
const;
3063 template <
class Archive>
3065 load(Archive &ar,
const unsigned int version);
3072 template <
class Archive>
3074 serialize(Archive &archive,
const unsigned int version);
3078 BOOST_SERIALIZATION_SPLIT_MEMBER()
3092 <<
"The number of partitions you gave is " << arg1
3093 <<
", but must be greater than zero.");
3099 <<
"The subdomain id " << arg1
3100 <<
" has no cells associated with it.");
3111 <<
"The scaling factor must be positive, but it is " << arg1
3119 <<
"The point <" << arg1
3120 <<
"> could not be found inside any of the "
3121 <<
"coarse grid cells.");
3128 <<
"The point <" << arg1
3129 <<
"> could not be found inside any of the "
3130 <<
"subcells of a coarse grid cell.");
3137 <<
"The given vertex with index " << arg1
3138 <<
" is not used in the given triangulation.");
3171 template <
int dim,
typename T>
3176 return std::numeric_limits<double>::quiet_NaN();
3194 template <
int dim,
typename Predicate,
int spacedim>
3199 std::vector<bool> treated_vertices(
triangulation.n_vertices(),
false);
3210 for (; cell != endc; ++cell)
3212 if (treated_vertices[cell->vertex_index(v)] ==
false)
3215 cell->vertex(v) = predicate(cell->vertex(v));
3217 treated_vertices[cell->vertex_index(v)] =
true;
3227 for (; cell != endc; ++cell)
3229 if (cell->face(face)->has_children() &&
3230 !cell->face(face)->at_boundary())
3233 cell->face(face)->child(0)->vertex(1) =
3234 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1)) /
3243 for (; cell != endc; ++cell)
3245 if (cell->face(face)->has_children() &&
3246 !cell->face(face)->at_boundary())
3249 cell->face(face)->child(0)->vertex(1) =
3250 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1)) /
3252 cell->face(face)->child(0)->vertex(2) =
3253 (cell->face(face)->vertex(0) + cell->face(face)->vertex(2)) /
3255 cell->face(face)->child(1)->vertex(3) =
3256 (cell->face(face)->vertex(1) + cell->face(face)->vertex(3)) /
3258 cell->face(face)->child(2)->vertex(3) =
3259 (cell->face(face)->vertex(2) + cell->face(face)->vertex(3)) /
3263 cell->face(face)->child(0)->vertex(3) =
3264 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1) +
3265 cell->face(face)->vertex(2) + cell->face(face)->vertex(3)) /
3276 template <
class MeshType>
3277 std::vector<typename MeshType::active_cell_iterator>
3280 std::vector<typename MeshType::active_cell_iterator> child_cells;
3282 if (cell->has_children())
3284 for (
unsigned int child = 0; child < cell->n_children(); ++child)
3285 if (cell->child(child)->has_children())
3287 const std::vector<typename MeshType::active_cell_iterator>
3288 children = get_active_child_cells<MeshType>(cell->child(child));
3289 child_cells.insert(child_cells.end(),
3294 child_cells.push_back(cell->child(child));
3302 template <
class MeshType>
3305 const typename MeshType::active_cell_iterator & cell,
3306 std::vector<typename MeshType::active_cell_iterator> &active_neighbors)
3308 active_neighbors.clear();
3309 for (
const unsigned int n :
3311 if (!cell->at_boundary(n))
3313 if (MeshType::dimension == 1)
3320 typename MeshType::cell_iterator neighbor_child =
3322 if (!neighbor_child->is_active())
3324 while (neighbor_child->has_children())
3325 neighbor_child = neighbor_child->child(n == 0 ? 1 : 0);
3327 Assert(neighbor_child->neighbor(n == 0 ? 1 : 0) == cell,
3330 active_neighbors.push_back(neighbor_child);
3334 if (cell->face(n)->has_children())
3338 for (
unsigned int c = 0;
3339 c < cell->face(n)->number_of_children();
3341 active_neighbors.push_back(
3342 cell->neighbor_child_on_subface(n, c));
3348 active_neighbors.push_back(cell->neighbor(n));
3358 namespace ProjectToObject
3372 struct CrossDerivative
3374 const unsigned int direction_0;
3375 const unsigned int direction_1;
3377 CrossDerivative(
const unsigned int d0,
const unsigned int d1);
3380 inline CrossDerivative::CrossDerivative(
const unsigned int d0,
3381 const unsigned int d1)
3392 template <
typename F>
3394 centered_first_difference(
const double center,
3398 return (f(
center + step) - f(
center - step)) / (2.0 * step);
3407 template <
typename F>
3409 centered_second_difference(
const double center,
3428 template <
int structdim,
typename F>
3431 const CrossDerivative cross_derivative,
3437 simplex_vector[cross_derivative.direction_0] = 0.5 * step;
3438 simplex_vector[cross_derivative.direction_1] = -0.5 * step;
3439 return (-4.0 * f(
center) - 1.0 * f(
center + simplex_vector) -
3440 1.0 / 3.0 * f(
center - simplex_vector) +
3441 16.0 / 3.0 * f(
center + 0.5 * simplex_vector)) /
3453 template <
int spacedim,
int structdim,
typename F>
3456 const unsigned int row_n,
3457 const unsigned int dependent_direction,
3464 dependent_direction <
3466 ExcMessage(
"This function assumes that the last weight is a "
3467 "dependent variable (and hence we cannot take its "
3468 "derivative directly)."));
3469 Assert(row_n != dependent_direction,
3471 "We cannot differentiate with respect to the variable "
3472 "that is assumed to be dependent."));
3476 {row_n, dependent_direction},
center, step, f);
3478 for (
unsigned int dim_n = 0; dim_n < spacedim; ++dim_n)
3480 -2.0 * (p0[dim_n] - manifold_point[dim_n]) * stencil_value[dim_n];
3489 template <
typename Iterator,
int spacedim,
int structdim>
3491 project_to_d_linear_object(
const Iterator &
object,
3526 for (
unsigned int d = 0;
d < structdim; ++
d)
3531 x_k +=
object->vertex(i) *
3537 for (
const unsigned int i :
3540 (x_k - trial_point) * object->vertex(i) *
3545 for (
const unsigned int i :
3547 for (
const unsigned int j :
3555 H_k += (
object->vertex(i) *
object->vertex(j)) * tmp;
3562 for (
const unsigned int i :
3564 x_k +=
object->vertex(i) *
3567 if (delta_xi.
norm() < 1
e-7)
3583 template <
int structdim>
3590 static const std::size_t n_vertices_per_cell =
3592 n_independent_components;
3593 std::array<double, n_vertices_per_cell> copied_weights;
3594 for (
unsigned int i = 0; i < n_vertices_per_cell; ++i)
3596 copied_weights[i] = v[i];
3597 if (v[i] < 0.0 || v[i] > 1.0)
3602 std::sort(copied_weights.begin(), copied_weights.end());
3604 std::accumulate(copied_weights.begin(), copied_weights.end(), 0.0);
3609 template <
typename Iterator>
3615 const int spacedim = Iterator::AccessorType::space_dimension;
3616 const int structdim = Iterator::AccessorType::structure_dimension;
3620 if (structdim >= spacedim)
3621 return projected_point;
3622 else if (structdim == 1 || structdim == 2)
3624 using namespace internal::ProjectToObject;
3629 const int dim = Iterator::AccessorType::dimension;
3632 &manifold) !=
nullptr)
3635 project_to_d_linear_object<Iterator, spacedim, structdim>(
3636 object, trial_point);
3681 const double step_size =
object->diameter() / 64.0;
3683 constexpr
unsigned int n_vertices_per_cell =
3686 std::array<Point<spacedim>, n_vertices_per_cell>
vertices;
3687 for (
unsigned int vertex_n = 0; vertex_n < n_vertices_per_cell;
3689 vertices[vertex_n] = object->vertex(vertex_n);
3691 auto get_point_from_weights =
3694 return object->get_manifold().get_new_point(
3702 double guess_weights_sum = 0.0;
3703 for (
unsigned int vertex_n = 0; vertex_n < n_vertices_per_cell;
3706 const double distance =
3708 if (distance == 0.0)
3710 guess_weights = 0.0;
3711 guess_weights[vertex_n] = 1.0;
3712 guess_weights_sum = 1.0;
3717 guess_weights[vertex_n] = 1.0 / distance;
3718 guess_weights_sum += guess_weights[vertex_n];
3721 guess_weights /= guess_weights_sum;
3722 Assert(internal::weights_are_ok<structdim>(guess_weights),
3731 for (
unsigned int outer_n = 0; outer_n < 40; ++outer_n)
3733 const unsigned int dependent_direction =
3734 n_vertices_per_cell - 1;
3736 for (
unsigned int row_n = 0; row_n < n_vertices_per_cell;
3739 if (row_n != dependent_direction)
3741 current_gradient[row_n] =
3742 gradient_entry<spacedim, structdim>(
3744 dependent_direction,
3748 get_point_from_weights);
3750 current_gradient[dependent_direction] -=
3751 current_gradient[row_n];
3769 double gradient_weight = -0.5;
3770 auto gradient_weight_objective_function =
3771 [&](
const double gradient_weight_guess) ->
double {
3772 return (trial_point -
3773 get_point_from_weights(guess_weights +
3774 gradient_weight_guess *
3779 for (
unsigned int inner_n = 0; inner_n < 10; ++inner_n)
3781 const double update_numerator = centered_first_difference(
3784 gradient_weight_objective_function);
3785 const double update_denominator =
3786 centered_second_difference(
3789 gradient_weight_objective_function);
3793 if (std::abs(update_denominator) == 0.0)
3796 gradient_weight - update_numerator / update_denominator;
3801 if (std::abs(gradient_weight) > 10)
3803 gradient_weight = -10.0;
3813 guess_weights + gradient_weight * current_gradient;
3815 double new_gradient_weight = gradient_weight;
3816 for (
unsigned int iteration_count = 0; iteration_count < 40;
3819 if (internal::weights_are_ok<structdim>(tentative_weights))
3822 for (
unsigned int i = 0; i < n_vertices_per_cell; ++i)
3824 if (tentative_weights[i] < 0.0)
3826 tentative_weights -=
3827 (tentative_weights[i] / current_gradient[i]) *
3830 if (tentative_weights[i] < 0.0 ||
3831 1.0 < tentative_weights[i])
3833 new_gradient_weight /= 2.0;
3836 new_gradient_weight * current_gradient;
3843 if (!internal::weights_are_ok<structdim>(tentative_weights))
3848 if (get_point_from_weights(tentative_weights)
3849 .distance(trial_point) <
3850 get_point_from_weights(guess_weights).distance(trial_point))
3851 guess_weights = tentative_weights;
3854 Assert(internal::weights_are_ok<structdim>(guess_weights),
3857 Assert(internal::weights_are_ok<structdim>(guess_weights),
3859 projected_point = get_point_from_weights(guess_weights);
3869 for (
unsigned int line_n = 0;
3870 line_n < GeometryInfo<structdim>::lines_per_cell;
3873 line_projections[line_n] =
3876 std::sort(line_projections.begin(),
3877 line_projections.end(),
3879 return a.distance(trial_point) <
3880 b.distance(trial_point);
3882 if (line_projections[0].distance(trial_point) <
3883 projected_point.
distance(trial_point))
3884 projected_point = line_projections[0];
3890 return projected_point;
3893 return projected_point;
3898 template <
int dim,
typename T>
3899 template <
class Archive>
3901 CellDataTransferBuffer<dim, T>::save(Archive &ar,
3902 const unsigned int )
const
3904 Assert(cell_ids.size() == data.size(),
3907 const std::size_t
n_cells = cell_ids.size();
3909 for (
const auto &
id : cell_ids)
3912 ar & binary_cell_id;
3920 template <
int dim,
typename T>
3921 template <
class Archive>
3923 CellDataTransferBuffer<dim, T>::load(Archive &ar,
3924 const unsigned int )
3930 for (
unsigned int c = 0; c <
n_cells; ++c)
3934 cell_ids.emplace_back(
value);
3941 template <
typename DataType,
typename MeshType>
3944 const MeshType & mesh,
3945 const std::function<std_cxx17::optional<DataType>(
3946 const typename MeshType::active_cell_iterator &)> &
pack,
3947 const std::function<
void(
const typename MeshType::active_cell_iterator &,
3948 const DataType &)> &
unpack)
3950 # ifndef DEAL_II_WITH_MPI
3956 "GridTools::exchange_cell_data_to_ghosts() requires MPI."));
3958 constexpr
int dim = MeshType::dimension;
3959 constexpr
int spacedim = MeshType::space_dimension;
3962 &mesh.get_triangulation());
3966 "The function exchange_cell_data_to_ghosts() only works with parallel triangulations."));
3969 using DestinationToBufferMap =
3971 CellDataTransferBuffer<dim, DataType>>;
3972 DestinationToBufferMap destination_to_data_buffer_map;
3974 std::map<unsigned int, std::set<::types::subdomain_id>>
3979 if (cell->is_locally_owned())
3981 std::set<::types::subdomain_id> send_to;
3984 const std::map<
unsigned int,
3985 std::set<::types::subdomain_id>>::
3986 const_iterator neighbor_subdomains_of_vertex =
3989 if (neighbor_subdomains_of_vertex ==
3993 Assert(neighbor_subdomains_of_vertex->second.size() != 0,
3996 send_to.insert(neighbor_subdomains_of_vertex->second.begin(),
3997 neighbor_subdomains_of_vertex->second.end());
4000 if (send_to.size() > 0)
4003 typename MeshType::active_cell_iterator mesh_it(tria,
4008 const std_cxx17::optional<DataType> data =
pack(mesh_it);
4012 const CellId cellid = cell->id();
4014 for (
const auto subdomain : send_to)
4018 typename DestinationToBufferMap::iterator p =
4019 destination_to_data_buffer_map
4020 .insert(std::make_pair(
4021 subdomain, CellDataTransferBuffer<dim, DataType>()))
4024 p->second.cell_ids.emplace_back(cellid);
4025 p->second.data.emplace_back(*data);
4034 tria->get_communicator());
4040 std::set<::types::subdomain_id> ghost_owners = tria->ghost_owners();
4041 const unsigned int n_ghost_owners = ghost_owners.size();
4042 std::vector<std::vector<char>> sendbuffers(n_ghost_owners);
4043 std::vector<MPI_Request> requests(n_ghost_owners);
4045 unsigned int idx = 0;
4046 for (
auto it = ghost_owners.begin(); it != ghost_owners.end(); ++it, ++idx)
4048 CellDataTransferBuffer<dim, DataType> &data =
4049 destination_to_data_buffer_map[*it];
4055 const int ierr = MPI_Isend(sendbuffers[idx].data(),
4056 sendbuffers[idx].size(),
4060 tria->get_communicator(),
4066 std::vector<char> receive;
4067 for (
unsigned int idx = 0; idx < n_ghost_owners; ++idx)
4071 MPI_Probe(MPI_ANY_SOURCE, mpi_tag, tria->get_communicator(), &status);
4075 ierr = MPI_Get_count(&status, MPI_BYTE, &len);
4078 receive.resize(len);
4080 char *ptr = receive.data();
4081 ierr = MPI_Recv(ptr,
4086 tria->get_communicator(),
4091 Utilities::unpack<CellDataTransferBuffer<dim, DataType>>(
4094 DataType *data = cellinfo.data.data();
4095 for (
unsigned int c = 0; c < cellinfo.cell_ids.size(); ++c, ++data)
4098 tria_cell = cellinfo.cell_ids[c].to_cell(*tria);
4100 const typename MeshType::active_cell_iterator cell(
4101 tria, tria_cell->level(), tria_cell->index(), &mesh);
4109 if (requests.size())
4112 MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
4115 # endif // DEAL_II_WITH_MPI