16#ifndef dealii_grid_tools_h
17#define dealii_grid_tools_h
53#include <boost/archive/binary_iarchive.hpp>
54#include <boost/archive/binary_oarchive.hpp>
55#include <boost/random/mersenne_twister.hpp>
56#include <boost/serialization/array.hpp>
57#include <boost/serialization/vector.hpp>
59#ifdef DEAL_II_WITH_ZLIB
60# include <boost/iostreams/device/back_inserter.hpp>
61# include <boost/iostreams/filter/gzip.hpp>
62# include <boost/iostreams/filtering_stream.hpp>
63# include <boost/iostreams/stream.hpp>
87 class MappingCollection;
94 template <
int dim,
int spacedim>
101 template <
int dim,
int spacedim,
class MeshType>
106 using type =
typename MeshType::active_cell_iterator;
113 template <
int dim,
int spacedim>
144 template <
int dim,
int spacedim>
174 template <
int dim,
int spacedim>
178 (ReferenceCells::get_hypercube<dim>()
180 .
template get_default_linear_mapping<dim, spacedim>()
182 .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
196 template <
int dim,
int spacedim>
201 (ReferenceCells::get_hypercube<dim>()
203 .
template get_default_linear_mapping<dim, spacedim>()
205 .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
219 template <
int dim,
int spacedim>
224 (ReferenceCells::get_hypercube<dim>()
226 .
template get_default_linear_mapping<dim, spacedim>()
228 .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
295 template <
int dim,
int spacedim>
357 template <
int dim,
int spacedim>
378 template <
typename Iterator>
381 const Iterator &
object,
396 template <
int dim,
int spacedim>
398 tuple<std::vector<Point<spacedim>>, std::vector<CellData<dim>>,
SubCellData>
423 template <
int dim,
int spacedim>
446 template <
int dim,
int spacedim>
451 std::vector<unsigned int> & considered_vertices,
452 const double tol = 1e-12);
472 template <
int dim,
int spacedim>
487 template <
int dim,
int spacedim>
575 template <
int dim,
typename Transformation,
int spacedim>
586 template <
int dim,
int spacedim>
640 DEAL_II_DEPRECATED_EARLY
void
641 rotate(
const double angle,
642 const unsigned int axis,
707 const bool solve_for_absolute_positions =
false);
714 template <
int dim,
int spacedim>
715 std::map<unsigned int, Point<spacedim>>
725 template <
int dim,
int spacedim>
727 scale(
const double scaling_factor,
744 template <
int dim,
int spacedim>
749 const bool keep_boundary =
true,
750 const unsigned int seed = boost::random::mt19937::default_seed);
785 template <
int dim,
int spacedim>
788 const bool isotropic =
false,
789 const unsigned int max_iterations = 100);
815 template <
int dim,
int spacedim>
818 const double max_ratio = 1.6180339887,
819 const unsigned int max_iterations = 5);
910 template <
int dim,
int spacedim>
913 const double limit_angle_fraction = .75);
975 template <
int dim,
int spacedim>
978 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
979 std::vector<std::vector<Point<dim>>>,
980 std::vector<std::vector<unsigned int>>>
1024 template <
int dim,
int spacedim>
1027 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
1028 std::vector<std::vector<Point<dim>>>,
1029 std::vector<std::vector<unsigned int>>,
1030 std::vector<unsigned int>>
1110 template <
int dim,
int spacedim>
1113 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
1114 std::vector<std::vector<Point<dim>>>,
1115 std::vector<std::vector<unsigned int>>,
1116 std::vector<std::vector<Point<spacedim>>>,
1117 std::vector<std::vector<unsigned int>>>
1125 const double tolerance = 1e-10);
1143 template <
int dim,
int spacedim>
1152 std::vector<std::tuple<std::pair<int, int>,
1182 std::vector<std::tuple<unsigned int, unsigned int, unsigned int>>
1207 template <
int dim,
int spacedim>
1213 const std::vector<bool> & marked_vertices,
1214 const double tolerance,
1215 const bool perform_handshake,
1216 const bool enforce_unique_mapping =
false);
1256 template <
int dim,
int spacedim>
1257 std::map<unsigned int, Point<spacedim>>
1261 (ReferenceCells::get_hypercube<dim>()
1263 .
template get_default_linear_mapping<dim, spacedim>()
1265 .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
1278 template <
int spacedim>
1306 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1310 const std::vector<bool> & marked_vertices = {});
1335 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1338 const MeshType<dim, spacedim> &mesh,
1340 const std::vector<bool> & marked_vertices = {});
1362 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1364 std::vector<typename MeshType<dim, spacedim>::active_cell_iterator>
1367 typename ::internal::
1368 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type>
1371 const unsigned int vertex_index);
1435 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1437 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim>>
1439 std::pair<typename ::internal::
1440 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1444 const MeshType<dim, spacedim> &mesh,
1446 const std::vector<bool> &marked_vertices = {},
1447 const double tolerance = 1.e-10);
1456 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1458 typename MeshType<dim, spacedim>::active_cell_iterator
1460 typename ::internal::
1461 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type
1465 const std::vector<bool> &marked_vertices = {},
1466 const double tolerance = 1.e-10);
1474 template <
int dim,
int spacedim>
1475 std::pair<typename DoFHandler<dim, spacedim>::active_cell_iterator,
1481 const double tolerance = 1.e-10);
1534 template <
int dim,
int spacedim>
1535 std::pair<typename Triangulation<dim, spacedim>::active_cell_iterator,
1538 const Cache<dim, spacedim> &cache,
1542 const std::vector<bool> &marked_vertices = {},
1543 const double tolerance = 1.e-10);
1558 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1560 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim>>
1562 std::pair<typename ::internal::
1563 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1568 const MeshType<dim, spacedim> &mesh,
1571 std::set<
typename MeshType<dim, spacedim>::active_cell_iterator>>
1574 const typename MeshType<dim, spacedim>::active_cell_iterator &cell_hint =
1575 typename MeshType<dim, spacedim>::active_cell_iterator(),
1576 const std::vector<bool> & marked_vertices = {},
1579 const double tolerance = 1.e-10,
1581 std::pair<BoundingBox<spacedim>,
1583 *relevant_cell_bounding_boxes_rtree =
nullptr);
1606 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1608 std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
1611 std::vector<std::pair<
1612 typename ::internal::
1613 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1618 const MeshType<dim, spacedim> &mesh,
1620 const double tolerance,
1621 const std::pair<
typename MeshType<dim, spacedim>::active_cell_iterator,
1630 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1632 std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
1635 std::vector<std::pair<
1636 typename ::internal::
1637 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1642 const MeshType<dim, spacedim> &mesh,
1644 const double tolerance = 1e-10,
1645 const std::vector<bool> & marked_vertices = {});
1668 template <
class MeshType>
1669 std::vector<typename MeshType::active_cell_iterator>
1696 template <
class MeshType>
1699 const typename MeshType::active_cell_iterator & cell,
1700 std::vector<typename MeshType::active_cell_iterator> &active_neighbors);
1751 template <
class MeshType>
1752 std::vector<typename MeshType::active_cell_iterator>
1754 const MeshType &mesh,
1755 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1766 template <
class MeshType>
1767 std::vector<typename MeshType::cell_iterator>
1769 const MeshType &mesh,
1770 const std::function<
bool(
const typename MeshType::cell_iterator &)>
1772 const unsigned int level);
1787 template <
class MeshType>
1788 std::vector<typename MeshType::active_cell_iterator>
1839 template <
class MeshType>
1840 std::vector<typename MeshType::active_cell_iterator>
1842 const MeshType &mesh,
1843 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1845 const double layer_thickness);
1869 template <
class MeshType>
1870 std::vector<typename MeshType::active_cell_iterator>
1872 const double layer_thickness);
1889 template <
class MeshType>
1892 const MeshType &mesh,
1893 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1948 template <
class MeshType>
1949 std::vector<BoundingBox<MeshType::space_dimension>>
1951 const MeshType &mesh,
1952 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1954 const unsigned int refinement_level = 0,
1955 const bool allow_merge =
false,
1985 template <
int spacedim>
1987 std::tuple<std::vector<std::vector<unsigned int>>,
1988 std::map<unsigned int, unsigned int>,
1989 std::map<unsigned int, std::vector<unsigned int>>>
2032 template <
int spacedim>
2034 std::tuple<std::map<unsigned int, std::vector<unsigned int>>,
2035 std::map<unsigned int, unsigned int>,
2036 std::map<unsigned int, std::vector<unsigned int>>>
2053 template <
int dim,
int spacedim>
2055 std::set<typename Triangulation<dim, spacedim>::active_cell_iterator>>
2070 template <
int dim,
int spacedim>
2071 std::vector<std::vector<Tensor<1, spacedim>>>
2086 template <
int dim,
int spacedim>
2092 (ReferenceCells::get_hypercube<dim>()
2094 .
template get_default_linear_mapping<dim, spacedim>()
2096 .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
2111 template <
int dim,
int spacedim>
2112 std::map<unsigned int, types::global_vertex_index>
2127 template <
int dim,
int spacedim>
2128 std::pair<unsigned int, double>
2146 template <
int dim,
int spacedim>
2160 template <
int dim,
int spacedim>
2174 template <
int dim,
int spacedim>
2178 const unsigned int level,
2201 template <
int dim,
int spacedim>
2218 template <
int dim,
int spacedim>
2221 const std::vector<unsigned int> &cell_weights,
2271 template <
int dim,
int spacedim>
2289 template <
int dim,
int spacedim>
2292 const std::vector<unsigned int> &cell_weights,
2312 template <
int dim,
int spacedim>
2316 const bool group_siblings =
true);
2329 template <
int dim,
int spacedim>
2340 template <
int dim,
int spacedim>
2341 std::vector<types::subdomain_id>
2343 const std::vector<CellId> & cell_ids);
2355 template <
int dim,
int spacedim>
2358 std::vector<types::subdomain_id> & subdomain);
2374 template <
int dim,
int spacedim>
2409 template <
int dim,
int spacedim>
2452 template <
typename MeshType>
2453 std::list<std::pair<
typename MeshType::cell_iterator,
2454 typename MeshType::cell_iterator>>
2466 template <
int dim,
int spacedim>
2480 template <
typename MeshType>
2505 template <
int dim,
int spacedim>
2561 template <
class MeshType>
2562 std::vector<typename MeshType::active_cell_iterator>
2587 template <
class Container>
2588 std::vector<typename Container::cell_iterator>
2590 const std::vector<typename Container::active_cell_iterator> &patch_cells);
2658 template <
class Container>
2661 const std::vector<typename Container::active_cell_iterator> &patch,
2663 &local_triangulation,
2666 Container::space_dimension>::active_cell_iterator,
2667 typename Container::active_cell_iterator> &patch_to_global_tria_map);
2700 template <
int dim,
int spacedim>
2703 std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>>
2719 template <
typename CellIterator>
2819 template <
typename FaceIterator>
2822 std::bitset<3> & orientation,
2823 const FaceIterator & face1,
2824 const FaceIterator & face2,
2825 const unsigned int direction,
2834 template <
typename FaceIterator>
2837 const FaceIterator & face1,
2838 const FaceIterator & face2,
2839 const unsigned int direction,
2901 template <
typename MeshType>
2904 const MeshType & mesh,
2907 const unsigned int direction,
2937 template <
typename MeshType>
2940 const MeshType & mesh,
2942 const unsigned int direction,
2945 const ::Tensor<1, MeshType::space_dimension> &offset =
2975 template <
int dim,
int spacedim>
2978 const bool reset_boundary_ids =
false);
3001 template <
int dim,
int spacedim>
3004 const std::vector<types::boundary_id> &src_boundary_ids,
3005 const std::vector<types::manifold_id> &dst_manifold_ids,
3007 const std::vector<types::boundary_id> &reset_boundary_ids = {});
3038 template <
int dim,
int spacedim>
3041 const bool compute_face_ids =
false);
3067 template <
int dim,
int spacedim>
3072 const std::set<types::manifold_id> &)> &disambiguation_function =
3073 [](
const std::set<types::manifold_id> &manifold_ids) {
3074 if (manifold_ids.size() == 1)
3075 return *manifold_ids.begin();
3079 bool overwrite_only_flat_manifold_ids =
true);
3166 template <
typename DataType,
typename MeshType>
3169 const MeshType & mesh,
3170 const std::function<std_cxx17::optional<DataType>(
3171 const typename MeshType::active_cell_iterator &)> &pack,
3172 const std::function<
void(
const typename MeshType::active_cell_iterator &,
3173 const DataType &)> & unpack,
3174 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
3188 template <
typename DataType,
typename MeshType>
3191 const MeshType & mesh,
3192 const std::function<std_cxx17::optional<DataType>(
3193 const typename MeshType::level_cell_iterator &)> &pack,
3194 const std::function<
void(
const typename MeshType::level_cell_iterator &,
3195 const DataType &)> & unpack,
3196 const std::function<
bool(
const typename MeshType::level_cell_iterator &)> &
3211 template <
int spacedim>
3212 std::vector<std::vector<BoundingBox<spacedim>>>
3215 const MPI_Comm & mpi_communicator);
3249 template <
int spacedim>
3253 const MPI_Comm & mpi_communicator);
3272 template <
int dim,
int spacedim>
3285 template <
int dim,
int spacedim>
3286 std::map<unsigned int, std::set<::types::subdomain_id>>
3305 template <
int dim,
typename T>
3326 template <
class Archive>
3328 save(Archive &ar,
const unsigned int)
const
3332 Assert(cell_ids.size() == data.size(),
3335 const std::size_t n_cells = cell_ids.size();
3337 for (
const auto &
id : cell_ids)
3340 ar & binary_cell_id;
3352 template <
class Archive>
3354 load(Archive &ar,
const unsigned int)
3358 std::size_t n_cells;
3361 cell_ids.reserve(n_cells);
3362 for (
unsigned int c = 0; c < n_cells; ++c)
3366 cell_ids.emplace_back(value);
3378 template <
class Archive>
3384 BOOST_SERIALIZATION_SPLIT_MEMBER()
3407 template <
int dim,
typename VectorType>
3430 const VectorType & ls_vector,
3431 const double iso_level,
3443 const VectorType & ls_vector,
3444 const double iso_level,
3462 const double iso_level,
3474 const std::vector<
Point<2>> & points,
3475 const std::vector<unsigned int> mask,
3476 const double iso_level,
3485 const std::vector<
Point<3>> & points,
3486 const std::vector<unsigned int> mask,
3487 const double iso_level,
3522 <<
"The number of partitions you gave is " << arg1
3523 <<
", but must be greater than zero.");
3529 <<
"The subdomain id " << arg1
3530 <<
" has no cells associated with it.");
3541 <<
"The scaling factor must be positive, but it is " << arg1
3549 <<
"The given vertex with index " << arg1
3550 <<
" is not used in the given triangulation.");
3567 "The edges of the mesh are not consistently orientable.");
3604 template <
int dim,
typename Predicate,
int spacedim>
3609 std::vector<bool> treated_vertices(
triangulation.n_vertices(),
false);
3620 for (; cell != endc; ++cell)
3621 for (
const unsigned int v : cell->vertex_indices())
3622 if (treated_vertices[cell->vertex_index(v)] ==
false)
3625 cell->vertex(v) = predicate(cell->vertex(v));
3627 treated_vertices[cell->vertex_index(v)] =
true;
3637 for (; cell != endc; ++cell)
3638 for (
const unsigned int face : cell->face_indices())
3639 if (cell->face(face)->has_children() &&
3640 !cell->face(face)->at_boundary())
3642 Assert(cell->reference_cell() ==
3643 ReferenceCells::get_hypercube<dim>(),
3647 cell->face(face)->child(0)->vertex(1) =
3648 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1)) /
3657 for (; cell != endc; ++cell)
3658 for (
const unsigned int face : cell->face_indices())
3659 if (cell->face(face)->has_children() &&
3660 !cell->face(face)->at_boundary())
3662 Assert(cell->reference_cell() ==
3663 ReferenceCells::get_hypercube<dim>(),
3667 cell->face(face)->child(0)->vertex(1) =
3668 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1)) /
3670 cell->face(face)->child(0)->vertex(2) =
3671 (cell->face(face)->vertex(0) + cell->face(face)->vertex(2)) /
3673 cell->face(face)->child(1)->vertex(3) =
3674 (cell->face(face)->vertex(1) + cell->face(face)->vertex(3)) /
3676 cell->face(face)->child(2)->vertex(3) =
3677 (cell->face(face)->vertex(2) + cell->face(face)->vertex(3)) /
3681 cell->face(face)->child(0)->vertex(3) =
3682 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1) +
3683 cell->face(face)->vertex(2) + cell->face(face)->vertex(3)) /
3694 template <
class MeshType>
3695 std::vector<typename MeshType::active_cell_iterator>
3698 std::vector<typename MeshType::active_cell_iterator> child_cells;
3700 if (cell->has_children())
3702 for (
unsigned int child = 0; child < cell->n_children(); ++child)
3703 if (cell->child(child)->has_children())
3705 const std::vector<typename MeshType::active_cell_iterator>
3706 children = get_active_child_cells<MeshType>(cell->child(child));
3707 child_cells.insert(child_cells.end(),
3712 child_cells.push_back(cell->child(child));
3720 template <
class MeshType>
3723 const typename MeshType::active_cell_iterator & cell,
3724 std::vector<typename MeshType::active_cell_iterator> &active_neighbors)
3726 active_neighbors.clear();
3727 for (
const unsigned int n : cell->face_indices())
3728 if (!cell->at_boundary(n))
3730 if (MeshType::dimension == 1)
3737 typename MeshType::cell_iterator neighbor_child =
3739 if (!neighbor_child->is_active())
3741 while (neighbor_child->has_children())
3742 neighbor_child = neighbor_child->child(n == 0 ? 1 : 0);
3744 Assert(neighbor_child->neighbor(n == 0 ? 1 : 0) == cell,
3747 active_neighbors.push_back(neighbor_child);
3751 if (cell->face(n)->has_children())
3755 for (
unsigned int c = 0;
3756 c < cell->face(n)->n_active_descendants();
3758 active_neighbors.push_back(
3759 cell->neighbor_child_on_subface(n, c));
3765 active_neighbors.push_back(cell->neighbor(n));
3775 namespace ProjectToObject
3789 struct CrossDerivative
3791 const unsigned int direction_0;
3792 const unsigned int direction_1;
3794 CrossDerivative(
const unsigned int d0,
const unsigned int d1);
3797 inline CrossDerivative::CrossDerivative(
const unsigned int d0,
3798 const unsigned int d1)
3809 template <
typename F>
3811 centered_first_difference(
const double center,
3815 return (f(
center + step) - f(
center - step)) / (2.0 * step);
3824 template <
typename F>
3826 centered_second_difference(
const double center,
3845 template <
int structdim,
typename F>
3848 const CrossDerivative cross_derivative,
3854 simplex_vector[cross_derivative.direction_0] = 0.5 * step;
3855 simplex_vector[cross_derivative.direction_1] = -0.5 * step;
3856 return (-4.0 * f(
center) - 1.0 * f(
center + simplex_vector) -
3857 1.0 / 3.0 * f(
center - simplex_vector) +
3858 16.0 / 3.0 * f(
center + 0.5 * simplex_vector)) /
3870 template <
int spacedim,
int structdim,
typename F>
3873 const unsigned int row_n,
3874 const unsigned int dependent_direction,
3881 dependent_direction <
3883 ExcMessage(
"This function assumes that the last weight is a "
3884 "dependent variable (and hence we cannot take its "
3885 "derivative directly)."));
3886 Assert(row_n != dependent_direction,
3888 "We cannot differentiate with respect to the variable "
3889 "that is assumed to be dependent."));
3893 {row_n, dependent_direction},
center, step, f);
3895 for (
unsigned int dim_n = 0; dim_n < spacedim; ++dim_n)
3897 -2.0 * (p0[dim_n] - manifold_point[dim_n]) * stencil_value[dim_n];
3906 template <
typename Iterator,
int spacedim,
int structdim>
3908 project_to_d_linear_object(
const Iterator &
object,
3943 for (
unsigned int d = 0;
d < structdim; ++
d)
3948 x_k +=
object->vertex(i) *
3954 for (
const unsigned int i :
3957 (x_k - trial_point) * object->vertex(i) *
3962 for (
const unsigned int i :
3964 for (
const unsigned int j :
3972 H_k += (
object->vertex(i) *
object->vertex(j)) * tmp;
3979 for (
const unsigned int i :
3981 x_k +=
object->vertex(i) *
3984 if (delta_xi.
norm() < 1e-7)
4000 template <
int structdim>
4007 static const std::size_t n_vertices_per_cell =
4009 n_independent_components;
4010 std::array<double, n_vertices_per_cell> copied_weights;
4011 for (
unsigned int i = 0; i < n_vertices_per_cell; ++i)
4013 copied_weights[i] = v[i];
4014 if (v[i] < 0.0 || v[i] > 1.0)
4019 std::sort(copied_weights.begin(), copied_weights.end());
4021 std::accumulate(copied_weights.begin(), copied_weights.end(), 0.0);
4026 template <
typename Iterator>
4029 const Iterator &
object,
4032 const int spacedim = Iterator::AccessorType::space_dimension;
4033 const int structdim = Iterator::AccessorType::structure_dimension;
4037 if (structdim >= spacedim)
4038 return projected_point;
4039 else if (structdim == 1 || structdim == 2)
4041 using namespace internal::ProjectToObject;
4046 const int dim = Iterator::AccessorType::dimension;
4049 &manifold) !=
nullptr)
4052 project_to_d_linear_object<Iterator, spacedim, structdim>(
4053 object, trial_point);
4098 const double step_size =
object->diameter() / 64.0;
4100 constexpr unsigned int n_vertices_per_cell =
4103 std::array<Point<spacedim>, n_vertices_per_cell>
vertices;
4104 for (
unsigned int vertex_n = 0; vertex_n < n_vertices_per_cell;
4106 vertices[vertex_n] = object->vertex(vertex_n);
4108 auto get_point_from_weights =
4111 return object->get_manifold().get_new_point(
4119 double guess_weights_sum = 0.0;
4120 for (
unsigned int vertex_n = 0; vertex_n < n_vertices_per_cell;
4123 const double distance =
4125 if (distance == 0.0)
4127 guess_weights = 0.0;
4128 guess_weights[vertex_n] = 1.0;
4129 guess_weights_sum = 1.0;
4134 guess_weights[vertex_n] = 1.0 / distance;
4135 guess_weights_sum += guess_weights[vertex_n];
4138 guess_weights /= guess_weights_sum;
4139 Assert(internal::weights_are_ok<structdim>(guess_weights),
4148 for (
unsigned int outer_n = 0; outer_n < 40; ++outer_n)
4150 const unsigned int dependent_direction =
4151 n_vertices_per_cell - 1;
4153 for (
unsigned int row_n = 0; row_n < n_vertices_per_cell;
4156 if (row_n != dependent_direction)
4158 current_gradient[row_n] =
4159 gradient_entry<spacedim, structdim>(
4161 dependent_direction,
4165 get_point_from_weights);
4167 current_gradient[dependent_direction] -=
4168 current_gradient[row_n];
4186 double gradient_weight = -0.5;
4187 auto gradient_weight_objective_function =
4188 [&](
const double gradient_weight_guess) ->
double {
4189 return (trial_point -
4190 get_point_from_weights(guess_weights +
4191 gradient_weight_guess *
4196 for (
unsigned int inner_n = 0; inner_n < 10; ++inner_n)
4198 const double update_numerator = centered_first_difference(
4201 gradient_weight_objective_function);
4202 const double update_denominator =
4203 centered_second_difference(
4206 gradient_weight_objective_function);
4210 if (
std::abs(update_denominator) == 0.0)
4213 gradient_weight - update_numerator / update_denominator;
4218 if (
std::abs(gradient_weight) > 10)
4220 gradient_weight = -10.0;
4230 guess_weights + gradient_weight * current_gradient;
4232 double new_gradient_weight = gradient_weight;
4233 for (
unsigned int iteration_count = 0; iteration_count < 40;
4236 if (internal::weights_are_ok<structdim>(tentative_weights))
4239 for (
unsigned int i = 0; i < n_vertices_per_cell; ++i)
4241 if (tentative_weights[i] < 0.0)
4243 tentative_weights -=
4244 (tentative_weights[i] / current_gradient[i]) *
4247 if (tentative_weights[i] < 0.0 ||
4248 1.0 < tentative_weights[i])
4250 new_gradient_weight /= 2.0;
4253 new_gradient_weight * current_gradient;
4260 if (!internal::weights_are_ok<structdim>(tentative_weights))
4265 if (get_point_from_weights(tentative_weights)
4266 .distance(trial_point) <
4267 get_point_from_weights(guess_weights).distance(trial_point))
4268 guess_weights = tentative_weights;
4271 Assert(internal::weights_are_ok<structdim>(guess_weights),
4274 Assert(internal::weights_are_ok<structdim>(guess_weights),
4276 projected_point = get_point_from_weights(guess_weights);
4286 for (
unsigned int line_n = 0;
4287 line_n < GeometryInfo<structdim>::lines_per_cell;
4290 line_projections[line_n] =
4293 std::sort(line_projections.begin(),
4294 line_projections.end(),
4296 return a.distance(trial_point) <
4297 b.distance(trial_point);
4299 if (line_projections[0].distance(trial_point) <
4300 projected_point.
distance(trial_point))
4301 projected_point = line_projections[0];
4307 return projected_point;
4310 return projected_point;
4317 template <
typename DataType,
4319 typename MeshCellIteratorType>
4322 const MeshType &mesh,
4323 const std::function<
4324 std_cxx17::optional<DataType>(
const MeshCellIteratorType &)> &pack,
4325 const std::function<
void(
const MeshCellIteratorType &,
const DataType &)>
4327 const std::function<
bool(
const MeshCellIteratorType &)> & cell_filter,
4328 const std::function<
void(
4329 const std::function<
void(
const MeshCellIteratorType &,
4331 const std::function<std::set<types::subdomain_id>(
4333 MeshType::space_dimension> &)>
4334 &compute_ghost_owners)
4336# ifndef DEAL_II_WITH_MPI
4341 (void)process_cells;
4342 (void)compute_ghost_owners;
4345 constexpr int dim = MeshType::dimension;
4346 constexpr int spacedim = MeshType::space_dimension;
4349 &mesh.get_triangulation());
4353 "The function exchange_cell_data_to_ghosts() only works with parallel triangulations."));
4355 if (
const auto tria =
dynamic_cast<
4357 &mesh.get_triangulation()))
4360 tria->with_artificial_cells(),
4362 "The functions GridTools::exchange_cell_data_to_ghosts() and "
4363 "GridTools::exchange_cell_data_to_level_ghosts() can only "
4364 "operate on a single layer of ghost cells. However, you have "
4365 "given a Triangulation object of type "
4366 "parallel::shared::Triangulation without artificial cells "
4367 "resulting in arbitrary numbers of ghost layers."));
4371 std::set<::types::subdomain_id> ghost_owners =
4372 compute_ghost_owners(*
tria);
4374 std::vector<typename CellId::binary_type>>
4377 for (
const auto ghost_owner : ghost_owners)
4378 neighbor_cell_list[ghost_owner] = {};
4380 process_cells([&](
const auto &cell,
const auto key) ->
void {
4381 if (cell_filter(cell))
4383 constexpr int spacedim = MeshType::space_dimension;
4384 neighbor_cell_list[key].emplace_back(
4385 cell->id().template to_binary<spacedim>());
4389 Assert(ghost_owners.size() == neighbor_cell_list.size(),
4401 const int mpi_tag_reply =
4405 std::vector<MPI_Request> requests(ghost_owners.size());
4407 unsigned int idx = 0;
4408 for (
const auto &it : neighbor_cell_list)
4411 const int ierr = MPI_Isend(it.second.data(),
4412 it.second.size() *
sizeof(it.second[0]),
4424 std::vector<MPI_Request> reply_requests(ghost_owners.size());
4425 std::vector<std::vector<char>> sendbuffers(ghost_owners.size());
4427 for (
unsigned int idx = 0; idx < ghost_owners.size(); ++idx)
4430 int ierr = MPI_Probe(MPI_ANY_SOURCE,
4437 ierr = MPI_Get_count(&status, MPI_BYTE, &len);
4444 std::vector<typename CellId::binary_type> cells_with_requests(
4446 std::vector<DataType> data_to_send;
4447 data_to_send.reserve(n_cells);
4448 std::vector<bool> cell_carries_data(n_cells,
false);
4450 ierr = MPI_Recv(cells_with_requests.data(),
4460 for (
unsigned int c = 0; c < static_cast<unsigned int>(n_cells); ++c)
4465 MeshCellIteratorType mesh_it(
tria,
4470 const std_cxx17::optional<DataType> data =
pack(mesh_it);
4473 data_to_send.emplace_back(*data);
4474 cell_carries_data[c] =
true;
4482 sendbuffers[idx].resize(
sizeof(std::size_t));
4487 std::size_t size_of_send =
4491 std::memcpy(sendbuffers[idx].data(),
4493 sizeof(std::size_t));
4497 if (data_to_send.size() < n_cells)
4503 ierr = MPI_Isend(sendbuffers[idx].data(),
4504 sendbuffers[idx].size(),
4509 &reply_requests[idx]);
4514 std::vector<char> receive;
4515 for (
unsigned int id = 0;
id < neighbor_cell_list.size(); ++id)
4518 int ierr = MPI_Probe(MPI_ANY_SOURCE,
4525 ierr = MPI_Get_count(&status, MPI_BYTE, &len);
4528 receive.resize(len);
4530 ierr = MPI_Recv(receive.data(),
4541 auto data_iterator = receive.begin();
4542 std::size_t size_of_received_data =
4543 Utilities::unpack<std::size_t>(data_iterator,
4544 data_iterator +
sizeof(std::size_t));
4545 data_iterator +=
sizeof(std::size_t);
4548 auto received_data = Utilities::unpack<std::vector<DataType>>(
4550 data_iterator + size_of_received_data,
4552 data_iterator += size_of_received_data;
4557 const std::vector<typename CellId::binary_type> &this_cell_list =
4558 neighbor_cell_list[status.MPI_SOURCE];
4560 std::vector<bool> cells_with_data;
4561 if (received_data.size() < this_cell_list.size())
4563 cells_with_data = Utilities::unpack<std::vector<bool>>(
4564 data_iterator, receive.end(),
false);
4570 auto received_data_iterator = received_data.begin();
4571 for (
unsigned int c = 0; c < this_cell_list.size(); ++c)
4572 if (cells_with_data.empty() || cells_with_data[c])
4577 MeshCellIteratorType cell(
tria,
4582 unpack(cell, *received_data_iterator);
4583 ++received_data_iterator;
4589 if (requests.size() > 0)
4592 MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
4595 if (reply_requests.size() > 0)
4597 const int ierr = MPI_Waitall(reply_requests.size(),
4598 reply_requests.data(),
4599 MPI_STATUSES_IGNORE);
4609 template <
typename DataType,
typename MeshType>
4612 const MeshType & mesh,
4613 const std::function<std_cxx17::optional<DataType>(
4614 const typename MeshType::active_cell_iterator &)> &pack,
4615 const std::function<
void(
const typename MeshType::active_cell_iterator &,
4616 const DataType &)> & unpack,
4617 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
4620# ifndef DEAL_II_WITH_MPI
4627 internal::exchange_cell_data<DataType,
4629 typename MeshType::active_cell_iterator>(
4634 [&](
const auto &process) {
4635 for (
const auto &cell : mesh.active_cell_iterators())
4636 if (cell->is_ghost())
4637 process(cell, cell->subdomain_id());
4639 [](
const auto &
tria) {
return tria.ghost_owners(); });
4645 template <
typename DataType,
typename MeshType>
4648 const MeshType & mesh,
4649 const std::function<std_cxx17::optional<DataType>(
4650 const typename MeshType::level_cell_iterator &)> &pack,
4651 const std::function<
void(
const typename MeshType::level_cell_iterator &,
4652 const DataType &)> & unpack,
4653 const std::function<
bool(
const typename MeshType::level_cell_iterator &)>
4656# ifndef DEAL_II_WITH_MPI
4663 internal::exchange_cell_data<DataType,
4665 typename MeshType::level_cell_iterator>(
4670 [&](
const auto &process) {
4671 for (
const auto &cell : mesh.cell_iterators())
4672 if (cell->level_subdomain_id() !=
4674 !cell->is_locally_owned_on_level())
4675 process(cell, cell->level_subdomain_id());
4677 [](
const auto &
tria) {
return tria.level_ghost_owners(); });
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
std::array< unsigned int, 4 > binary_type
Abstract base class for mapping classes.
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
numbers::NumberTraits< Number >::real_type norm() const
virtual MPI_Comm get_communicator() const
cell_iterator create_cell_iterator(const CellId &cell_id) const
typename MeshType::active_cell_iterator type
#define DEAL_II_DEPRECATED
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
#define DeclException0(Exception0)
static ::ExceptionBase & ExcMeshNotOrientable()
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)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcTriangulationHasBeenRefined()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcInvalidNumberOfPartitions(int arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcScalingFactorNotPositive(double arg1)
typename ActiveSelector::active_cell_iterator active_cell_iterator
RTree< std::pair< BoundingBox< spacedim >, unsigned int > > build_global_description_tree(const std::vector< BoundingBox< spacedim > > &local_description, const MPI_Comm &mpi_communicator)
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={})
std::vector< std::vector< BoundingBox< spacedim > > > exchange_local_bounding_boxes(const std::vector< BoundingBox< spacedim > > &local_bboxes, const MPI_Comm &mpi_communicator)
std::map< unsigned int, std::set<::types::subdomain_id > > compute_vertices_with_ghost_neighbors(const Triangulation< dim, spacedim > &tria)
void collect_coinciding_vertices(const Triangulation< dim, spacedim > &tria, std::map< unsigned int, std::vector< unsigned int > > &coinciding_vertex_groups, std::map< unsigned int, unsigned int > &vertex_to_coinciding_vertex_group)
void exchange_cell_data_to_ghosts(const MeshType &mesh, const std::function< std_cxx17::optional< DataType >(const typename MeshType::active_cell_iterator &)> &pack, const std::function< void(const typename MeshType::active_cell_iterator &, const DataType &)> &unpack, const std::function< bool(const typename MeshType::active_cell_iterator &)> &cell_filter=always_return< typename MeshType::active_cell_iterator, bool >{true})
void exchange_cell_data_to_level_ghosts(const MeshType &mesh, const std::function< std_cxx17::optional< DataType >(const typename MeshType::level_cell_iterator &)> &pack, const std::function< void(const typename MeshType::level_cell_iterator &, const DataType &)> &unpack, const std::function< bool(const typename MeshType::level_cell_iterator &)> &cell_filter=always_return< typename MeshType::level_cell_iterator, bool >{ true})
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)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
T sum(const T &t, const MPI_Comm &mpi_communicator)
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::subdomain_id artificial_subdomain_id
static const unsigned int invalid_unsigned_int
const types::manifold_id flat_manifold_id
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int global_dof_index
unsigned int subdomain_id
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
boost::geometry::index::rtree< LeafType, IndexType, IndexableGetter > RTree
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
static Tensor< 1, dim > d_linear_shape_function_gradient(const Point< dim > &xi, const unsigned int i)
constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
constexpr SymmetricTensor< 4, dim, Number > outer_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2)
std::map< unsigned int, unsigned int > vertex_to_coinciding_vertex_group
std::map< unsigned int, std::vector< unsigned int > > coinciding_vertex_groups
const ::Triangulation< dim, spacedim > & tria