16#ifndef dealii_grid_tools_h
17# define dealii_grid_tools_h
52# include <boost/archive/binary_iarchive.hpp>
53# include <boost/archive/binary_oarchive.hpp>
54# include <boost/random/mersenne_twister.hpp>
55# include <boost/serialization/array.hpp>
56# include <boost/serialization/vector.hpp>
58# ifdef DEAL_II_WITH_ZLIB
59# include <boost/iostreams/device/back_inserter.hpp>
60# include <boost/iostreams/filter/gzip.hpp>
61# include <boost/iostreams/filtering_stream.hpp>
62# include <boost/iostreams/stream.hpp>
86 class MappingCollection;
94 template <
int dim,
int spacedim,
class MeshType>
99 using type =
typename MeshType::active_cell_iterator;
106 template <
int dim,
int spacedim>
117 template <
int dim,
int spacedim>
118 class ActiveCellIterator<dim, spacedim, ::
hp::DoFHandler<dim, spacedim>>
137 template <
int dim,
int spacedim>
151 template <
int dim,
int spacedim>
181 template <
int dim,
int spacedim>
185 (ReferenceCells::get_hypercube<dim>()
186 .
template get_default_linear_mapping<dim, spacedim>()));
198 template <
int dim,
int spacedim>
203 (ReferenceCells::get_hypercube<dim>()
204 .
template get_default_linear_mapping<dim, spacedim>()));
216 template <
int dim,
int spacedim>
221 (ReferenceCells::get_hypercube<dim>()
222 .
template get_default_linear_mapping<dim, spacedim>()));
284 template <
int dim,
int spacedim>
346 template <
int dim,
int spacedim>
367 template <
typename Iterator>
370 const Iterator &
object,
385 template <
int dim,
int spacedim>
387 tuple<std::vector<Point<spacedim>>, std::vector<CellData<dim>>,
SubCellData>
412 template <
int dim,
int spacedim>
435 template <
int dim,
int spacedim>
440 std::vector<unsigned int> & considered_vertices,
441 const double tol = 1
e-12);
461 template <
int dim,
int spacedim>
539 template <
int dim,
typename Transformation,
int spacedim>
549 template <
int dim,
int spacedim>
583 const unsigned int axis,
648 const bool solve_for_absolute_positions =
false);
655 template <
int dim,
int spacedim>
656 std::map<unsigned int, Point<spacedim>>
666 template <
int dim,
int spacedim>
668 scale(
const double scaling_factor,
685 template <
int dim,
int spacedim>
690 const bool keep_boundary =
true,
691 const unsigned int seed = boost::random::mt19937::default_seed);
726 template <
int dim,
int spacedim>
729 const bool isotropic =
false,
730 const unsigned int max_iterations = 100);
756 template <
int dim,
int spacedim>
759 const double max_ratio = 1.6180339887,
760 const unsigned int max_iterations = 5);
851 template <
int dim,
int spacedim>
854 const double limit_angle_fraction = .75);
916 template <
int dim,
int spacedim>
919 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
920 std::vector<std::vector<Point<dim>>>,
921 std::vector<std::vector<unsigned int>>>
965 template <
int dim,
int spacedim>
968 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
969 std::vector<std::vector<Point<dim>>>,
970 std::vector<std::vector<unsigned int>>,
971 std::vector<unsigned int>>
1051 template <
int dim,
int spacedim>
1054 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
1055 std::vector<std::vector<Point<dim>>>,
1056 std::vector<std::vector<unsigned int>>,
1057 std::vector<std::vector<Point<spacedim>>>,
1058 std::vector<std::vector<unsigned int>>>
1066 const double tolerance = 1
e-10);
1084 template <
int dim,
int spacedim>
1093 std::vector<std::tuple<std::pair<int, int>,
1123 std::vector<std::tuple<unsigned int, unsigned int, unsigned int>>
1148 template <
int dim,
int spacedim>
1154 const double tolerance,
1155 const bool perform_handshake,
1156 const bool enforce_unique_mapping =
false);
1196 template <
int dim,
int spacedim>
1197 std::map<unsigned int, Point<spacedim>>
1201 (ReferenceCells::get_hypercube<dim>()
1202 .
template get_default_linear_mapping<dim, spacedim>()));
1213 template <
int spacedim>
1241 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1245 const std::vector<bool> & marked_vertices = {});
1270 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1273 const MeshType<dim, spacedim> &mesh,
1275 const std::vector<bool> & marked_vertices = {});
1297 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1299 std::vector<typename MeshType<dim, spacedim>::active_cell_iterator>
1302 typename ::internal::
1303 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type>
1306 const unsigned int vertex_index);
1370 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1372 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim>>
1374 std::pair<typename ::internal::
1375 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1379 const MeshType<dim, spacedim> &mesh,
1381 const std::vector<bool> &marked_vertices = {},
1382 const double tolerance = 1.e-10);
1391 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1393 typename MeshType<dim, spacedim>::active_cell_iterator
1395 typename ::internal::
1396 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type
1400 const std::vector<bool> &marked_vertices = {},
1401 const double tolerance = 1.e-10);
1409 template <
int dim,
int spacedim>
1410 std::pair<typename DoFHandler<dim, spacedim>::active_cell_iterator,
1416 const double tolerance = 1.e-10);
1469 template <
int dim,
int spacedim>
1470 std::pair<typename Triangulation<dim, spacedim>::active_cell_iterator,
1473 const Cache<dim, spacedim> &cache,
1477 const std::vector<bool> &marked_vertices = {},
1478 const double tolerance = 1.e-10);
1493 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1495 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim>>
1497 std::pair<typename ::internal::
1498 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1503 const MeshType<dim, spacedim> &mesh,
1506 std::set<
typename MeshType<dim, spacedim>::active_cell_iterator>>
1509 const typename MeshType<dim, spacedim>::active_cell_iterator &cell_hint =
1510 typename MeshType<dim, spacedim>::active_cell_iterator(),
1511 const std::vector<bool> & marked_vertices = {},
1514 const double tolerance = 1.e-10);
1537 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1539 std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
1542 std::vector<std::pair<
1543 typename ::internal::
1544 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1549 const MeshType<dim, spacedim> &mesh,
1551 const double tolerance,
1552 const std::pair<
typename MeshType<dim, spacedim>::active_cell_iterator,
1561 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1563 std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
1566 std::vector<std::pair<
1567 typename ::internal::
1568 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1573 const MeshType<dim, spacedim> &mesh,
1575 const double tolerance = 1
e-10,
1576 const std::vector<bool> & marked_vertices = {});
1599 template <
class MeshType>
1600 std::vector<typename MeshType::active_cell_iterator>
1627 template <
class MeshType>
1630 const typename MeshType::active_cell_iterator & cell,
1631 std::vector<typename MeshType::active_cell_iterator> &active_neighbors);
1682 template <
class MeshType>
1683 std::vector<typename MeshType::active_cell_iterator>
1685 const MeshType &mesh,
1686 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1697 template <
class MeshType>
1698 std::vector<typename MeshType::cell_iterator>
1700 const MeshType &mesh,
1701 const std::function<
bool(
const typename MeshType::cell_iterator &)>
1703 const unsigned int level);
1718 template <
class MeshType>
1719 std::vector<typename MeshType::active_cell_iterator>
1770 template <
class MeshType>
1771 std::vector<typename MeshType::active_cell_iterator>
1773 const MeshType &mesh,
1774 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1776 const double layer_thickness);
1800 template <
class MeshType>
1801 std::vector<typename MeshType::active_cell_iterator>
1803 const double layer_thickness);
1820 template <
class MeshType>
1823 const MeshType &mesh,
1824 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1879 template <
class MeshType>
1880 std::vector<BoundingBox<MeshType::space_dimension>>
1882 const MeshType &mesh,
1883 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1885 const unsigned int refinement_level = 0,
1886 const bool allow_merge =
false,
1916 template <
int spacedim>
1918 std::tuple<std::vector<std::vector<unsigned int>>,
1919 std::map<unsigned int, unsigned int>,
1920 std::map<unsigned int, std::vector<unsigned int>>>
1963 template <
int spacedim>
1965 std::tuple<std::map<unsigned int, std::vector<unsigned int>>,
1966 std::map<unsigned int, unsigned int>,
1967 std::map<unsigned int, std::vector<unsigned int>>>
1984 template <
int dim,
int spacedim>
1986 std::set<typename Triangulation<dim, spacedim>::active_cell_iterator>>
2001 template <
int dim,
int spacedim>
2002 std::vector<std::vector<Tensor<1, spacedim>>>
2017 template <
int dim,
int spacedim>
2023 (ReferenceCells::get_hypercube<dim>()
2024 .
template get_default_linear_mapping<dim, spacedim>()));
2037 template <
int dim,
int spacedim>
2038 std::map<unsigned int, types::global_vertex_index>
2053 template <
int dim,
int spacedim>
2054 std::pair<unsigned int, double>
2072 template <
int dim,
int spacedim>
2086 template <
int dim,
int spacedim>
2100 template <
int dim,
int spacedim>
2104 const unsigned int level,
2127 template <
int dim,
int spacedim>
2144 template <
int dim,
int spacedim>
2147 const std::vector<unsigned int> &cell_weights,
2197 template <
int dim,
int spacedim>
2215 template <
int dim,
int spacedim>
2218 const std::vector<unsigned int> &cell_weights,
2238 template <
int dim,
int spacedim>
2242 const bool group_siblings =
true);
2255 template <
int dim,
int spacedim>
2266 template <
int dim,
int spacedim>
2267 std::vector<types::subdomain_id>
2269 const std::vector<CellId> & cell_ids);
2281 template <
int dim,
int spacedim>
2284 std::vector<types::subdomain_id> & subdomain);
2300 template <
int dim,
int spacedim>
2335 template <
int dim,
int spacedim>
2378 template <
typename MeshType>
2379 std::list<std::pair<
typename MeshType::cell_iterator,
2380 typename MeshType::cell_iterator>>
2392 template <
int dim,
int spacedim>
2406 template <
typename MeshType>
2431 template <
int dim,
int spacedim>
2487 template <
class MeshType>
2488 std::vector<typename MeshType::active_cell_iterator>
2513 template <
class Container>
2514 std::vector<typename Container::cell_iterator>
2516 const std::vector<typename Container::active_cell_iterator> &patch_cells);
2584 template <
class Container>
2587 const std::vector<typename Container::active_cell_iterator> &patch,
2589 &local_triangulation,
2592 Container::space_dimension>::active_cell_iterator,
2593 typename Container::active_cell_iterator> &patch_to_global_tria_map);
2626 template <
int dim,
int spacedim>
2629 std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>>
2645 template <
typename CellIterator>
2745 template <
typename FaceIterator>
2747 std::bitset<3> & orientation,
2748 const FaceIterator & face1,
2749 const FaceIterator & face2,
2750 const int direction,
2759 template <
typename FaceIterator>
2762 const FaceIterator & face1,
2763 const FaceIterator & face2,
2764 const int direction,
2826 template <
typename MeshType>
2829 const MeshType & mesh,
2832 const int direction,
2862 template <
typename MeshType>
2865 const MeshType & mesh,
2867 const int direction,
2870 const ::Tensor<1, MeshType::space_dimension> &offset =
2900 template <
int dim,
int spacedim>
2903 const bool reset_boundary_ids =
false);
2926 template <
int dim,
int spacedim>
2929 const std::vector<types::boundary_id> &src_boundary_ids,
2930 const std::vector<types::manifold_id> &dst_manifold_ids,
2932 const std::vector<types::boundary_id> &reset_boundary_ids = {});
2963 template <
int dim,
int spacedim>
2966 const bool compute_face_ids =
false);
2992 template <
int dim,
int spacedim>
2997 const std::set<types::manifold_id> &)> &disambiguation_function =
2998 [](
const std::set<types::manifold_id> &manifold_ids) {
2999 if (manifold_ids.size() == 1)
3000 return *manifold_ids.begin();
3004 bool overwrite_only_flat_manifold_ids =
true);
3091 template <
typename DataType,
typename MeshType>
3094 const MeshType & mesh,
3095 const std::function<std_cxx17::optional<DataType>(
3096 const typename MeshType::active_cell_iterator &)> &
pack,
3097 const std::function<
void(
const typename MeshType::active_cell_iterator &,
3098 const DataType &)> &
unpack,
3099 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
3113 template <
typename DataType,
typename MeshType>
3116 const MeshType & mesh,
3117 const std::function<std_cxx17::optional<DataType>(
3118 const typename MeshType::level_cell_iterator &)> &
pack,
3119 const std::function<
void(
const typename MeshType::level_cell_iterator &,
3120 const DataType &)> &
unpack,
3121 const std::function<
bool(
const typename MeshType::level_cell_iterator &)> &
3136 template <
int spacedim>
3137 std::vector<std::vector<BoundingBox<spacedim>>>
3140 const MPI_Comm & mpi_communicator);
3174 template <
int spacedim>
3178 const MPI_Comm & mpi_communicator);
3197 template <
int dim,
int spacedim>
3201 std::map<
unsigned int, std::vector<unsigned int>> &coinciding_vertex_groups,
3202 std::map<unsigned int, unsigned int> &vertex_to_coinciding_vertex_group);
3210 template <
int dim,
int spacedim>
3211 std::map<unsigned int, std::set<::types::subdomain_id>>
3227 template <
int dim,
typename T>
3248 template <
class Archive>
3250 save(Archive &ar,
const unsigned int version)
const;
3258 template <
class Archive>
3260 load(Archive &ar,
const unsigned int version);
3268 template <
class Archive>
3274 BOOST_SERIALIZATION_SPLIT_MEMBER()
3297 template <
int dim,
typename VectorType>
3320 const VectorType & ls_vector,
3321 const double iso_level,
3333 const VectorType & ls_vector,
3334 const double iso_level,
3352 const double iso_level,
3364 const std::vector<
Point<2>> & points,
3365 const std::vector<unsigned int> mask,
3366 const double iso_level,
3375 const std::vector<
Point<3>> & points,
3376 const std::vector<unsigned int> mask,
3377 const double iso_level,
3412 <<
"The number of partitions you gave is " << arg1
3413 <<
", but must be greater than zero.");
3419 <<
"The subdomain id " << arg1
3420 <<
" has no cells associated with it.");
3431 <<
"The scaling factor must be positive, but it is " << arg1
3439 <<
"The given vertex with index " << arg1
3440 <<
" is not used in the given triangulation.");
3457 "The edges of the mesh are not consistently orientable.");
3494 template <
int dim,
typename Predicate,
int spacedim>
3499 std::vector<bool> treated_vertices(
triangulation.n_vertices(),
false);
3510 for (; cell != endc; ++cell)
3511 for (
const unsigned int v : cell->vertex_indices())
3512 if (treated_vertices[cell->vertex_index(v)] ==
false)
3515 cell->vertex(v) = predicate(cell->vertex(v));
3517 treated_vertices[cell->vertex_index(v)] =
true;
3527 for (; cell != endc; ++cell)
3528 for (
const unsigned int face : cell->face_indices())
3529 if (cell->face(face)->has_children() &&
3530 !cell->face(face)->at_boundary())
3532 Assert(cell->reference_cell() ==
3533 ReferenceCells::get_hypercube<dim>(),
3537 cell->face(face)->child(0)->vertex(1) =
3538 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1)) /
3547 for (; cell != endc; ++cell)
3548 for (
const unsigned int face : cell->face_indices())
3549 if (cell->face(face)->has_children() &&
3550 !cell->face(face)->at_boundary())
3552 Assert(cell->reference_cell() ==
3553 ReferenceCells::get_hypercube<dim>(),
3557 cell->face(face)->child(0)->vertex(1) =
3558 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1)) /
3560 cell->face(face)->child(0)->vertex(2) =
3561 (cell->face(face)->vertex(0) + cell->face(face)->vertex(2)) /
3563 cell->face(face)->child(1)->vertex(3) =
3564 (cell->face(face)->vertex(1) + cell->face(face)->vertex(3)) /
3566 cell->face(face)->child(2)->vertex(3) =
3567 (cell->face(face)->vertex(2) + cell->face(face)->vertex(3)) /
3571 cell->face(face)->child(0)->vertex(3) =
3572 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1) +
3573 cell->face(face)->vertex(2) + cell->face(face)->vertex(3)) /
3584 template <
class MeshType>
3585 std::vector<typename MeshType::active_cell_iterator>
3588 std::vector<typename MeshType::active_cell_iterator> child_cells;
3590 if (cell->has_children())
3592 for (
unsigned int child = 0; child < cell->n_children(); ++child)
3593 if (cell->child(child)->has_children())
3595 const std::vector<typename MeshType::active_cell_iterator>
3596 children = get_active_child_cells<MeshType>(cell->child(child));
3597 child_cells.insert(child_cells.end(),
3602 child_cells.push_back(cell->child(child));
3610 template <
class MeshType>
3613 const typename MeshType::active_cell_iterator & cell,
3614 std::vector<typename MeshType::active_cell_iterator> &active_neighbors)
3616 active_neighbors.clear();
3617 for (
const unsigned int n : cell->face_indices())
3618 if (!cell->at_boundary(n))
3620 if (MeshType::dimension == 1)
3627 typename MeshType::cell_iterator neighbor_child =
3629 if (!neighbor_child->is_active())
3631 while (neighbor_child->has_children())
3632 neighbor_child = neighbor_child->child(n == 0 ? 1 : 0);
3634 Assert(neighbor_child->neighbor(n == 0 ? 1 : 0) == cell,
3637 active_neighbors.push_back(neighbor_child);
3641 if (cell->face(n)->has_children())
3645 for (
unsigned int c = 0;
3646 c < cell->face(n)->n_active_descendants();
3648 active_neighbors.push_back(
3649 cell->neighbor_child_on_subface(n, c));
3655 active_neighbors.push_back(cell->neighbor(n));
3665 namespace ProjectToObject
3679 struct CrossDerivative
3681 const unsigned int direction_0;
3682 const unsigned int direction_1;
3684 CrossDerivative(
const unsigned int d0,
const unsigned int d1);
3687 inline CrossDerivative::CrossDerivative(
const unsigned int d0,
3688 const unsigned int d1)
3699 template <
typename F>
3701 centered_first_difference(
const double center,
3705 return (f(
center + step) - f(
center - step)) / (2.0 * step);
3714 template <
typename F>
3716 centered_second_difference(
const double center,
3735 template <
int structdim,
typename F>
3738 const CrossDerivative cross_derivative,
3744 simplex_vector[cross_derivative.direction_0] = 0.5 * step;
3745 simplex_vector[cross_derivative.direction_1] = -0.5 * step;
3746 return (-4.0 * f(
center) - 1.0 * f(
center + simplex_vector) -
3747 1.0 / 3.0 * f(
center - simplex_vector) +
3748 16.0 / 3.0 * f(
center + 0.5 * simplex_vector)) /
3760 template <
int spacedim,
int structdim,
typename F>
3763 const unsigned int row_n,
3764 const unsigned int dependent_direction,
3771 dependent_direction <
3773 ExcMessage(
"This function assumes that the last weight is a "
3774 "dependent variable (and hence we cannot take its "
3775 "derivative directly)."));
3776 Assert(row_n != dependent_direction,
3778 "We cannot differentiate with respect to the variable "
3779 "that is assumed to be dependent."));
3783 {row_n, dependent_direction},
center, step, f);
3785 for (
unsigned int dim_n = 0; dim_n < spacedim; ++dim_n)
3787 -2.0 * (p0[dim_n] - manifold_point[dim_n]) * stencil_value[dim_n];
3796 template <
typename Iterator,
int spacedim,
int structdim>
3798 project_to_d_linear_object(
const Iterator &
object,
3833 for (
unsigned int d = 0;
d < structdim; ++
d)
3838 x_k +=
object->vertex(i) *
3844 for (
const unsigned int i :
3847 (x_k - trial_point) * object->vertex(i) *
3852 for (
const unsigned int i :
3854 for (
const unsigned int j :
3862 H_k += (
object->vertex(i) *
object->vertex(j)) * tmp;
3869 for (
const unsigned int i :
3871 x_k +=
object->vertex(i) *
3874 if (delta_xi.
norm() < 1
e-7)
3890 template <
int structdim>
3897 static const std::size_t n_vertices_per_cell =
3899 n_independent_components;
3900 std::array<double, n_vertices_per_cell> copied_weights;
3901 for (
unsigned int i = 0; i < n_vertices_per_cell; ++i)
3903 copied_weights[i] = v[i];
3904 if (v[i] < 0.0 || v[i] > 1.0)
3909 std::sort(copied_weights.begin(), copied_weights.end());
3911 std::accumulate(copied_weights.begin(), copied_weights.end(), 0.0);
3916 template <
typename Iterator>
3919 const Iterator &
object,
3922 const int spacedim = Iterator::AccessorType::space_dimension;
3923 const int structdim = Iterator::AccessorType::structure_dimension;
3927 if (structdim >= spacedim)
3928 return projected_point;
3929 else if (structdim == 1 || structdim == 2)
3931 using namespace internal::ProjectToObject;
3936 const int dim = Iterator::AccessorType::dimension;
3939 &manifold) !=
nullptr)
3942 project_to_d_linear_object<Iterator, spacedim, structdim>(
3943 object, trial_point);
3988 const double step_size =
object->diameter() / 64.0;
3990 constexpr unsigned int n_vertices_per_cell =
3993 std::array<Point<spacedim>, n_vertices_per_cell>
vertices;
3994 for (
unsigned int vertex_n = 0; vertex_n < n_vertices_per_cell;
3996 vertices[vertex_n] = object->vertex(vertex_n);
3998 auto get_point_from_weights =
4001 return object->get_manifold().get_new_point(
4009 double guess_weights_sum = 0.0;
4010 for (
unsigned int vertex_n = 0; vertex_n < n_vertices_per_cell;
4013 const double distance =
4015 if (distance == 0.0)
4017 guess_weights = 0.0;
4018 guess_weights[vertex_n] = 1.0;
4019 guess_weights_sum = 1.0;
4024 guess_weights[vertex_n] = 1.0 / distance;
4025 guess_weights_sum += guess_weights[vertex_n];
4028 guess_weights /= guess_weights_sum;
4029 Assert(internal::weights_are_ok<structdim>(guess_weights),
4038 for (
unsigned int outer_n = 0; outer_n < 40; ++outer_n)
4040 const unsigned int dependent_direction =
4041 n_vertices_per_cell - 1;
4043 for (
unsigned int row_n = 0; row_n < n_vertices_per_cell;
4046 if (row_n != dependent_direction)
4048 current_gradient[row_n] =
4049 gradient_entry<spacedim, structdim>(
4051 dependent_direction,
4055 get_point_from_weights);
4057 current_gradient[dependent_direction] -=
4058 current_gradient[row_n];
4076 double gradient_weight = -0.5;
4077 auto gradient_weight_objective_function =
4078 [&](
const double gradient_weight_guess) ->
double {
4079 return (trial_point -
4080 get_point_from_weights(guess_weights +
4081 gradient_weight_guess *
4086 for (
unsigned int inner_n = 0; inner_n < 10; ++inner_n)
4088 const double update_numerator = centered_first_difference(
4091 gradient_weight_objective_function);
4092 const double update_denominator =
4093 centered_second_difference(
4096 gradient_weight_objective_function);
4100 if (
std::abs(update_denominator) == 0.0)
4103 gradient_weight - update_numerator / update_denominator;
4108 if (
std::abs(gradient_weight) > 10)
4110 gradient_weight = -10.0;
4120 guess_weights + gradient_weight * current_gradient;
4122 double new_gradient_weight = gradient_weight;
4123 for (
unsigned int iteration_count = 0; iteration_count < 40;
4126 if (internal::weights_are_ok<structdim>(tentative_weights))
4129 for (
unsigned int i = 0; i < n_vertices_per_cell; ++i)
4131 if (tentative_weights[i] < 0.0)
4133 tentative_weights -=
4134 (tentative_weights[i] / current_gradient[i]) *
4137 if (tentative_weights[i] < 0.0 ||
4138 1.0 < tentative_weights[i])
4140 new_gradient_weight /= 2.0;
4143 new_gradient_weight * current_gradient;
4150 if (!internal::weights_are_ok<structdim>(tentative_weights))
4155 if (get_point_from_weights(tentative_weights)
4156 .distance(trial_point) <
4157 get_point_from_weights(guess_weights).distance(trial_point))
4158 guess_weights = tentative_weights;
4161 Assert(internal::weights_are_ok<structdim>(guess_weights),
4164 Assert(internal::weights_are_ok<structdim>(guess_weights),
4166 projected_point = get_point_from_weights(guess_weights);
4176 for (
unsigned int line_n = 0;
4177 line_n < GeometryInfo<structdim>::lines_per_cell;
4180 line_projections[line_n] =
4183 std::sort(line_projections.begin(),
4184 line_projections.end(),
4186 return a.distance(trial_point) <
4187 b.distance(trial_point);
4189 if (line_projections[0].distance(trial_point) <
4190 projected_point.
distance(trial_point))
4191 projected_point = line_projections[0];
4197 return projected_point;
4200 return projected_point;
4205 template <
int dim,
typename T>
4206 template <
class Archive>
4208 CellDataTransferBuffer<dim, T>::save(Archive &ar,
4209 const unsigned int )
const
4211 Assert(cell_ids.size() == data.size(),
4214 const std::size_t
n_cells = cell_ids.size();
4216 for (
const auto &
id : cell_ids)
4219 ar & binary_cell_id;
4227 template <
int dim,
typename T>
4228 template <
class Archive>
4230 CellDataTransferBuffer<dim, T>::load(Archive &ar,
4231 const unsigned int )
4237 for (
unsigned int c = 0; c <
n_cells; ++c)
4241 cell_ids.emplace_back(value);
4249 template <
typename DataType,
4251 typename MeshCellIteratorType>
4254 const MeshType &mesh,
4255 const std::function<
4256 std_cxx17::optional<DataType>(
const MeshCellIteratorType &)> &
pack,
4257 const std::function<
void(
const MeshCellIteratorType &,
const DataType &)>
4259 const std::function<
bool(
const MeshCellIteratorType &)> & cell_filter,
4260 const std::function<
void(
4261 const std::function<
void(
const MeshCellIteratorType &,
4263 const std::function<std::set<types::subdomain_id>(
4265 MeshType::space_dimension> &)>
4266 &compute_ghost_owners)
4268# ifndef DEAL_II_WITH_MPI
4273 (void)process_cells;
4274 (void)compute_ghost_owners;
4277 constexpr int dim = MeshType::dimension;
4278 constexpr int spacedim = MeshType::space_dimension;
4281 &mesh.get_triangulation());
4285 "The function exchange_cell_data_to_ghosts() only works with parallel triangulations."));
4287 if (
const auto tria =
dynamic_cast<
4289 &mesh.get_triangulation()))
4292 tria->with_artificial_cells(),
4294 "The functions GridTools::exchange_cell_data_to_ghosts() and "
4295 "GridTools::exchange_cell_data_to_level_ghosts() can only "
4296 "operate on a single layer ghost cells. However, you have "
4297 "given a Triangulation object of type "
4298 "parallel::shared::Triangulation without artificial cells "
4299 "resulting in arbitrary numbers of ghost layers."));
4303 std::set<::types::subdomain_id> ghost_owners =
4304 compute_ghost_owners(*tria);
4306 std::vector<typename CellId::binary_type>>
4309 for (
const auto ghost_owner : ghost_owners)
4310 neighbor_cell_list[ghost_owner] = {};
4312 process_cells([&](
const auto &cell,
const auto key) {
4313 if (cell_filter(cell))
4314 neighbor_cell_list[key].emplace_back(
4315 cell->id().template to_binary<spacedim>());
4318 Assert(ghost_owners.size() == neighbor_cell_list.size(),
4330 const int mpi_tag_reply =
4334 std::vector<MPI_Request> requests(ghost_owners.size());
4336 unsigned int idx = 0;
4337 for (
const auto &it : neighbor_cell_list)
4340 const int ierr = MPI_Isend(it.second.data(),
4341 it.second.size() *
sizeof(it.second[0]),
4352 using DestinationToBufferMap =
4355 DestinationToBufferMap destination_to_data_buffer_map;
4358 std::vector<std::vector<typename CellId::binary_type>> cell_data_to_send(
4359 ghost_owners.size());
4360 std::vector<std::vector<types::global_dof_index>>
4361 send_dof_numbers_and_indices(ghost_owners.size());
4362 std::vector<MPI_Request> reply_requests(ghost_owners.size());
4363 std::vector<std::vector<char>> sendbuffers(ghost_owners.size());
4365 for (
unsigned int idx = 0; idx < ghost_owners.size(); ++idx)
4368 int ierr = MPI_Probe(MPI_ANY_SOURCE,
4375 ierr = MPI_Get_count(&status, MPI_BYTE, &len);
4377 Assert(len %
sizeof(cell_data_to_send[idx][0]) == 0,
4382 cell_data_to_send[idx].resize(
n_cells);
4384 ierr = MPI_Recv(cell_data_to_send[idx].data(),
4394 for (
unsigned int c = 0; c < static_cast<unsigned int>(
n_cells); ++c)
4399 MeshCellIteratorType mesh_it(tria,
4403 const std_cxx17::optional<DataType> data =
pack(mesh_it);
4407 typename DestinationToBufferMap::iterator p =
4408 destination_to_data_buffer_map
4409 .insert(std::make_pair(
4414 p->second.cell_ids.emplace_back(cell->id());
4415 p->second.data.emplace_back(*data);
4421 destination_to_data_buffer_map[idx];
4425 ierr = MPI_Isend(sendbuffers[idx].data(),
4426 sendbuffers[idx].size(),
4431 &reply_requests[idx]);
4436 std::vector<char> receive;
4437 for (
unsigned int idx = 0; idx < ghost_owners.size(); ++idx)
4440 int ierr = MPI_Probe(MPI_ANY_SOURCE,
4447 ierr = MPI_Get_count(&status, MPI_BYTE, &len);
4450 receive.resize(len);
4452 char *ptr = receive.data();
4453 ierr = MPI_Recv(ptr,
4463 Utilities::unpack<CellDataTransferBuffer<dim, DataType>>(
4466 DataType *data = cellinfo.
data.data();
4467 for (
unsigned int c = 0; c < cellinfo.cell_ids.size(); ++c, ++data)
4472 MeshCellIteratorType cell(tria,
4483 if (requests.size() > 0)
4486 MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
4489 if (reply_requests.size() > 0)
4491 const int ierr = MPI_Waitall(reply_requests.size(),
4492 reply_requests.data(),
4493 MPI_STATUSES_IGNORE);
4503 template <
typename DataType,
typename MeshType>
4506 const MeshType & mesh,
4507 const std::function<std_cxx17::optional<DataType>(
4508 const typename MeshType::active_cell_iterator &)> &
pack,
4509 const std::function<
void(
const typename MeshType::active_cell_iterator &,
4510 const DataType &)> &
unpack,
4511 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
4514# ifndef DEAL_II_WITH_MPI
4521 internal::exchange_cell_data<DataType,
4523 typename MeshType::active_cell_iterator>(
4528 [&](
const auto &process) {
4529 for (
const auto &cell : mesh.active_cell_iterators())
4530 if (cell->is_ghost())
4531 process(cell, cell->subdomain_id());
4533 [](
const auto &tria) {
return tria.ghost_owners(); });
4539 template <
typename DataType,
typename MeshType>
4542 const MeshType & mesh,
4543 const std::function<std_cxx17::optional<DataType>(
4544 const typename MeshType::level_cell_iterator &)> &
pack,
4545 const std::function<
void(
const typename MeshType::level_cell_iterator &,
4546 const DataType &)> &
unpack,
4547 const std::function<
bool(
const typename MeshType::level_cell_iterator &)>
4550# ifndef DEAL_II_WITH_MPI
4557 internal::exchange_cell_data<DataType,
4559 typename MeshType::level_cell_iterator>(
4564 [&](
const auto &process) {
4565 for (
const auto &cell : mesh.cell_iterators())
4566 if (cell->level_subdomain_id() !=
4568 !cell->is_locally_owned_on_level())
4569 process(cell, cell->level_subdomain_id());
4571 [](
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
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &t)
constexpr SymmetricTensor< 4, dim, Number > outer_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2)
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_NAMESPACE_OPEN
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_DEPRECATED_EARLY
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
__global__ void set(Number *val, const Number s, const size_type N)
#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 AssertThrowMPI(error_code)
#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)
@ matrix
Contents is actually a matrix.
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Tensor< 2, dim, Number > F(const Tensor< 2, dim, Number > &Grad_u)
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)