16 #ifndef dealii_grid_tools_h 17 # define dealii_grid_tools_h 20 # include <deal.II/base/config.h> 22 # include <deal.II/base/bounding_box.h> 23 # include <deal.II/base/geometry_info.h> 25 # include <deal.II/boost_adaptors/bounding_box.h> 27 # include <deal.II/dofs/dof_handler.h> 29 # include <deal.II/fe/mapping.h> 30 # include <deal.II/fe/mapping_q1.h> 32 # include <deal.II/grid/manifold.h> 33 # include <deal.II/grid/tria.h> 34 # include <deal.II/grid/tria_accessor.h> 35 # include <deal.II/grid/tria_iterator.h> 37 # include <deal.II/hp/dof_handler.h> 39 # include <deal.II/lac/sparsity_tools.h> 41 # include <deal.II/numerics/rtree.h> 43 # include <boost/archive/binary_iarchive.hpp> 44 # include <boost/archive/binary_oarchive.hpp> 45 # include <boost/geometry/index/detail/serialization.hpp> 46 # include <boost/geometry/index/rtree.hpp> 47 # include <boost/optional.hpp> 48 # include <boost/serialization/array.hpp> 49 # include <boost/serialization/vector.hpp> 51 # ifdef DEAL_II_WITH_ZLIB 52 # include <boost/iostreams/device/back_inserter.hpp> 53 # include <boost/iostreams/filter/gzip.hpp> 54 # include <boost/iostreams/filtering_stream.hpp> 55 # include <boost/iostreams/stream.hpp> 63 DEAL_II_NAMESPACE_OPEN
77 class MappingCollection;
84 template <
int dim,
int spacedim,
class MeshType>
85 class ActiveCellIterator
89 using type =
typename MeshType::active_cell_iterator;
96 template <
int dim,
int spacedim>
97 class ActiveCellIterator<dim, spacedim, ::
DoFHandler<dim, spacedim>>
104 template <
int dim,
int spacedim>
105 class ActiveCellIterator<dim, spacedim, ::
hp::DoFHandler<dim, spacedim>>
124 template <
int dim,
int spacedim>
138 template <
int dim,
int spacedim>
168 template <
int dim,
int spacedim>
184 template <
int dim,
int spacedim>
200 template <
int dim,
int spacedim>
226 template <
int dim,
typename T>
243 template <
int dim,
int spacedim>
266 template <
typename Iterator>
269 const Iterator &
object,
280 template <
int dim,
int spacedim>
282 tuple<std::vector<Point<spacedim>>, std::vector<CellData<dim>>,
SubCellData>
307 template <
int dim,
int spacedim>
329 template <
int dim,
int spacedim>
334 std::vector<unsigned int> & considered_vertices,
335 const double tol = 1e-12);
377 template <
int dim,
typename Transformation,
int spacedim>
379 transform(
const Transformation & transformation,
387 template <
int dim,
int spacedim>
417 rotate(
const double angle,
418 const unsigned int axis,
483 const bool solve_for_absolute_positions =
false);
490 template <
int dim,
int spacedim>
491 std::map<unsigned int, Point<spacedim>>
501 template <
int dim,
int spacedim>
503 scale(
const double scaling_factor,
516 template <
int dim,
int spacedim>
520 const bool keep_boundary =
true);
557 template <
int dim,
int spacedim>
560 const bool isotropic =
false,
561 const unsigned int max_iterations = 100);
589 template <
int dim,
int spacedim>
592 const double max_ratio = 1.6180339887,
593 const unsigned int max_iterations = 5);
686 template <
int dim,
int spacedim>
689 const double limit_angle_fraction = .75);
748 template <
int dim,
int spacedim>
751 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
752 std::vector<std::vector<Point<dim>>>,
753 std::vector<std::vector<unsigned int>>>
792 template <
int dim,
int spacedim>
795 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
796 std::vector<std::vector<Point<dim>>>,
797 std::vector<std::vector<unsigned int>>,
798 std::vector<unsigned int>>
873 template <
int dim,
int spacedim>
876 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
877 std::vector<std::vector<Point<dim>>>,
878 std::vector<std::vector<unsigned int>>,
879 std::vector<std::vector<Point<spacedim>>>,
880 std::vector<std::vector<unsigned int>>>
927 template <
int dim,
int spacedim>
928 std::map<unsigned int, Point<spacedim>>
944 template <
int spacedim>
975 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
979 const std::vector<bool> & marked_vertices = {});
1006 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1009 const MeshType<dim, spacedim> &mesh,
1011 const std::vector<bool> & marked_vertices = {});
1038 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1040 std::vector<typename MeshType<dim, spacedim>::active_cell_iterator>
1043 typename ::internal::
1044 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type>
1047 const unsigned int vertex_index);
1074 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1076 typename MeshType<dim, spacedim>::active_cell_iterator
1078 typename ::internal::ActiveCellIterator<dim,
1080 MeshType<dim, spacedim>>::type
1084 const std::vector<bool> &marked_vertices = {});
1173 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1175 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim>>
1177 std::pair<typename ::internal::
1178 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1182 const MeshType<dim, spacedim> &mesh,
1184 const std::vector<bool> &marked_vertices = {});
1197 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1199 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim>>
1201 std::pair<typename ::internal::
1202 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1207 const MeshType<dim, spacedim> &mesh,
1210 std::set<
typename MeshType<dim, spacedim>::active_cell_iterator>>
1213 const typename MeshType<dim, spacedim>::active_cell_iterator &cell_hint =
1214 typename MeshType<dim, spacedim>::active_cell_iterator(),
1215 const std::vector<bool> & marked_vertices = {},
1216 const RTree<std::pair<Point<spacedim>,
unsigned int>> &used_vertices_rtree =
1217 RTree<std::pair<Point<spacedim>,
unsigned int>>{});
1240 template <
int dim,
int spacedim>
1241 std::pair<typename hp::DoFHandler<dim, spacedim>::active_cell_iterator,
1254 template <
int dim,
int spacedim>
1255 std::pair<typename Triangulation<dim, spacedim>::active_cell_iterator,
1258 const Cache<dim, spacedim> &cache,
1262 const std::vector<bool> &marked_vertices = {});
1281 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1283 std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
1286 std::vector<std::pair<
1287 typename ::internal::
1288 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1293 const MeshType<dim, spacedim> &mesh,
1295 const double tolerance = 1e-12,
1296 const std::vector<bool> & marked_vertices = {});
1319 template <
class MeshType>
1320 std::vector<typename MeshType::active_cell_iterator>
1340 template <
class MeshType>
1343 const typename MeshType::active_cell_iterator & cell,
1344 std::vector<typename MeshType::active_cell_iterator> &active_neighbors);
1398 template <
class MeshType>
1399 std::vector<typename MeshType::active_cell_iterator>
1401 const MeshType &mesh,
1402 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1413 template <
class MeshType>
1414 std::vector<typename MeshType::cell_iterator>
1416 const MeshType &mesh,
1417 const std::function<
bool(
const typename MeshType::cell_iterator &)>
1419 const unsigned int level);
1437 template <
class MeshType>
1438 std::vector<typename MeshType::active_cell_iterator>
1492 template <
class MeshType>
1493 std::vector<typename MeshType::active_cell_iterator>
1495 const MeshType &mesh,
1496 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1498 const double layer_thickness);
1525 template <
class MeshType>
1526 std::vector<typename MeshType::active_cell_iterator>
1528 const double layer_thickness);
1545 template <
class MeshType>
1548 const MeshType &mesh,
1549 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1604 template <
class MeshType>
1605 std::vector<BoundingBox<MeshType::space_dimension>>
1607 const MeshType &mesh,
1608 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1610 const unsigned int refinement_level = 0,
1611 const bool allow_merge =
false,
1643 template <
int spacedim>
1645 std::tuple<std::vector<std::vector<unsigned int>>,
1646 std::map<unsigned int, unsigned int>,
1647 std::map<unsigned int, std::vector<unsigned int>>>
1690 template <
int spacedim>
1692 std::tuple<std::map<unsigned int, std::vector<unsigned int>>,
1693 std::map<unsigned int, unsigned int>,
1694 std::map<unsigned int, std::vector<unsigned int>>>
1711 template <
int dim,
int spacedim>
1713 std::set<typename Triangulation<dim, spacedim>::active_cell_iterator>>
1730 template <
int dim,
int spacedim>
1731 std::vector<std::vector<Tensor<1, spacedim>>>
1745 template <
int dim,
int spacedim>
1762 template <
int dim,
int spacedim>
1763 std::map<unsigned int, types::global_vertex_index>
1780 template <
int dim,
int spacedim>
1781 std::pair<unsigned int, double>
1799 template <
int dim,
int spacedim>
1813 template <
int dim,
int spacedim>
1827 template <
int dim,
int spacedim>
1831 const unsigned int level,
1854 template <
int dim,
int spacedim>
1871 template <
int dim,
int spacedim>
1874 const std::vector<unsigned int> &cell_weights,
1924 template <
int dim,
int spacedim>
1942 template <
int dim,
int spacedim>
1945 const std::vector<unsigned int> &cell_weights,
1958 template <
int dim,
int spacedim>
1974 template <
int dim,
int spacedim>
1988 template <
int dim,
int spacedim>
1991 std::vector<types::subdomain_id> & subdomain);
2007 template <
int dim,
int spacedim>
2043 template <
int dim,
int spacedim>
2086 template <
typename MeshType>
2087 std::list<std::pair<
typename MeshType::cell_iterator,
2088 typename MeshType::cell_iterator>>
2100 template <
int dim,
int spacedim>
2114 template <
typename MeshType>
2139 template <
int dim,
int spacedim>
2197 template <
class MeshType>
2198 std::vector<typename MeshType::active_cell_iterator>
2225 template <
class Container>
2226 std::vector<typename Container::cell_iterator>
2228 const std::vector<typename Container::active_cell_iterator> &patch_cells);
2298 template <
class Container>
2301 const std::vector<typename Container::active_cell_iterator> &patch,
2303 &local_triangulation,
2306 Container::space_dimension>::active_cell_iterator,
2307 typename Container::active_cell_iterator> &patch_to_global_tria_map);
2345 template <
class DoFHandlerType>
2347 std::vector<typename DoFHandlerType::active_cell_iterator>>
2363 template <
typename CellIterator>
2364 struct PeriodicFacePair
2465 template <
typename FaceIterator>
2467 std::bitset<3> & orientation,
2468 const FaceIterator & face1,
2469 const FaceIterator & face2,
2470 const int direction,
2479 template <
typename FaceIterator>
2482 const FaceIterator & face1,
2483 const FaceIterator & face2,
2484 const int direction,
2548 template <
typename MeshType>
2551 const MeshType & mesh,
2554 const int direction,
2584 template <
typename MeshType>
2587 const MeshType & mesh,
2589 const int direction,
2592 const ::Tensor<1, MeshType::space_dimension> &offset =
2624 template <
int dim,
int spacedim>
2627 const bool reset_boundary_ids =
false);
2652 template <
int dim,
int spacedim>
2655 const std::vector<types::boundary_id> &src_boundary_ids,
2656 const std::vector<types::manifold_id> &dst_manifold_ids,
2658 const std::vector<types::boundary_id> &reset_boundary_ids = {});
2691 template <
int dim,
int spacedim>
2694 const bool compute_face_ids =
false);
2722 template <
int dim,
int spacedim>
2727 const std::set<types::manifold_id> &)> &disambiguation_function =
2728 [](
const std::set<types::manifold_id> &manifold_ids) {
2729 if (manifold_ids.size() == 1)
2730 return *manifold_ids.begin();
2734 bool overwrite_only_flat_manifold_ids =
true);
2819 template <
typename DataType,
typename MeshType>
2822 const MeshType & mesh,
2823 const std::function<boost::optional<DataType>(
2824 const typename MeshType::active_cell_iterator &)> &pack,
2825 const std::function<
void(
const typename MeshType::active_cell_iterator &,
2826 const DataType &)> & unpack);
2839 template <
int spacedim>
2840 std::vector<std::vector<BoundingBox<spacedim>>>
2841 exchange_local_bounding_boxes(
2843 MPI_Comm mpi_communicator);
2879 template <
int spacedim>
2880 RTree<std::pair<BoundingBox<spacedim>,
unsigned int>>
2883 MPI_Comm mpi_communicator);
2897 template <
int dim,
typename T>
2917 template <
class Archive>
2919 save(Archive &ar,
const unsigned int version)
const;
2925 template <
class Archive>
2927 load(Archive &ar,
const unsigned int version);
2929 BOOST_SERIALIZATION_SPLIT_MEMBER()
2942 <<
"The number of partitions you gave is " << arg1
2943 <<
", but must be greater than zero.");
2949 <<
"The subdomain id " << arg1
2950 <<
" has no cells associated with it.");
2961 <<
"The scaling factor must be positive, but it is " << arg1
2969 <<
"The point <" << arg1
2970 <<
"> could not be found inside any of the " 2971 <<
"coarse grid cells.");
2978 <<
"The point <" << arg1
2979 <<
"> could not be found inside any of the " 2980 <<
"subcells of a coarse grid cell.");
2987 <<
"The given vertex with index " << arg1
2988 <<
" is not used in the given triangulation.");
3003 template <
int dim,
typename T>
3008 return std::numeric_limits<double>::quiet_NaN();
3011 template <
int dim,
typename Predicate,
int spacedim>
3016 std::vector<bool> treated_vertices(triangulation.n_vertices(),
false);
3025 cell = triangulation.begin_active(),
3026 endc = triangulation.end();
3027 for (; cell != endc; ++cell)
3028 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell; ++v)
3029 if (treated_vertices[cell->vertex_index(v)] ==
false)
3032 cell->vertex(v) = predicate(cell->vertex(v));
3034 treated_vertices[cell->vertex_index(v)] =
true;
3042 cell = triangulation.begin_active(),
3043 endc = triangulation.end();
3044 for (; cell != endc; ++cell)
3045 for (
unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell;
3047 if (cell->face(face)->has_children() &&
3048 !cell->face(face)->at_boundary())
3051 cell->face(face)->child(0)->vertex(1) =
3052 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1)) /
3059 cell = triangulation.begin_active(),
3060 endc = triangulation.end();
3061 for (; cell != endc; ++cell)
3062 for (
unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell;
3064 if (cell->face(face)->has_children() &&
3065 !cell->face(face)->at_boundary())
3068 cell->face(face)->child(0)->vertex(1) =
3069 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1)) /
3071 cell->face(face)->child(0)->vertex(2) =
3072 (cell->face(face)->vertex(0) + cell->face(face)->vertex(2)) /
3074 cell->face(face)->child(1)->vertex(3) =
3075 (cell->face(face)->vertex(1) + cell->face(face)->vertex(3)) /
3077 cell->face(face)->child(2)->vertex(3) =
3078 (cell->face(face)->vertex(2) + cell->face(face)->vertex(3)) /
3082 cell->face(face)->child(0)->vertex(3) =
3083 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1) +
3084 cell->face(face)->vertex(2) + cell->face(face)->vertex(3)) /
3090 triangulation.signals.mesh_movement();
3095 template <
class MeshType>
3096 std::vector<typename MeshType::active_cell_iterator>
3099 std::vector<typename MeshType::active_cell_iterator> child_cells;
3101 if (cell->has_children())
3103 for (
unsigned int child = 0; child < cell->n_children(); ++child)
3104 if (cell->child(child)->has_children())
3106 const std::vector<typename MeshType::active_cell_iterator>
3107 children = get_active_child_cells<MeshType>(cell->child(child));
3108 child_cells.insert(child_cells.end(),
3113 child_cells.push_back(cell->child(child));
3121 template <
class MeshType>
3124 const typename MeshType::active_cell_iterator & cell,
3125 std::vector<typename MeshType::active_cell_iterator> &active_neighbors)
3127 active_neighbors.clear();
3128 for (
unsigned int n = 0;
3129 n < GeometryInfo<MeshType::dimension>::faces_per_cell;
3131 if (!cell->at_boundary(n))
3133 if (MeshType::dimension == 1)
3140 typename MeshType::cell_iterator neighbor_child =
3142 if (!neighbor_child->active())
3144 while (neighbor_child->has_children())
3145 neighbor_child = neighbor_child->child(n == 0 ? 1 : 0);
3147 Assert(neighbor_child->neighbor(n == 0 ? 1 : 0) == cell,
3150 active_neighbors.push_back(neighbor_child);
3154 if (cell->face(n)->has_children())
3158 for (
unsigned int c = 0;
3159 c < cell->face(n)->number_of_children();
3161 active_neighbors.push_back(
3162 cell->neighbor_child_on_subface(n, c));
3168 active_neighbors.push_back(cell->neighbor(n));
3178 namespace ProjectToObject
3192 struct CrossDerivative
3194 const unsigned int direction_0;
3195 const unsigned int direction_1;
3197 CrossDerivative(
const unsigned int d0,
const unsigned int d1);
3200 inline CrossDerivative::CrossDerivative(
const unsigned int d0,
3201 const unsigned int d1)
3212 template <
typename F>
3214 centered_first_difference(
const double center,
3216 const F &f) -> decltype(f(center) - f(center))
3218 return (f(center + step) - f(center - step)) / (2.0 * step);
3227 template <
typename F>
3229 centered_second_difference(
const double center,
3231 const F &f) -> decltype(f(center) - f(center))
3233 return (f(center + step) - 2.0 * f(center) + f(center - step)) /
3248 template <
int structdim,
typename F>
3251 const CrossDerivative cross_derivative,
3254 const F &f) -> decltype(f(center) - f(center))
3257 simplex_vector[cross_derivative.direction_0] = 0.5 * step;
3258 simplex_vector[cross_derivative.direction_1] = -0.5 * step;
3259 return (-4.0 * f(center) - 1.0 * f(center + simplex_vector) -
3260 1.0 / 3.0 * f(center - simplex_vector) +
3261 16.0 / 3.0 * f(center + 0.5 * simplex_vector)) /
3273 template <
int spacedim,
int structdim,
typename F>
3276 const unsigned int row_n,
3277 const unsigned int dependent_direction,
3284 dependent_direction <
3286 ExcMessage(
"This function assumes that the last weight is a " 3287 "dependent variable (and hence we cannot take its " 3288 "derivative directly)."));
3289 Assert(row_n != dependent_direction,
3291 "We cannot differentiate with respect to the variable " 3292 "that is assumed to be dependent."));
3296 {row_n, dependent_direction}, center, step, f);
3298 for (
unsigned int dim_n = 0; dim_n < spacedim; ++dim_n)
3300 -2.0 * (p0[dim_n] - manifold_point[dim_n]) * stencil_value[dim_n];
3309 template <
typename Iterator,
int spacedim,
int structdim>
3311 project_to_d_linear_object(
const Iterator &
object,
3346 for (
unsigned int d = 0;
d < structdim; ++
d)
3350 for (
unsigned int i = 0; i < GeometryInfo<structdim>::vertices_per_cell;
3352 x_k += object->vertex(i) *
3358 for (
unsigned int i = 0;
3359 i < GeometryInfo<structdim>::vertices_per_cell;
3362 (x_k - trial_point) *
object->vertex(i) *
3367 for (
unsigned int i = 0;
3368 i < GeometryInfo<structdim>::vertices_per_cell;
3370 for (
unsigned int j = 0;
3371 j < GeometryInfo<structdim>::vertices_per_cell;
3379 H_k += (
object->vertex(i) *
object->vertex(j)) * tmp;
3386 for (
unsigned int i = 0;
3387 i < GeometryInfo<structdim>::vertices_per_cell;
3389 x_k += object->vertex(i) *
3392 if (delta_xi.
norm() < 1
e-7)
3408 template <
int structdim>
3415 static const std::size_t n_vertices_per_cell =
3417 n_independent_components;
3418 std::array<double, n_vertices_per_cell> copied_weights;
3419 for (
unsigned int i = 0; i < n_vertices_per_cell; ++i)
3421 copied_weights[i] = v[i];
3422 if (v[i] < 0.0 || v[i] > 1.0)
3427 std::sort(copied_weights.begin(), copied_weights.end());
3429 std::accumulate(copied_weights.begin(), copied_weights.end(), 0.0);
3430 return std::abs(sum - 1.0) < 1
e-10;
3434 template <
typename Iterator>
3437 const Iterator &
object,
3440 const int spacedim = Iterator::AccessorType::space_dimension;
3441 const int structdim = Iterator::AccessorType::structure_dimension;
3445 if (structdim >= spacedim)
3446 return projected_point;
3447 else if (structdim == 1 || structdim == 2)
3449 using namespace internal::ProjectToObject;
3454 const int dim = Iterator::AccessorType::dimension;
3457 &manifold) !=
nullptr)
3460 project_to_d_linear_object<Iterator, spacedim, structdim>(
3461 object, trial_point);
3506 const double step_size =
object->diameter() / 64.0;
3508 constexpr
unsigned int n_vertices_per_cell =
3511 std::array<Point<spacedim>, n_vertices_per_cell> vertices;
3512 for (
unsigned int vertex_n = 0; vertex_n < n_vertices_per_cell;
3514 vertices[vertex_n] = object->vertex(vertex_n);
3516 auto get_point_from_weights =
3519 return object->get_manifold().get_new_point(
3520 make_array_view(vertices.begin(), vertices.end()),
3521 make_array_view(weights.begin_raw(), weights.end_raw()));
3527 double guess_weights_sum = 0.0;
3528 for (
unsigned int vertex_n = 0; vertex_n < n_vertices_per_cell;
3531 const double distance =
3532 vertices[vertex_n].
distance(trial_point);
3533 if (distance == 0.0)
3535 guess_weights = 0.0;
3536 guess_weights[vertex_n] = 1.0;
3537 guess_weights_sum = 1.0;
3542 guess_weights[vertex_n] = 1.0 / distance;
3543 guess_weights_sum += guess_weights[vertex_n];
3546 guess_weights /= guess_weights_sum;
3547 Assert(internal::weights_are_ok<structdim>(guess_weights),
3556 for (
unsigned int outer_n = 0; outer_n < 40; ++outer_n)
3558 const unsigned int dependent_direction =
3559 n_vertices_per_cell - 1;
3561 for (
unsigned int row_n = 0; row_n < n_vertices_per_cell;
3564 if (row_n != dependent_direction)
3566 current_gradient[row_n] =
3567 gradient_entry<spacedim, structdim>(
3569 dependent_direction,
3573 get_point_from_weights);
3575 current_gradient[dependent_direction] -=
3576 current_gradient[row_n];
3594 double gradient_weight = -0.5;
3595 auto gradient_weight_objective_function =
3596 [&](
const double gradient_weight_guess) ->
double {
3597 return (trial_point -
3598 get_point_from_weights(guess_weights +
3599 gradient_weight_guess *
3604 for (
unsigned int inner_n = 0; inner_n < 10; ++inner_n)
3606 const double update_numerator = centered_first_difference(
3609 gradient_weight_objective_function);
3610 const double update_denominator =
3611 centered_second_difference(
3614 gradient_weight_objective_function);
3618 if (std::abs(update_denominator) == 0.0)
3621 gradient_weight - update_numerator / update_denominator;
3626 if (std::abs(gradient_weight) > 10)
3628 gradient_weight = -10.0;
3638 guess_weights + gradient_weight * current_gradient;
3640 double new_gradient_weight = gradient_weight;
3641 for (
unsigned int iteration_count = 0; iteration_count < 40;
3644 if (internal::weights_are_ok<structdim>(tentative_weights))
3647 for (
unsigned int i = 0; i < n_vertices_per_cell; ++i)
3649 if (tentative_weights[i] < 0.0)
3651 tentative_weights -=
3652 (tentative_weights[i] / current_gradient[i]) *
3655 if (tentative_weights[i] < 0.0 ||
3656 1.0 < tentative_weights[i])
3658 new_gradient_weight /= 2.0;
3661 new_gradient_weight * current_gradient;
3668 if (!internal::weights_are_ok<structdim>(tentative_weights))
3673 if (get_point_from_weights(tentative_weights)
3674 .distance(trial_point) <
3675 get_point_from_weights(guess_weights).distance(trial_point))
3676 guess_weights = tentative_weights;
3679 Assert(internal::weights_are_ok<structdim>(guess_weights),
3682 Assert(internal::weights_are_ok<structdim>(guess_weights),
3684 projected_point = get_point_from_weights(guess_weights);
3694 for (
unsigned int line_n = 0;
3695 line_n < GeometryInfo<structdim>::lines_per_cell;
3698 line_projections[line_n] =
3701 std::sort(line_projections.begin(),
3702 line_projections.end(),
3705 b.distance(trial_point);
3707 if (line_projections[0].distance(trial_point) <
3708 projected_point.
distance(trial_point))
3709 projected_point = line_projections[0];
3715 return projected_point;
3718 return projected_point;
3723 template <
int dim,
typename T>
3724 template <
class Archive>
3726 CellDataTransferBuffer<dim, T>::save(Archive &ar,
3727 const unsigned int )
const 3729 Assert(cell_ids.size() == data.size(),
3732 const std::size_t n_cells = cell_ids.size();
3734 for (
const auto &
id : cell_ids)
3737 ar & binary_cell_id;
3745 template <
int dim,
typename T>
3746 template <
class Archive>
3748 CellDataTransferBuffer<dim, T>::load(Archive &ar,
3749 const unsigned int )
3751 std::size_t n_cells;
3754 cell_ids.reserve(n_cells);
3755 for (
unsigned int c = 0; c < n_cells; ++c)
3759 cell_ids.emplace_back(value);
3766 template <
typename DataType,
typename MeshType>
3769 const MeshType & mesh,
3770 const std::function<boost::optional<DataType>(
3771 const typename MeshType::active_cell_iterator &)> &pack,
3772 const std::function<
void(
const typename MeshType::active_cell_iterator &,
3773 const DataType &)> & unpack)
3775 # ifndef DEAL_II_WITH_MPI 3781 "GridTools::exchange_cell_data_to_ghosts() requires MPI."));
3783 constexpr
int dim = MeshType::dimension;
3784 constexpr
int spacedim = MeshType::space_dimension;
3786 &mesh.get_triangulation());
3790 "The function exchange_cell_data_to_ghosts() only works with parallel triangulations."));
3793 using DestinationToBufferMap =
3795 CellDataTransferBuffer<dim, DataType>>;
3796 DestinationToBufferMap destination_to_data_buffer_map;
3798 std::map<unsigned int, std::set<::types::subdomain_id>>
3799 vertices_with_ghost_neighbors =
3800 tria->compute_vertices_with_ghost_neighbors();
3803 if (cell->is_locally_owned())
3805 std::set<::types::subdomain_id> send_to;
3806 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell;
3809 const std::map<
unsigned int,
3810 std::set<::types::subdomain_id>>::
3811 const_iterator neighbor_subdomains_of_vertex =
3812 vertices_with_ghost_neighbors.find(cell->vertex_index(v));
3814 if (neighbor_subdomains_of_vertex ==
3815 vertices_with_ghost_neighbors.end())
3818 Assert(neighbor_subdomains_of_vertex->second.size() != 0,
3821 send_to.insert(neighbor_subdomains_of_vertex->second.begin(),
3822 neighbor_subdomains_of_vertex->second.end());
3825 if (send_to.size() > 0)
3828 typename MeshType::active_cell_iterator mesh_it(tria,
3833 const boost::optional<DataType> data =
pack(mesh_it);
3837 const CellId cellid = cell->id();
3839 for (
const auto subdomain : send_to)
3843 typename DestinationToBufferMap::iterator p =
3844 destination_to_data_buffer_map
3845 .insert(std::make_pair(
3846 subdomain, CellDataTransferBuffer<dim, DataType>()))
3849 p->second.cell_ids.emplace_back(cellid);
3850 p->second.data.emplace_back(data.get());
3858 std::set<::types::subdomain_id> ghost_owners = tria->ghost_owners();
3859 const unsigned int n_ghost_owners = ghost_owners.size();
3860 std::vector<std::vector<char>> sendbuffers(n_ghost_owners);
3861 std::vector<MPI_Request> requests(n_ghost_owners);
3863 unsigned int idx = 0;
3864 for (
auto it = ghost_owners.begin(); it != ghost_owners.end(); ++it, ++idx)
3866 CellDataTransferBuffer<dim, DataType> &data =
3867 destination_to_data_buffer_map[*it];
3873 const int ierr = MPI_Isend(sendbuffers[idx].data(),
3874 sendbuffers[idx].size(),
3878 tria->get_communicator(),
3884 std::vector<char> receive;
3885 for (
unsigned int idx = 0; idx < n_ghost_owners; ++idx)
3890 MPI_Probe(MPI_ANY_SOURCE, 786, tria->get_communicator(), &status);
3892 ierr = MPI_Get_count(&status, MPI_BYTE, &len);
3895 receive.resize(len);
3897 char *ptr = receive.data();
3898 ierr = MPI_Recv(ptr,
3903 tria->get_communicator(),
3908 Utilities::unpack<CellDataTransferBuffer<dim, DataType>>(receive);
3910 DataType *data = cellinfo.data.data();
3911 for (
unsigned int c = 0; c < cellinfo.cell_ids.size(); ++c, ++data)
3914 tria_cell = cellinfo.cell_ids[c].to_cell(*tria);
3916 const typename MeshType::active_cell_iterator cell(
3917 tria, tria_cell->level(), tria_cell->index(), &mesh);
3925 if (requests.size())
3928 MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
3931 # endif // DEAL_II_WITH_MPI 3937 DEAL_II_NAMESPACE_CLOSE
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={})
static ::ExceptionBase & ExcScalingFactorNotPositive(double arg1)
const types::manifold_id flat_manifold_id
static const unsigned int invalid_unsigned_int
void copy_boundary_to_manifold_id(Triangulation< dim, spacedim > &tria, const bool reset_boundary_ids=false)
static ::ExceptionBase & ExcPointNotFoundInCoarseGrid(Point< N > arg1)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
IteratorRange< active_cell_iterator > active_cell_iterators() const
SymmetricTensor< 4, dim, Number > outer_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2)
static ::ExceptionBase & ExcNonExistentSubdomain(int arg1)
static Tensor< 1, dim > d_linear_shape_function_gradient(const Point< dim > &xi, const unsigned int i)
numbers::NumberTraits< Number >::real_type norm() const
std::array< unsigned int, 4 > binary_type
SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
static ::ExceptionBase & ExcTriangulationHasBeenRefined()
static ::ExceptionBase & ExcInvalidNumberOfPartitions(int arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
unsigned int subdomain_id
#define DeclException1(Exception1, type1, outsequence)
T sum(const T &t, const MPI_Comm &mpi_communicator)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
#define DeclException0(Exception0)
void copy_material_to_manifold_id(Triangulation< dim, spacedim > &tria, const bool compute_face_ids=false)
size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
unsigned int global_dof_index
__global__ void set(Number *val, const Number s, const size_type N)
#define AssertThrowMPI(error_code)
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
T unpack(const std::vector< char > &buffer, const bool allow_compression=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)
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcVertexNotUsed(unsigned int arg1)
static ::ExceptionBase & ExcPointNotFound(Point< N > arg1)
static ::ExceptionBase & ExcInternalError()