15#ifndef dealii_grid_tools_h
16#define dealii_grid_tools_h
48#include <boost/archive/binary_iarchive.hpp>
49#include <boost/archive/binary_oarchive.hpp>
50#include <boost/random/mersenne_twister.hpp>
51#include <boost/serialization/array.hpp>
52#include <boost/serialization/vector.hpp>
54#ifdef DEAL_II_WITH_ZLIB
55# include <boost/iostreams/device/back_inserter.hpp>
56# include <boost/iostreams/filter/gzip.hpp>
57# include <boost/iostreams/filtering_stream.hpp>
58# include <boost/iostreams/stream.hpp>
64#ifdef DEAL_II_HAVE_CXX20
77 template <
int dim,
int spacedim>
86 class MappingCollection;
93 template <
int dim,
int spacedim>
100 template <
int dim,
int spacedim,
typename MeshType>
106 using type =
typename MeshType::active_cell_iterator;
113 template <
int dim,
int spacedim>
219 template <
int dim,
typename Transformation,
int spacedim>
222 std::assignable_from<
234 template <
int dim,
int spacedim>
250 template <
int dim,
int spacedim>
289 rotate(const
double angle,
290 const
unsigned int axis,
354 const
Function<dim,
double> *coefficient =
nullptr,
355 const
bool solve_for_absolute_positions = false);
364 template <
int dim,
int spacedim>
366 scale(const
double scaling_factor,
398 template <
int dim,
int spacedim>
403 const
bool keep_boundary = true,
404 const
unsigned int seed =
boost::random::mt19937::default_seed);
495 template <
int dim,
int spacedim>
498 const
double limit_angle_fraction = .75);
563 template <
int dim,
int spacedim>
566 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
567 std::vector<std::vector<Point<dim>>>,
568 std::vector<std::vector<unsigned int>>>
612 template <
int dim,
int spacedim>
615 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
616 std::vector<std::vector<Point<dim>>>,
617 std::vector<std::vector<unsigned int>>,
618 std::vector<unsigned int>>
708 template <
int dim,
int spacedim>
711 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
712 std::vector<std::vector<Point<dim>>>,
713 std::vector<std::vector<unsigned int>>,
714 std::vector<std::vector<Point<spacedim>>>,
715 std::vector<std::vector<unsigned int>>>
723 const double tolerance = 1e-10,
724 const std::vector<bool> &marked_vertices = {},
725 const bool enforce_unique_mapping =
true);
743 template <
int dim,
int spacedim>
768 std::vector<std::tuple<std::pair<int, int>,
798 std::vector<std::tuple<unsigned int, unsigned int, unsigned int>>
823 template <
int dim,
int spacedim>
829 const std::vector<bool> &marked_vertices,
830 const double tolerance,
831 const bool perform_handshake,
832 const bool enforce_unique_mapping =
false);
842 template <
int structdim,
int spacedim>
849 std::array<::Point<spacedim>, structdim + 1>;
859 std::vector<std::tuple<std::pair<int, int>,
873 std::vector<std::tuple<unsigned int, unsigned int, IntersectionType>>
905 const unsigned int n_points_1D,
910 const bool consistent_numbering_of_sender_and_receiver =
false)
const;
920 std::map<unsigned int, std::vector<unsigned int>>
922 const std::vector<std::tuple<unsigned int, unsigned int, unsigned int>>
923 &point_recv_components,
934 template <
int structdim,
int dim,
int spacedim>
938 const std::vector<std::vector<
Point<spacedim>>> &intersection_requests,
940 const std::vector<bool> &marked_vertices,
941 const double tolerance);
954 template <
int spacedim>
985 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
989 const MeshType<dim, spacedim> &mesh,
990 const
Point<spacedim> &p,
991 const
std::vector<
bool> &marked_vertices = {});
1019 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1023 const
Mapping<dim, spacedim> &mapping,
1024 const MeshType<dim, spacedim> &mesh,
1025 const
Point<spacedim> &p,
1026 const
std::vector<
bool> &marked_vertices = {});
1051 template <
typename MeshType>
1054 std::vector<typename MeshType::active_cell_iterator>
1057 typename ::internal::ActiveCellIterator<MeshType::dimension,
1058 MeshType::space_dimension,
1062 const unsigned int vertex_index);
1129 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1133 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim>>
1135 std::pair<typename ::internal::
1136 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1140 const MeshType<dim, spacedim> &mesh,
1142 const std::vector<bool> &marked_vertices = {},
1143 const double tolerance = 1.e-10);
1155 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1159 typename MeshType<dim, spacedim>::active_cell_iterator
1161 typename ::internal::
1162 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type
1166 const std::vector<bool> &marked_vertices = {},
1167 const double tolerance = 1.e-10);
1175 template <
int dim,
int spacedim>
1176 std::pair<typename DoFHandler<dim, spacedim>::active_cell_iterator,
1182 const double tolerance = 1.e-10);
1235 template <
int dim,
int spacedim>
1236 std::pair<typename Triangulation<dim, spacedim>::active_cell_iterator,
1239 const Cache<dim, spacedim> &cache,
1243 const std::vector<bool> &marked_vertices = {},
1244 const double tolerance = 1.e-10);
1262 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1266 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim>>
1268 std::pair<typename ::internal::
1269 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1274 const MeshType<dim, spacedim> &mesh,
1277 std::set<
typename MeshType<dim, spacedim>::active_cell_iterator>>
1280 &vertex_to_cell_centers,
1281 const typename MeshType<dim, spacedim>::active_cell_iterator &cell_hint =
1282 typename MeshType<dim, spacedim>::active_cell_iterator(),
1283 const std::vector<bool> &marked_vertices = {},
1286 const double tolerance = 1.e-10,
1288 std::pair<BoundingBox<spacedim>,
1290 *relevant_cell_bounding_boxes_rtree =
nullptr);
1321 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1325 std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
1328 std::vector<std::pair<
1329 typename ::internal::
1330 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1335 const MeshType<dim, spacedim> &mesh,
1337 const double tolerance,
1338 const std::pair<
typename MeshType<dim, spacedim>::active_cell_iterator,
1341 std::set<
typename MeshType<dim, spacedim>::active_cell_iterator>>
1342 *vertex_to_cells =
nullptr);
1353 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1357 std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
1360 std::vector<std::pair<
1361 typename ::internal::
1362 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1367 const MeshType<dim, spacedim> &mesh,
1369 const double tolerance = 1e-10,
1370 const std::vector<bool> &marked_vertices = {});
1395 template <
typename MeshType>
1398 const typename MeshType::cell_iterator &cell);
1426 template <typename MeshType>
1429 const typename MeshType::active_cell_iterator &cell,
1430 std::vector<typename MeshType::active_cell_iterator> &active_neighbors);
1483 template <typename MeshType>
1487 const MeshType &mesh,
1488 const
std::function<
bool(const typename MeshType::active_cell_iterator &)>
1501 template <typename MeshType>
1505 const MeshType &mesh,
1506 const
std::function<
bool(const typename MeshType::cell_iterator &)>
1508 const
unsigned int level);
1525 template <typename MeshType>
1581 template <typename MeshType>
1585 const MeshType &mesh,
1586 const
std::function<
bool(const typename MeshType::active_cell_iterator &)>
1588 const
double layer_thickness);
1614 template <typename MeshType>
1618 const MeshType &mesh,
1619 const
double layer_thickness);
1693 template <typename MeshType>
1697 const MeshType &mesh,
1698 const
std::function<
bool(const typename MeshType::active_cell_iterator &)>
1700 const
unsigned int refinement_level = 0,
1701 const
bool allow_merge = false,
1702 const
unsigned int max_boxes =
numbers::invalid_unsigned_int);
1731 template <
int spacedim>
1733 std::tuple<std::vector<std::vector<unsigned int>>,
1734 std::map<unsigned int, unsigned int>,
1735 std::map<unsigned int, std::vector<unsigned int>>>
1778 template <
int spacedim>
1780 std::tuple<std::map<unsigned int, std::vector<unsigned int>>,
1781 std::map<unsigned int, unsigned int>,
1782 std::map<unsigned int, std::vector<unsigned int>>>
1802 template <
int dim,
int spacedim>
1803 std::vector<std::vector<Tensor<1, spacedim>>>
1818 template <
int dim,
int spacedim>
1843 template <
int dim,
int spacedim>
1844 std::map<unsigned int, types::global_vertex_index>
1874 template <
int dim,
int spacedim>
1891 template <
int dim,
int spacedim>
1894 const std::vector<unsigned int> &cell_weights,
1944 template <
int dim,
int spacedim>
1962 template <
int dim,
int spacedim>
1965 const std::vector<unsigned int> &cell_weights,
1985 template <
int dim,
int spacedim>
1989 const bool group_siblings =
true);
2002 template <
int dim,
int spacedim>
2013 template <
int dim,
int spacedim>
2014 std::vector<types::subdomain_id>
2016 const std::vector<CellId> &cell_ids);
2028 template <
int dim,
int spacedim>
2031 std::vector<types::subdomain_id> &subdomain);
2047 template <
int dim,
int spacedim>
2082 template <
int dim,
int spacedim>
2112 template <
int dim,
int spacedim>
2170 template <
typename MeshType>
2173 const typename MeshType::active_cell_iterator &cell);
2197 template <
class Container>
2198 std::vector<typename Container::cell_iterator>
2200 const std::vector<typename Container::active_cell_iterator> &patch_cells);
2268 template <
class Container>
2271 const std::vector<typename Container::active_cell_iterator> &patch,
2273 &local_triangulation,
2276 Container::space_dimension>::active_cell_iterator,
2277 typename Container::active_cell_iterator> &patch_to_global_tria_map);
2310 template <
int dim,
int spacedim>
2313 std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>>
2331 template <
typename CellIterator>
2395 template <
typename FaceIterator>
2396 std::optional<unsigned char>
2398 const FaceIterator &face1,
2399 const FaceIterator &face2,
2400 const unsigned int direction,
2463 template <
typename MeshType>
2466 const MeshType &mesh,
2469 const unsigned int direction,
2501 template <
typename MeshType>
2504 const MeshType &mesh,
2505 const
types::boundary_id b_id,
2506 const
unsigned int direction,
2509 const ::
Tensor<1, MeshType::space_dimension> &offset =
2510 ::
Tensor<1, MeshType::space_dimension>(),
2539 template <
int dim,
int spacedim>
2542 const
bool reset_boundary_ids = false);
2565 template <
int dim,
int spacedim>
2568 const
std::vector<
types::boundary_id> &src_boundary_ids,
2569 const
std::vector<
types::manifold_id> &dst_manifold_ids,
2571 const
std::vector<
types::boundary_id> &reset_boundary_ids = {});
2602 template <
int dim,
int spacedim>
2605 const bool compute_face_ids =
false);
2631 template <
int dim,
int spacedim>
2636 const std::set<types::manifold_id> &)> &disambiguation_function =
2637 [](
const std::set<types::manifold_id> &manifold_ids) {
2638 if (manifold_ids.size() == 1)
2639 return *manifold_ids.begin();
2643 bool overwrite_only_flat_manifold_ids =
true);
2732 template <
typename DataType,
typename MeshType>
2735 const MeshType &mesh,
2736 const
std::function<
std::optional<DataType>(
2737 const typename MeshType::active_cell_iterator &)> &pack,
2738 const
std::function<
void(const typename MeshType::active_cell_iterator &,
2739 const DataType &)> &unpack,
2740 const
std::function<
bool(const typename MeshType::active_cell_iterator &)>
2742 always_return<typename MeshType::active_cell_iterator,
bool>{
true});
2756 template <
typename DataType,
typename MeshType>
2759 const MeshType &mesh,
2760 const
std::function<
std::optional<DataType>(
2761 const typename MeshType::level_cell_iterator &)> &pack,
2762 const
std::function<
void(const typename MeshType::level_cell_iterator &,
2763 const DataType &)> &unpack,
2764 const
std::function<
bool(const typename MeshType::level_cell_iterator &)> &
2765 cell_filter =
always_return<typename MeshType::level_cell_iterator,
bool>{
2779 template <
int spacedim>
2780 std::vector<std::vector<BoundingBox<spacedim>>>
2817 template <
int spacedim>
2840 template <
int dim,
int spacedim>
2844 std::map<
unsigned int, std::vector<unsigned int>> &coinciding_vertex_groups,
2845 std::map<unsigned int, unsigned int> &vertex_to_coinciding_vertex_group);
2866 template <
int dim,
int spacedim>
2867 std::map<unsigned int, std::set<::types::subdomain_id>>
2891 template <
int dim,
typename VectorType>
2920 const VectorType &ls_vector,
2921 const double iso_level,
2931 const VectorType &ls_vector,
2932 const double iso_level,
2946 const VectorType &ls_vector,
2947 const double iso_level,
2957 const VectorType &ls_vector,
2958 const double iso_level,
2975 const double iso_level,
2978 const bool write_back_cell_data =
true)
const;
2986 const std::vector<unsigned int> &,
3003 const std::vector<
Point<2>> &points,
3004 const std::vector<unsigned int> &mask,
3005 const double iso_level,
3008 const bool write_back_cell_data)
const;
3015 const std::vector<
Point<3>> &points,
3016 const std::vector<unsigned int> &mask,
3017 const double iso_level,
3020 const bool write_back_cell_data)
const;
3053 <<
"The number of partitions you gave is " << arg1
3054 <<
", but must be greater than zero.");
3060 <<
"The subdomain id " << arg1
3061 <<
" has no cells associated with it.");
3072 <<
"The scaling factor must be positive, but it is " << arg1
3080 <<
"The given vertex with index " << arg1
3081 <<
" is not used in the given triangulation.");
3120 template <
int dim,
typename Transformation,
int spacedim>
3123 std::assignable_from<
3126 void
transform(const Transformation &transformation,
3140 for (; cell != endc; ++cell)
3142 if (treated_vertices[cell->vertex_index(v)] == false)
3145 cell->vertex(v) = transformation(cell->vertex(v));
3147 treated_vertices[cell->vertex_index(v)] =
true;
3157 for (; cell != endc; ++cell)
3158 for (
const unsigned int face : cell->face_indices())
3159 if (cell->face(face)->has_children() &&
3160 !cell->face(face)->at_boundary())
3162 Assert(cell->reference_cell() ==
3167 cell->face(face)->child(0)->vertex(1) =
3168 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1)) /
3177 for (; cell != endc; ++cell)
3178 for (
const unsigned int face : cell->face_indices())
3179 if (cell->face(face)->has_children() &&
3180 !cell->face(face)->at_boundary())
3182 Assert(cell->reference_cell() ==
3187 cell->face(face)->child(0)->vertex(1) =
3188 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1)) /
3190 cell->face(face)->child(0)->vertex(2) =
3191 (cell->face(face)->vertex(0) + cell->face(face)->vertex(2)) /
3193 cell->face(face)->child(1)->vertex(3) =
3194 (cell->face(face)->vertex(1) + cell->face(face)->vertex(3)) /
3196 cell->face(face)->child(2)->vertex(3) =
3197 (cell->face(face)->vertex(2) + cell->face(face)->vertex(3)) /
3201 cell->face(face)->child(0)->vertex(3) =
3202 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1) +
3203 cell->face(face)->vertex(2) + cell->face(face)->vertex(3)) /
3214 template <
typename MeshType>
3217 const typename MeshType::cell_iterator &cell)
3219 std::vector<typename MeshType::active_cell_iterator> child_cells;
3221 if (cell->has_children())
3223 for (
unsigned int child = 0; child < cell->n_children(); ++child)
3224 if (cell->child(child)->has_children())
3226 const std::vector<typename MeshType::active_cell_iterator>
3228 child_cells.insert(child_cells.end(),
3233 child_cells.push_back(cell->child(child));
3241 template <
typename MeshType>
3244 const typename MeshType::active_cell_iterator &cell,
3245 std::vector<typename MeshType::active_cell_iterator> &active_neighbors)
3247 active_neighbors.clear();
3248 for (
const unsigned int n : cell->face_indices())
3249 if (!cell->at_boundary(n))
3251 if (MeshType::dimension == 1)
3258 typename MeshType::cell_iterator neighbor_child =
3260 if (!neighbor_child->is_active())
3262 while (neighbor_child->has_children())
3263 neighbor_child = neighbor_child->child(n == 0 ? 1 : 0);
3265 Assert(neighbor_child->neighbor(n == 0 ? 1 : 0) == cell,
3268 active_neighbors.push_back(neighbor_child);
3272 if (cell->face(n)->has_children())
3276 for (
unsigned int c = 0;
3277 c < cell->face(n)->n_active_descendants();
3279 active_neighbors.push_back(
3280 cell->neighbor_child_on_subface(n, c));
3286 active_neighbors.push_back(cell->neighbor(n));
3294 template <
typename CellIterator>
3298 return sizeof(*this) +
matrix.memory_consumption();
3305 template <
typename DataType,
3307 typename MeshCellIteratorType>
3309 inline void exchange_cell_data(
3310 const MeshType &mesh,
3311 const std::function<std::optional<DataType>(
const MeshCellIteratorType &)>
3313 const std::function<
void(
const MeshCellIteratorType &,
const DataType &)>
3315 const std::function<
bool(
const MeshCellIteratorType &)> &cell_filter,
3316 const std::function<
void(
3317 const std::function<
void(
const MeshCellIteratorType &,
3319 const std::function<std::set<types::subdomain_id>(
3321 MeshType::space_dimension> &)>
3322 &compute_ghost_owners)
3324# ifndef DEAL_II_WITH_MPI
3329 (void)process_cells;
3330 (void)compute_ghost_owners;
3333 constexpr int dim = MeshType::dimension;
3334 constexpr int spacedim = MeshType::space_dimension;
3337 &mesh.get_triangulation());
3341 "The function exchange_cell_data_to_ghosts() only works with parallel triangulations."));
3343 if (
const auto tria =
dynamic_cast<
3345 &mesh.get_triangulation()))
3348 tria->with_artificial_cells(),
3350 "The functions GridTools::exchange_cell_data_to_ghosts() and "
3351 "GridTools::exchange_cell_data_to_level_ghosts() can only "
3352 "operate on a single layer of ghost cells. However, you have "
3353 "given a Triangulation object of type "
3354 "parallel::shared::Triangulation without artificial cells "
3355 "resulting in an arbitrary number of ghost layers. "
3356 "To use this function for a Triangulation object of type "
3357 "parallel::shared::Triangulation, make sure to create the "
3358 "Triangulation object with allow_artificial_cells set to true. "
3359 "This results in a parallel::shared::Triangulation with only "
3360 "a single layer of ghost cells."));
3364 std::set<::types::subdomain_id> ghost_owners =
3365 compute_ghost_owners(*tria);
3367 std::vector<typename CellId::binary_type>>
3370 for (
const auto ghost_owner : ghost_owners)
3371 neighbor_cell_list[ghost_owner] = {};
3373 process_cells([&](
const auto &cell,
const auto key) ->
void {
3374 if (cell_filter(cell))
3376 constexpr int spacedim = MeshType::space_dimension;
3377 neighbor_cell_list[key].emplace_back(
3378 cell->id().template to_binary<spacedim>());
3382 Assert(ghost_owners.size() == neighbor_cell_list.size(),
3394 const int mpi_tag_reply =
3398 std::vector<MPI_Request> requests(ghost_owners.size());
3400 unsigned int idx = 0;
3401 for (
const auto &it : neighbor_cell_list)
3404 const int ierr = MPI_Isend(it.second.data(),
3405 it.second.size() *
sizeof(it.second[0]),
3417 std::vector<MPI_Request> reply_requests(ghost_owners.size());
3418 std::vector<std::vector<char>> sendbuffers(ghost_owners.size());
3420 for (
unsigned int idx = 0; idx < ghost_owners.size(); ++idx)
3423 int ierr = MPI_Probe(MPI_ANY_SOURCE,
3430 ierr = MPI_Get_count(&status, MPI_BYTE, &len);
3437 std::vector<typename CellId::binary_type> cells_with_requests(
3439 std::vector<DataType> data_to_send;
3440 data_to_send.reserve(n_cells);
3441 std::vector<bool> cell_carries_data(n_cells,
false);
3443 ierr = MPI_Recv(cells_with_requests.data(),
3453 for (
unsigned int c = 0; c < static_cast<unsigned int>(n_cells); ++c)
3458 MeshCellIteratorType mesh_it(tria,
3463 const std::optional<DataType> data =
pack(mesh_it);
3466 data_to_send.emplace_back(*data);
3467 cell_carries_data[c] =
true;
3475 sendbuffers[idx].resize(
sizeof(std::size_t));
3480 std::size_t size_of_send =
3484 std::memcpy(sendbuffers[idx].data(),
3486 sizeof(std::size_t));
3490 if (data_to_send.size() < n_cells)
3496 ierr = MPI_Isend(sendbuffers[idx].data(),
3497 sendbuffers[idx].size(),
3502 &reply_requests[idx]);
3507 std::vector<char> receive;
3508 for (
unsigned int id = 0;
id < neighbor_cell_list.size(); ++id)
3511 int ierr = MPI_Probe(MPI_ANY_SOURCE,
3518 ierr = MPI_Get_count(&status, MPI_BYTE, &len);
3521 receive.resize(len);
3523 ierr = MPI_Recv(receive.data(),
3534 auto data_iterator = receive.begin();
3535 std::size_t size_of_received_data =
3537 data_iterator +
sizeof(std::size_t));
3538 data_iterator +=
sizeof(std::size_t);
3543 data_iterator + size_of_received_data,
3545 data_iterator += size_of_received_data;
3550 const std::vector<typename CellId::binary_type> &this_cell_list =
3551 neighbor_cell_list[status.MPI_SOURCE];
3553 std::vector<bool> cells_with_data;
3554 if (received_data.size() < this_cell_list.size())
3557 data_iterator, receive.end(),
false);
3563 auto received_data_iterator = received_data.begin();
3564 for (
unsigned int c = 0; c < this_cell_list.size(); ++c)
3565 if (cells_with_data.empty() || cells_with_data[c])
3570 MeshCellIteratorType cell(tria,
3575 unpack(cell, *received_data_iterator);
3576 ++received_data_iterator;
3582 if (requests.size() > 0)
3585 MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
3588 if (reply_requests.size() > 0)
3590 const int ierr = MPI_Waitall(reply_requests.size(),
3591 reply_requests.data(),
3592 MPI_STATUSES_IGNORE);
3602 template <
typename DataType,
typename MeshType>
3605 const MeshType &mesh,
3606 const std::function<std::optional<DataType>(
3607 const typename MeshType::active_cell_iterator &)> &pack,
3608 const std::function<
void(
const typename MeshType::active_cell_iterator &,
3609 const DataType &)> &unpack,
3610 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
3613# ifndef DEAL_II_WITH_MPI
3620 internal::exchange_cell_data<DataType,
3622 typename MeshType::active_cell_iterator>(
3627 [&](
const auto &process) {
3628 for (
const auto &cell : mesh.active_cell_iterators())
3629 if (cell->is_ghost())
3632 [](
const auto &tria) {
return tria.ghost_owners(); });
3638 template <
typename DataType,
typename MeshType>
3641 const MeshType &mesh,
3642 const std::function<std::optional<DataType>(
3643 const typename MeshType::level_cell_iterator &)> &pack,
3644 const std::function<
void(
const typename MeshType::level_cell_iterator &,
3645 const DataType &)> &unpack,
3646 const std::function<
bool(
const typename MeshType::level_cell_iterator &)>
3649# ifndef DEAL_II_WITH_MPI
3656 internal::exchange_cell_data<DataType,
3658 typename MeshType::level_cell_iterator>(
3663 [&](
const auto &process) {
3664 for (
const auto &cell : mesh.cell_iterators())
3665 if (cell->is_ghost_on_level())
3666 process(cell, cell->level_subdomain_id());
3668 [](
const auto &tria) {
return tria.level_ghost_owners(); });
std::array< unsigned int, 4 > binary_type
Abstract base class for mapping classes.
const Mapping< dim, spacedim > & get_default_linear_mapping() const
virtual MPI_Comm get_communicator() const
cell_iterator create_cell_iterator(const CellId &cell_id) const
cell_iterator end() const
unsigned int n_vertices() const
active_cell_iterator begin_active(const unsigned int level=0) const
typename MeshType::active_cell_iterator type
#define DEAL_II_DEPRECATED
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
#define DeclException0(Exception0)
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcNeedsMPI()
static ::ExceptionBase & ExcNonExistentSubdomain(int arg1)
#define Assert(cond, exc)
static ::ExceptionBase & ExcVertexNotUsed(unsigned int arg1)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcTriangulationHasBeenRefined()
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcInvalidNumberOfPartitions(int arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcScalingFactorNotPositive(double arg1)
#define AssertThrow(cond, exc)
typename ActiveSelector::active_cell_iterator active_cell_iterator
void copy_boundary_to_manifold_id(Triangulation< dim, spacedim > &tria, const bool reset_boundary_ids=false)
void copy_material_to_manifold_id(Triangulation< dim, spacedim > &tria, const bool compute_face_ids=false)
void map_boundary_to_manifold_ids(const std::vector< types::boundary_id > &src_boundary_ids, const std::vector< types::manifold_id > &dst_manifold_ids, Triangulation< dim, spacedim > &tria, const std::vector< types::boundary_id > &reset_boundary_ids={})
void assign_co_dimensional_manifold_indicators(Triangulation< dim, spacedim > &tria, const std::function< types::manifold_id(const std::set< types::manifold_id > &)> &disambiguation_function=[](const std::set< types::manifold_id > &manifold_ids) { if(manifold_ids.size()==1) return *manifold_ids.begin();else return numbers::flat_manifold_id;}, bool overwrite_only_flat_manifold_ids=true)
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
@ matrix
Contents is actually a matrix.
constexpr const ReferenceCell & get_hypercube()
size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
T unpack(const std::vector< char > &buffer, const bool allow_compression=true)
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
const types::manifold_id flat_manifold_id
unsigned int subdomain_id
unsigned int global_dof_index
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
*braid_SplitCommworld & comm
boost::geometry::index::rtree< LeafType, IndexType, IndexableGetter > RTree