Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
grid_tools.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2001 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_grid_tools_h
17# define dealii_grid_tools_h
18
19
20# include <deal.II/base/config.h>
21
25
27
29
31
32# include <deal.II/fe/fe_values.h>
33# include <deal.II/fe/mapping.h>
35
37# include <deal.II/grid/tria.h>
40
42
48
50
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>
57
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>
63# endif
65
66# include <bitset>
67# include <list>
68# include <set>
69
71
72// Forward declarations
73# ifndef DOXYGEN
74namespace parallel
75{
76 namespace distributed
77 {
78 template <int, int>
79 class Triangulation;
80 }
81} // namespace parallel
82
83namespace hp
84{
85 template <int, int>
86 class MappingCollection;
87}
88
89class SparsityPattern;
90# endif
91
92namespace internal
93{
94 template <int dim, int spacedim, class MeshType>
96 {
97 public:
98# ifndef _MSC_VER
99 using type = typename MeshType::active_cell_iterator;
100# else
102# endif
103 };
104
105# ifdef _MSC_VER
106 template <int dim, int spacedim>
107 class ActiveCellIterator<dim, spacedim, ::DoFHandler<dim, spacedim>>
108 {
109 public:
110 using type =
112 };
113# endif
114
115
116# ifdef _MSC_VER
117 template <int dim, int spacedim>
118 class ActiveCellIterator<dim, spacedim, ::hp::DoFHandler<dim, spacedim>>
119 {
120 public:
121 using type =
123 };
124# endif
125} // namespace internal
126
135namespace GridTools
136{
137 template <int dim, int spacedim>
138 class Cache;
139
144
151 template <int dim, int spacedim>
152 double
154
181 template <int dim, int spacedim>
182 double
184 const Mapping<dim, spacedim> & mapping =
185 (ReferenceCells::get_hypercube<dim>()
186 .template get_default_linear_mapping<dim, spacedim>()));
187
198 template <int dim, int spacedim>
199 double
202 const Mapping<dim, spacedim> & mapping =
203 (ReferenceCells::get_hypercube<dim>()
204 .template get_default_linear_mapping<dim, spacedim>()));
205
216 template <int dim, int spacedim>
217 double
220 const Mapping<dim, spacedim> & mapping =
221 (ReferenceCells::get_hypercube<dim>()
222 .template get_default_linear_mapping<dim, spacedim>()));
223
235 template <int dim>
238 const std::vector<Point<dim>> &all_vertices,
240
257 template <int dim>
258 double
259 cell_measure(const std::vector<Point<dim>> & all_vertices,
261
284 template <int dim, int spacedim>
285 std::pair<DerivativeForm<1, dim, spacedim>, Tensor<1, spacedim>>
287
314 template <int dim>
318 const Quadrature<dim> & quadrature);
319
327 template <int dim>
328 double
331 const Quadrature<dim> & quadrature);
332
346 template <int dim, int spacedim>
349
367 template <typename Iterator>
370 const Iterator & object,
372
385 template <int dim, int spacedim>
386 std::
387 tuple<std::vector<Point<spacedim>>, std::vector<CellData<dim>>, SubCellData>
389
395
412 template <int dim, int spacedim>
413 void
415 std::vector<CellData<dim>> & cells,
416 SubCellData & subcelldata);
417
435 template <int dim, int spacedim>
436 void
437 delete_duplicated_vertices(std::vector<Point<spacedim>> &all_vertices,
438 std::vector<CellData<dim>> & cells,
439 SubCellData & subcelldata,
440 std::vector<unsigned int> & considered_vertices,
441 const double tol = 1e-12);
442
461 template <int dim, int spacedim>
462 void
464 const std::vector<Point<spacedim>> &all_vertices,
465 std::vector<CellData<dim>> & cells);
466
476 template <int dim>
477 void
478 consistently_order_cells(std::vector<CellData<dim>> &cells);
479
485
539 template <int dim, typename Transformation, int spacedim>
540 void
541 transform(const Transformation & transformation,
543
549 template <int dim, int spacedim>
550 void
551 shift(const Tensor<1, spacedim> & shift_vector,
553
554
564 template <int dim>
565 void
567
580 template <int dim>
581 void
582 rotate(const double angle,
583 const unsigned int axis,
585
643 template <int dim>
644 void
645 laplace_transform(const std::map<unsigned int, Point<dim>> &new_points,
646 Triangulation<dim> & tria,
647 const Function<dim, double> *coefficient = nullptr,
648 const bool solve_for_absolute_positions = false);
649
655 template <int dim, int spacedim>
656 std::map<unsigned int, Point<spacedim>>
658
666 template <int dim, int spacedim>
667 void
668 scale(const double scaling_factor,
670
685 template <int dim, int spacedim>
686 void
688 const double factor,
690 const bool keep_boundary = true,
691 const unsigned int seed = boost::random::mt19937::default_seed);
692
726 template <int dim, int spacedim>
727 void
729 const bool isotropic = false,
730 const unsigned int max_iterations = 100);
731
756 template <int dim, int spacedim>
757 void
759 const double max_ratio = 1.6180339887,
760 const unsigned int max_iterations = 5);
761
851 template <int dim, int spacedim>
852 void
854 const double limit_angle_fraction = .75);
855
861
916 template <int dim, int spacedim>
917# ifndef DOXYGEN
918 std::tuple<
919 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
920 std::vector<std::vector<Point<dim>>>,
921 std::vector<std::vector<unsigned int>>>
922# else
923 return_type
924# endif
926 const Cache<dim, spacedim> & cache,
927 const std::vector<Point<spacedim>> &points,
929 &cell_hint =
931
965 template <int dim, int spacedim>
966# ifndef DOXYGEN
967 std::tuple<
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>>
972# else
973 return_type
974# endif
976 const Cache<dim, spacedim> & cache,
977 const std::vector<Point<spacedim>> &points,
979 &cell_hint =
981
1051 template <int dim, int spacedim>
1052# ifndef DOXYGEN
1053 std::tuple<
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>>>
1059# else
1060 return_type
1061# endif
1063 const GridTools::Cache<dim, spacedim> & cache,
1064 const std::vector<Point<spacedim>> & local_points,
1065 const std::vector<std::vector<BoundingBox<spacedim>>> &global_bboxes,
1066 const double tolerance = 1e-10);
1067
1068 namespace internal
1069 {
1084 template <int dim, int spacedim>
1086 {
1093 std::vector<std::tuple<std::pair<int, int>,
1094 unsigned int,
1095 unsigned int,
1096 Point<dim>,
1098 unsigned int>>
1100
1104 std::vector<unsigned int> send_ranks;
1105
1111 std::vector<unsigned int> send_ptrs;
1112
1123 std::vector<std::tuple<unsigned int, unsigned int, unsigned int>>
1125
1129 std::vector<unsigned int> recv_ranks;
1130
1136 std::vector<unsigned int> recv_ptrs;
1137 };
1138
1148 template <int dim, int spacedim>
1151 const GridTools::Cache<dim, spacedim> & cache,
1152 const std::vector<Point<spacedim>> & points,
1153 const std::vector<std::vector<BoundingBox<spacedim>>> &global_bboxes,
1154 const double tolerance,
1155 const bool perform_handshake,
1156 const bool enforce_unique_mapping = false);
1157
1158 } // namespace internal
1159
1196 template <int dim, int spacedim>
1197 std::map<unsigned int, Point<spacedim>>
1199 const Triangulation<dim, spacedim> &container,
1200 const Mapping<dim, spacedim> & mapping =
1201 (ReferenceCells::get_hypercube<dim>()
1202 .template get_default_linear_mapping<dim, spacedim>()));
1203
1213 template <int spacedim>
1214 unsigned int
1215 find_closest_vertex(const std::map<unsigned int, Point<spacedim>> &vertices,
1216 const Point<spacedim> & p);
1217
1241 template <int dim, template <int, int> class MeshType, int spacedim>
1242 unsigned int
1243 find_closest_vertex(const MeshType<dim, spacedim> &mesh,
1244 const Point<spacedim> & p,
1245 const std::vector<bool> & marked_vertices = {});
1246
1270 template <int dim, template <int, int> class MeshType, int spacedim>
1271 unsigned int
1273 const MeshType<dim, spacedim> &mesh,
1274 const Point<spacedim> & p,
1275 const std::vector<bool> & marked_vertices = {});
1276
1277
1297 template <int dim, template <int, int> class MeshType, int spacedim>
1298# ifndef _MSC_VER
1299 std::vector<typename MeshType<dim, spacedim>::active_cell_iterator>
1300# else
1301 std::vector<
1302 typename ::internal::
1303 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type>
1304# endif
1305 find_cells_adjacent_to_vertex(const MeshType<dim, spacedim> &container,
1306 const unsigned int vertex_index);
1307
1370 template <int dim, template <int, int> class MeshType, int spacedim>
1371# ifndef _MSC_VER
1372 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator, Point<dim>>
1373# else
1374 std::pair<typename ::internal::
1375 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1376 Point<dim>>
1377# endif
1379 const MeshType<dim, spacedim> &mesh,
1380 const Point<spacedim> & p,
1381 const std::vector<bool> &marked_vertices = {},
1382 const double tolerance = 1.e-10);
1383
1391 template <int dim, template <int, int> class MeshType, int spacedim>
1392# ifndef _MSC_VER
1393 typename MeshType<dim, spacedim>::active_cell_iterator
1394# else
1395 typename ::internal::
1396 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type
1397# endif
1398 find_active_cell_around_point(const MeshType<dim, spacedim> &mesh,
1399 const Point<spacedim> & p,
1400 const std::vector<bool> &marked_vertices = {},
1401 const double tolerance = 1.e-10);
1402
1409 template <int dim, int spacedim>
1410 std::pair<typename DoFHandler<dim, spacedim>::active_cell_iterator,
1411 Point<dim>>
1414 const DoFHandler<dim, spacedim> & mesh,
1415 const Point<spacedim> & p,
1416 const double tolerance = 1.e-10);
1417
1469 template <int dim, int spacedim>
1470 std::pair<typename Triangulation<dim, spacedim>::active_cell_iterator,
1471 Point<dim>>
1473 const Cache<dim, spacedim> &cache,
1474 const Point<spacedim> & p,
1477 const std::vector<bool> &marked_vertices = {},
1478 const double tolerance = 1.e-10);
1479
1493 template <int dim, template <int, int> class MeshType, int spacedim>
1494# ifndef _MSC_VER
1495 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator, Point<dim>>
1496# else
1497 std::pair<typename ::internal::
1498 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1499 Point<dim>>
1500# endif
1502 const Mapping<dim, spacedim> & mapping,
1503 const MeshType<dim, spacedim> &mesh,
1504 const Point<spacedim> & p,
1505 const std::vector<
1506 std::set<typename MeshType<dim, spacedim>::active_cell_iterator>>
1508 const std::vector<std::vector<Tensor<1, spacedim>>> &vertex_to_cell_centers,
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 = {},
1512 const RTree<std::pair<Point<spacedim>, unsigned int>> &used_vertices_rtree =
1513 RTree<std::pair<Point<spacedim>, unsigned int>>{},
1514 const double tolerance = 1.e-10);
1515
1537 template <int dim, template <int, int> class MeshType, int spacedim>
1538# ifndef _MSC_VER
1539 std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
1540 Point<dim>>>
1541# else
1542 std::vector<std::pair<
1543 typename ::internal::
1544 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1545 Point<dim>>>
1546# endif
1548 const Mapping<dim, spacedim> & mapping,
1549 const MeshType<dim, spacedim> &mesh,
1550 const Point<spacedim> & p,
1551 const double tolerance,
1552 const std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
1553 Point<dim>> & first_cell);
1554
1561 template <int dim, template <int, int> class MeshType, int spacedim>
1562# ifndef _MSC_VER
1563 std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
1564 Point<dim>>>
1565# else
1566 std::vector<std::pair<
1567 typename ::internal::
1568 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1569 Point<dim>>>
1570# endif
1572 const Mapping<dim, spacedim> & mapping,
1573 const MeshType<dim, spacedim> &mesh,
1574 const Point<spacedim> & p,
1575 const double tolerance = 1e-10,
1576 const std::vector<bool> & marked_vertices = {});
1577
1599 template <class MeshType>
1600 std::vector<typename MeshType::active_cell_iterator>
1601 get_active_child_cells(const typename MeshType::cell_iterator &cell);
1602
1627 template <class MeshType>
1628 void
1630 const typename MeshType::active_cell_iterator & cell,
1631 std::vector<typename MeshType::active_cell_iterator> &active_neighbors);
1632
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 &)>
1687 &predicate);
1688
1689
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 &)>
1702 & predicate,
1703 const unsigned int level);
1704
1705
1718 template <class MeshType>
1719 std::vector<typename MeshType::active_cell_iterator>
1720 compute_ghost_cell_halo_layer(const MeshType &mesh);
1721
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 &)>
1775 & predicate,
1776 const double layer_thickness);
1777
1800 template <class MeshType>
1801 std::vector<typename MeshType::active_cell_iterator>
1802 compute_ghost_cell_layer_within_distance(const MeshType &mesh,
1803 const double layer_thickness);
1804
1820 template <class MeshType>
1821 std::pair<Point<MeshType::space_dimension>, Point<MeshType::space_dimension>>
1823 const MeshType &mesh,
1824 const std::function<bool(const typename MeshType::active_cell_iterator &)>
1825 &predicate);
1826
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 &)>
1884 & predicate,
1885 const unsigned int refinement_level = 0,
1886 const bool allow_merge = false,
1887 const unsigned int max_boxes = numbers::invalid_unsigned_int);
1888
1916 template <int spacedim>
1917# ifndef DOXYGEN
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>>>
1921# else
1922 return_type
1923# endif
1925 const std::vector<std::vector<BoundingBox<spacedim>>> &global_bboxes,
1926 const std::vector<Point<spacedim>> & points);
1927
1928
1963 template <int spacedim>
1964# ifndef DOXYGEN
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>>>
1968# else
1969 return_type
1970# endif
1972 const RTree<std::pair<BoundingBox<spacedim>, unsigned int>> &covering_rtree,
1973 const std::vector<Point<spacedim>> & points);
1974
1975
1984 template <int dim, int spacedim>
1985 std::vector<
1986 std::set<typename Triangulation<dim, spacedim>::active_cell_iterator>>
1988
2001 template <int dim, int spacedim>
2002 std::vector<std::vector<Tensor<1, spacedim>>>
2004 const Triangulation<dim, spacedim> &mesh,
2005 const std::vector<
2007 &vertex_to_cells);
2008
2009
2017 template <int dim, int spacedim>
2018 unsigned int
2021 const Point<spacedim> & position,
2022 const Mapping<dim, spacedim> & mapping =
2023 (ReferenceCells::get_hypercube<dim>()
2024 .template get_default_linear_mapping<dim, spacedim>()));
2025
2037 template <int dim, int spacedim>
2038 std::map<unsigned int, types::global_vertex_index>
2041
2053 template <int dim, int spacedim>
2054 std::pair<unsigned int, double>
2057
2063
2072 template <int dim, int spacedim>
2073 void
2076 DynamicSparsityPattern & connectivity);
2077
2086 template <int dim, int spacedim>
2087 void
2090 DynamicSparsityPattern & connectivity);
2091
2100 template <int dim, int spacedim>
2101 void
2104 const unsigned int level,
2105 DynamicSparsityPattern & connectivity);
2106
2127 template <int dim, int spacedim>
2128 void
2129 partition_triangulation(const unsigned int n_partitions,
2131 const SparsityTools::Partitioner partitioner =
2133
2144 template <int dim, int spacedim>
2145 void
2146 partition_triangulation(const unsigned int n_partitions,
2147 const std::vector<unsigned int> &cell_weights,
2149 const SparsityTools::Partitioner partitioner =
2151
2197 template <int dim, int spacedim>
2198 void
2199 partition_triangulation(const unsigned int n_partitions,
2200 const SparsityPattern & cell_connection_graph,
2202 const SparsityTools::Partitioner partitioner =
2204
2215 template <int dim, int spacedim>
2216 void
2217 partition_triangulation(const unsigned int n_partitions,
2218 const std::vector<unsigned int> &cell_weights,
2219 const SparsityPattern & cell_connection_graph,
2221 const SparsityTools::Partitioner partitioner =
2223
2238 template <int dim, int spacedim>
2239 void
2240 partition_triangulation_zorder(const unsigned int n_partitions,
2242 const bool group_siblings = true);
2243
2255 template <int dim, int spacedim>
2256 void
2258
2266 template <int dim, int spacedim>
2267 std::vector<types::subdomain_id>
2269 const std::vector<CellId> & cell_ids);
2270
2281 template <int dim, int spacedim>
2282 void
2284 std::vector<types::subdomain_id> & subdomain);
2285
2300 template <int dim, int spacedim>
2301 unsigned int
2304 const types::subdomain_id subdomain);
2305
2335 template <int dim, int spacedim>
2336 std::vector<bool>
2338
2344
2378 template <typename MeshType>
2379 std::list<std::pair<typename MeshType::cell_iterator,
2380 typename MeshType::cell_iterator>>
2381 get_finest_common_cells(const MeshType &mesh_1, const MeshType &mesh_2);
2382
2392 template <int dim, int spacedim>
2393 bool
2395 const Triangulation<dim, spacedim> &mesh_2);
2396
2406 template <typename MeshType>
2407 bool
2408 have_same_coarse_mesh(const MeshType &mesh_1, const MeshType &mesh_2);
2409
2415
2431 template <int dim, int spacedim>
2435 & distorted_cells,
2437
2438
2439
2448
2449
2487 template <class MeshType>
2488 std::vector<typename MeshType::active_cell_iterator>
2489 get_patch_around_cell(const typename MeshType::active_cell_iterator &cell);
2490
2491
2513 template <class Container>
2514 std::vector<typename Container::cell_iterator>
2516 const std::vector<typename Container::active_cell_iterator> &patch_cells);
2517
2584 template <class Container>
2585 void
2587 const std::vector<typename Container::active_cell_iterator> &patch,
2589 &local_triangulation,
2590 std::map<
2591 typename Triangulation<Container::dimension,
2592 Container::space_dimension>::active_cell_iterator,
2593 typename Container::active_cell_iterator> &patch_to_global_tria_map);
2594
2626 template <int dim, int spacedim>
2627 std::map<
2629 std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>>
2631
2632
2639
2645 template <typename CellIterator>
2647 {
2651 CellIterator cell[2];
2652
2657 unsigned int face_idx[2];
2658
2664 std::bitset<3> orientation;
2665
2679 };
2680
2681
2745 template <typename FaceIterator>
2747 std::bitset<3> & orientation,
2748 const FaceIterator & face1,
2749 const FaceIterator & face2,
2750 const int direction,
2754
2755
2759 template <typename FaceIterator>
2760 bool
2762 const FaceIterator & face1,
2763 const FaceIterator & face2,
2764 const int direction,
2768
2769
2826 template <typename MeshType>
2827 void
2829 const MeshType & mesh,
2830 const types::boundary_id b_id1,
2831 const types::boundary_id b_id2,
2832 const int direction,
2834 & matched_pairs,
2838
2839
2862 template <typename MeshType>
2863 void
2865 const MeshType & mesh,
2866 const types::boundary_id b_id,
2867 const int direction,
2869 & matched_pairs,
2870 const ::Tensor<1, MeshType::space_dimension> &offset =
2873
2879
2900 template <int dim, int spacedim>
2901 void
2903 const bool reset_boundary_ids = false);
2904
2926 template <int dim, int spacedim>
2927 void
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 = {});
2933
2963 template <int dim, int spacedim>
2964 void
2966 const bool compute_face_ids = false);
2967
2992 template <int dim, int spacedim>
2993 void
2996 const std::function<types::manifold_id(
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();
3001 else
3003 },
3004 bool overwrite_only_flat_manifold_ids = true);
3091 template <typename DataType, typename MeshType>
3092 void
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 &)>
3100 &cell_filter =
3102
3113 template <typename DataType, typename MeshType>
3114 void
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 &)> &
3123 true});
3124
3125 /* Exchange with all processors of the MPI communicator @p mpi_communicator the vector of bounding
3126 * boxes @p local_bboxes.
3127 *
3128 * This function is meant to exchange bounding boxes describing the locally
3129 * owned cells in a distributed triangulation obtained with the function
3130 * GridTools::compute_mesh_predicate_bounding_box .
3131 *
3132 * The output vector's size is the number of processes of the MPI
3133 * communicator:
3134 * its i-th entry contains the vector @p local_bboxes of the i-th process.
3135 */
3136 template <int spacedim>
3137 std::vector<std::vector<BoundingBox<spacedim>>>
3139 const std::vector<BoundingBox<spacedim>> &local_bboxes,
3140 const MPI_Comm & mpi_communicator);
3141
3174 template <int spacedim>
3177 const std::vector<BoundingBox<spacedim>> &local_description,
3178 const MPI_Comm & mpi_communicator);
3179
3197 template <int dim, int spacedim>
3198 void
3200 const Triangulation<dim, spacedim> & tria,
3201 std::map<unsigned int, std::vector<unsigned int>> &coinciding_vertex_groups,
3202 std::map<unsigned int, unsigned int> &vertex_to_coinciding_vertex_group);
3203
3210 template <int dim, int spacedim>
3211 std::map<unsigned int, std::set<::types::subdomain_id>>
3213 const Triangulation<dim, spacedim> &tria);
3214
3227 template <int dim, typename T>
3229 {
3233 std::vector<CellId> cell_ids;
3234
3238 std::vector<T> data;
3239
3248 template <class Archive>
3249 void
3250 save(Archive &ar, const unsigned int version) const;
3251
3258 template <class Archive>
3259 void
3260 load(Archive &ar, const unsigned int version);
3261
3262# ifdef DOXYGEN
3268 template <class Archive>
3269 void
3270 serialize(Archive &archive, const unsigned int version);
3271# else
3272 // This macro defines the serialize() method that is compatible with
3273 // the templated save() and load() method that have been implemented.
3274 BOOST_SERIALIZATION_SPLIT_MEMBER()
3275# endif
3276 };
3277
3278
3279
3297 template <int dim, typename VectorType>
3299 {
3300 public:
3304 using value_type = typename VectorType::value_type;
3305
3310 const FiniteElement<dim, dim> &fe,
3311 const unsigned int n_subdivisions = 1,
3312 const double tolerance = 1e-10);
3313
3318 void
3319 process(const DoFHandler<dim> & background_dof_handler,
3320 const VectorType & ls_vector,
3321 const double iso_level,
3322 std::vector<Point<dim>> & vertices,
3323 std::vector<CellData<dim - 1>> &cells) const;
3324
3331 void
3333 const VectorType & ls_vector,
3334 const double iso_level,
3335 std::vector<Point<dim>> & vertices,
3336 std::vector<CellData<dim - 1>> &cells) const;
3337
3338 private:
3343 static Quadrature<dim>
3344 create_quadrature_rule(const unsigned int n_subdivisions);
3345
3349 void
3350 process_cell(std::vector<value_type> & ls_values,
3351 const std::vector<Point<dim>> & points,
3352 const double iso_level,
3353 std::vector<Point<dim>> & vertices,
3354 std::vector<CellData<dim - 1>> &cells) const;
3355
3362 void
3363 process_sub_cell(const std::vector<value_type> & ls_values,
3364 const std::vector<Point<2>> & points,
3365 const std::vector<unsigned int> mask,
3366 const double iso_level,
3367 std::vector<Point<2>> & vertices,
3368 std::vector<CellData<1>> & cells) const;
3369
3373 void
3374 process_sub_cell(const std::vector<value_type> & ls_values,
3375 const std::vector<Point<3>> & points,
3376 const std::vector<unsigned int> mask,
3377 const double iso_level,
3378 std::vector<Point<3>> & vertices,
3379 std::vector<CellData<2>> & cells) const;
3380
3385 const unsigned int n_subdivisions;
3386
3391 const double tolerance;
3392
3398 };
3399
3400
3401
3406
3411 int,
3412 << "The number of partitions you gave is " << arg1
3413 << ", but must be greater than zero.");
3418 int,
3419 << "The subdomain id " << arg1
3420 << " has no cells associated with it.");
3425
3430 double,
3431 << "The scaling factor must be positive, but it is " << arg1
3432 << ".");
3433
3438 unsigned int,
3439 << "The given vertex with index " << arg1
3440 << " is not used in the given triangulation.");
3441
3444} /*namespace GridTools*/
3445
3446
3457 "The edges of the mesh are not consistently orientable.");
3458
3459
3460
3461/* ----------------- Template function --------------- */
3462
3463# ifndef DOXYGEN
3464
3465namespace GridTools
3466{
3467 template <int dim>
3468 double
3470 const std::vector<Point<dim>> &all_vertices,
3471 const unsigned int (&indices)[GeometryInfo<dim>::vertices_per_cell])
3472 {
3473 // We forward call to the ArrayView version:
3476 return cell_measure(all_vertices, view);
3477 }
3478
3479
3480
3481 // This specialization is defined here so that the general template in the
3482 // source file doesn't need to have further 1D overloads for the internal
3483 // functions it calls.
3484 template <>
3488 {
3489 return {};
3490 }
3491
3492
3493
3494 template <int dim, typename Predicate, int spacedim>
3495 void
3496 transform(const Predicate & predicate,
3498 {
3499 std::vector<bool> treated_vertices(triangulation.n_vertices(), false);
3500
3501 // loop over all active cells, and
3502 // transform those vertices that
3503 // have not yet been touched. note
3504 // that we get to all vertices in
3505 // the triangulation by only
3506 // visiting the active cells.
3508 cell = triangulation.begin_active(),
3509 endc = triangulation.end();
3510 for (; cell != endc; ++cell)
3511 for (const unsigned int v : cell->vertex_indices())
3512 if (treated_vertices[cell->vertex_index(v)] == false)
3513 {
3514 // transform this vertex
3515 cell->vertex(v) = predicate(cell->vertex(v));
3516 // and mark it as treated
3517 treated_vertices[cell->vertex_index(v)] = true;
3518 };
3519
3520
3521 // now fix any vertices on hanging nodes so that we don't create any holes
3522 if (dim == 2)
3523 {
3525 cell = triangulation.begin_active(),
3526 endc = triangulation.end();
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())
3531 {
3532 Assert(cell->reference_cell() ==
3533 ReferenceCells::get_hypercube<dim>(),
3535
3536 // this line has children
3537 cell->face(face)->child(0)->vertex(1) =
3538 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1)) /
3539 2;
3540 }
3541 }
3542 else if (dim == 3)
3543 {
3545 cell = triangulation.begin_active(),
3546 endc = triangulation.end();
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())
3551 {
3552 Assert(cell->reference_cell() ==
3553 ReferenceCells::get_hypercube<dim>(),
3555
3556 // this face has hanging nodes
3557 cell->face(face)->child(0)->vertex(1) =
3558 (cell->face(face)->vertex(0) + cell->face(face)->vertex(1)) /
3559 2.0;
3560 cell->face(face)->child(0)->vertex(2) =
3561 (cell->face(face)->vertex(0) + cell->face(face)->vertex(2)) /
3562 2.0;
3563 cell->face(face)->child(1)->vertex(3) =
3564 (cell->face(face)->vertex(1) + cell->face(face)->vertex(3)) /
3565 2.0;
3566 cell->face(face)->child(2)->vertex(3) =
3567 (cell->face(face)->vertex(2) + cell->face(face)->vertex(3)) /
3568 2.0;
3569
3570 // center of the face
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)) /
3574 4.0;
3575 }
3576 }
3577
3578 // Make sure FEValues notices that the mesh has changed
3579 triangulation.signals.mesh_movement();
3580 }
3581
3582
3583
3584 template <class MeshType>
3585 std::vector<typename MeshType::active_cell_iterator>
3586 get_active_child_cells(const typename MeshType::cell_iterator &cell)
3587 {
3588 std::vector<typename MeshType::active_cell_iterator> child_cells;
3589
3590 if (cell->has_children())
3591 {
3592 for (unsigned int child = 0; child < cell->n_children(); ++child)
3593 if (cell->child(child)->has_children())
3594 {
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(),
3598 children.begin(),
3599 children.end());
3600 }
3601 else
3602 child_cells.push_back(cell->child(child));
3603 }
3604
3605 return child_cells;
3606 }
3607
3608
3609
3610 template <class MeshType>
3611 void
3613 const typename MeshType::active_cell_iterator & cell,
3614 std::vector<typename MeshType::active_cell_iterator> &active_neighbors)
3615 {
3616 active_neighbors.clear();
3617 for (const unsigned int n : cell->face_indices())
3618 if (!cell->at_boundary(n))
3619 {
3620 if (MeshType::dimension == 1)
3621 {
3622 // check children of neighbor. note
3623 // that in 1d children of the neighbor
3624 // may be further refined. In 1d the
3625 // case is simple since we know what
3626 // children bound to the present cell
3627 typename MeshType::cell_iterator neighbor_child =
3628 cell->neighbor(n);
3629 if (!neighbor_child->is_active())
3630 {
3631 while (neighbor_child->has_children())
3632 neighbor_child = neighbor_child->child(n == 0 ? 1 : 0);
3633
3634 Assert(neighbor_child->neighbor(n == 0 ? 1 : 0) == cell,
3636 }
3637 active_neighbors.push_back(neighbor_child);
3638 }
3639 else
3640 {
3641 if (cell->face(n)->has_children())
3642 // this neighbor has children. find
3643 // out which border to the present
3644 // cell
3645 for (unsigned int c = 0;
3646 c < cell->face(n)->n_active_descendants();
3647 ++c)
3648 active_neighbors.push_back(
3649 cell->neighbor_child_on_subface(n, c));
3650 else
3651 {
3652 // the neighbor must be active
3653 // himself
3654 Assert(cell->neighbor(n)->is_active(), ExcInternalError());
3655 active_neighbors.push_back(cell->neighbor(n));
3656 }
3657 }
3658 }
3659 }
3660
3661
3662
3663 namespace internal
3664 {
3665 namespace ProjectToObject
3666 {
3679 struct CrossDerivative
3680 {
3681 const unsigned int direction_0;
3682 const unsigned int direction_1;
3683
3684 CrossDerivative(const unsigned int d0, const unsigned int d1);
3685 };
3686
3687 inline CrossDerivative::CrossDerivative(const unsigned int d0,
3688 const unsigned int d1)
3689 : direction_0(d0)
3690 , direction_1(d1)
3691 {}
3692
3693
3694
3699 template <typename F>
3700 inline auto
3701 centered_first_difference(const double center,
3702 const double step,
3703 const F &f) -> decltype(f(center) - f(center))
3704 {
3705 return (f(center + step) - f(center - step)) / (2.0 * step);
3706 }
3707
3708
3709
3714 template <typename F>
3715 inline auto
3716 centered_second_difference(const double center,
3717 const double step,
3718 const F &f) -> decltype(f(center) - f(center))
3719 {
3720 return (f(center + step) - 2.0 * f(center) + f(center - step)) /
3721 (step * step);
3722 }
3723
3724
3725
3735 template <int structdim, typename F>
3736 inline auto
3737 cross_stencil(
3738 const CrossDerivative cross_derivative,
3740 const double step,
3741 const F &f) -> decltype(f(center) - f(center))
3742 {
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)) /
3749 step;
3750 }
3751
3752
3753
3760 template <int spacedim, int structdim, typename F>
3761 inline double
3762 gradient_entry(
3763 const unsigned int row_n,
3764 const unsigned int dependent_direction,
3765 const Point<spacedim> &p0,
3767 const double step,
3768 const F & f)
3769 {
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,
3777 ExcMessage(
3778 "We cannot differentiate with respect to the variable "
3779 "that is assumed to be dependent."));
3780
3781 const Point<spacedim> manifold_point = f(center);
3782 const Tensor<1, spacedim> stencil_value = cross_stencil<structdim>(
3783 {row_n, dependent_direction}, center, step, f);
3784 double entry = 0.0;
3785 for (unsigned int dim_n = 0; dim_n < spacedim; ++dim_n)
3786 entry +=
3787 -2.0 * (p0[dim_n] - manifold_point[dim_n]) * stencil_value[dim_n];
3788 return entry;
3789 }
3790
3796 template <typename Iterator, int spacedim, int structdim>
3798 project_to_d_linear_object(const Iterator & object,
3799 const Point<spacedim> &trial_point)
3800 {
3801 // let's look at this for simplicity for a quad (structdim==2) in a
3802 // space with spacedim>2 (notate trial_point by y): all points on the
3803 // surface are given by
3804 // x(\xi) = sum_i v_i phi_x(\xi)
3805 // where v_i are the vertices of the quad, and \xi=(\xi_1,\xi_2) are the
3806 // reference coordinates of the quad. so what we are trying to do is
3807 // find a point x on the surface that is closest to the point y. there
3808 // are different ways to solve this problem, but in the end it's a
3809 // nonlinear problem and we have to find reference coordinates \xi so
3810 // that J(\xi) = 1/2 || x(\xi)-y ||^2 is minimal. x(\xi) is a function
3811 // that is structdim-linear in \xi, so J(\xi) is a polynomial of degree
3812 // 2*structdim that we'd like to minimize. unless structdim==1, we'll
3813 // have to use a Newton method to find the answer. This leads to the
3814 // following formulation of Newton steps:
3815 //
3816 // Given \xi_k, find \delta\xi_k so that
3817 // H_k \delta\xi_k = - F_k
3818 // where H_k is an approximation to the second derivatives of J at
3819 // \xi_k, and F_k is the first derivative of J. We'll iterate this a
3820 // number of times until the right hand side is small enough. As a
3821 // stopping criterion, we terminate if ||\delta\xi||<eps.
3822 //
3823 // As for the Hessian, the best choice would be
3824 // H_k = J''(\xi_k)
3825 // but we'll opt for the simpler Gauss-Newton form
3826 // H_k = A^T A
3827 // i.e.
3828 // (H_k)_{nm} = \sum_{i,j} v_i*v_j *
3829 // \partial_n phi_i *
3830 // \partial_m phi_j
3831 // we start at xi=(0.5, 0.5).
3833 for (unsigned int d = 0; d < structdim; ++d)
3834 xi[d] = 0.5;
3835
3836 Point<spacedim> x_k;
3837 for (const unsigned int i : GeometryInfo<structdim>::vertex_indices())
3838 x_k += object->vertex(i) *
3840
3841 do
3842 {
3844 for (const unsigned int i :
3846 F_k +=
3847 (x_k - trial_point) * object->vertex(i) *
3849 i);
3850
3852 for (const unsigned int i :
3854 for (const unsigned int j :
3856 {
3859 xi, i),
3861 xi, j));
3862 H_k += (object->vertex(i) * object->vertex(j)) * tmp;
3863 }
3864
3865 const Tensor<1, structdim> delta_xi = -invert(H_k) * F_k;
3866 xi += delta_xi;
3867
3868 x_k = Point<spacedim>();
3869 for (const unsigned int i :
3871 x_k += object->vertex(i) *
3873
3874 if (delta_xi.norm() < 1e-7)
3875 break;
3876 }
3877 while (true);
3878
3879 return x_k;
3880 }
3881 } // namespace ProjectToObject
3882 } // namespace internal
3883
3884
3885
3886 namespace internal
3887 {
3888 // We hit an internal compiler error in ICC 15 if we define this as a lambda
3889 // inside the project_to_object function below.
3890 template <int structdim>
3891 inline bool
3892 weights_are_ok(
3894 {
3895 // clang has trouble figuring out structdim here, so define it
3896 // again:
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)
3902 {
3903 copied_weights[i] = v[i];
3904 if (v[i] < 0.0 || v[i] > 1.0)
3905 return false;
3906 }
3907
3908 // check the sum: try to avoid some roundoff errors by summing in order
3909 std::sort(copied_weights.begin(), copied_weights.end());
3910 const double sum =
3911 std::accumulate(copied_weights.begin(), copied_weights.end(), 0.0);
3912 return std::abs(sum - 1.0) < 1e-10; // same tolerance used in manifold.cc
3913 }
3914 } // namespace internal
3915
3916 template <typename Iterator>
3919 const Iterator & object,
3921 {
3922 const int spacedim = Iterator::AccessorType::space_dimension;
3923 const int structdim = Iterator::AccessorType::structure_dimension;
3924
3925 Point<spacedim> projected_point = trial_point;
3926
3927 if (structdim >= spacedim)
3928 return projected_point;
3929 else if (structdim == 1 || structdim == 2)
3930 {
3931 using namespace internal::ProjectToObject;
3932 // Try to use the special flat algorithm for quads (this is better
3933 // than the general algorithm in 3D). This does not take into account
3934 // whether projected_point is outside the quad, but we optimize along
3935 // lines below anyway:
3936 const int dim = Iterator::AccessorType::dimension;
3937 const Manifold<dim, spacedim> &manifold = object->get_manifold();
3938 if (structdim == 2 && dynamic_cast<const FlatManifold<dim, spacedim> *>(
3939 &manifold) != nullptr)
3940 {
3941 projected_point =
3942 project_to_d_linear_object<Iterator, spacedim, structdim>(
3943 object, trial_point);
3944 }
3945 else
3946 {
3947 // We want to find a point on the convex hull (defined by the
3948 // vertices of the object and the manifold description) that is
3949 // relatively close to the trial point. This has a few issues:
3950 //
3951 // 1. For a general convex hull we are not guaranteed that a unique
3952 // minimum exists.
3953 // 2. The independent variables in the optimization process are the
3954 // weights given to Manifold::get_new_point, which must sum to 1,
3955 // so we cannot use standard finite differences to approximate a
3956 // gradient.
3957 //
3958 // There is not much we can do about 1., but for 2. we can derive
3959 // finite difference stencils that work on a structdim-dimensional
3960 // simplex and rewrite the optimization problem to use those
3961 // instead. Consider the structdim 2 case and let
3962 //
3963 // F(c0, c1, c2, c3) = Manifold::get_new_point(vertices, {c0, c1,
3964 // c2, c3})
3965 //
3966 // where {c0, c1, c2, c3} are the weights for the four vertices on
3967 // the quadrilateral. We seek to minimize the Euclidean distance
3968 // between F(...) and trial_point. We can solve for c3 in terms of
3969 // the other weights and get, for one coordinate direction
3970 //
3971 // d/dc0 ((x0 - F(c0, c1, c2, 1 - c0 - c1 - c2))^2)
3972 // = -2(x0 - F(...)) (d/dc0 F(...) - d/dc3 F(...))
3973 //
3974 // where we substitute back in for c3 after taking the
3975 // derivative. We can compute a stencil for the cross derivative
3976 // d/dc0 - d/dc3: this is exactly what cross_stencil approximates
3977 // (and gradient_entry computes the sum over the independent
3978 // variables). Below, we somewhat arbitrarily pick the last
3979 // component as the dependent one.
3980 //
3981 // Since we can now calculate derivatives of the objective
3982 // function we can use gradient descent to minimize it.
3983 //
3984 // Of course, this is much simpler in the structdim = 1 case (we
3985 // could rewrite the projection as a 1D optimization problem), but
3986 // to reduce the potential for bugs we use the same code in both
3987 // cases.
3988 const double step_size = object->diameter() / 64.0;
3989
3990 constexpr unsigned int n_vertices_per_cell =
3992
3993 std::array<Point<spacedim>, n_vertices_per_cell> vertices;
3994 for (unsigned int vertex_n = 0; vertex_n < n_vertices_per_cell;
3995 ++vertex_n)
3996 vertices[vertex_n] = object->vertex(vertex_n);
3997
3998 auto get_point_from_weights =
3999 [&](const Tensor<1, n_vertices_per_cell> &weights)
4000 -> Point<spacedim> {
4001 return object->get_manifold().get_new_point(
4002 make_array_view(vertices.begin(), vertices.end()),
4003 make_array_view(weights.begin_raw(), weights.end_raw()));
4004 };
4005
4006 // pick the initial weights as (normalized) inverse distances from
4007 // the trial point:
4008 Tensor<1, n_vertices_per_cell> guess_weights;
4009 double guess_weights_sum = 0.0;
4010 for (unsigned int vertex_n = 0; vertex_n < n_vertices_per_cell;
4011 ++vertex_n)
4012 {
4013 const double distance =
4014 vertices[vertex_n].distance(trial_point);
4015 if (distance == 0.0)
4016 {
4017 guess_weights = 0.0;
4018 guess_weights[vertex_n] = 1.0;
4019 guess_weights_sum = 1.0;
4020 break;
4021 }
4022 else
4023 {
4024 guess_weights[vertex_n] = 1.0 / distance;
4025 guess_weights_sum += guess_weights[vertex_n];
4026 }
4027 }
4028 guess_weights /= guess_weights_sum;
4029 Assert(internal::weights_are_ok<structdim>(guess_weights),
4031
4032 // The optimization algorithm consists of two parts:
4033 //
4034 // 1. An outer loop where we apply the gradient descent algorithm.
4035 // 2. An inner loop where we do a line search to find the optimal
4036 // length of the step one should take in the gradient direction.
4037 //
4038 for (unsigned int outer_n = 0; outer_n < 40; ++outer_n)
4039 {
4040 const unsigned int dependent_direction =
4041 n_vertices_per_cell - 1;
4042 Tensor<1, n_vertices_per_cell> current_gradient;
4043 for (unsigned int row_n = 0; row_n < n_vertices_per_cell;
4044 ++row_n)
4045 {
4046 if (row_n != dependent_direction)
4047 {
4048 current_gradient[row_n] =
4049 gradient_entry<spacedim, structdim>(
4050 row_n,
4051 dependent_direction,
4052 trial_point,
4053 guess_weights,
4054 step_size,
4055 get_point_from_weights);
4056
4057 current_gradient[dependent_direction] -=
4058 current_gradient[row_n];
4059 }
4060 }
4061
4062 // We need to travel in the -gradient direction, as noted
4063 // above, but we may not want to take a full step in that
4064 // direction; instead, guess that we will go -0.5*gradient and
4065 // do quasi-Newton iteration to pick the best multiplier. The
4066 // goal is to find a scalar alpha such that
4067 //
4068 // F(x - alpha g)
4069 //
4070 // is minimized, where g is the gradient and F is the
4071 // objective function. To find the optimal value we find roots
4072 // of the derivative of the objective function with respect to
4073 // alpha by Newton iteration, where we approximate the first
4074 // and second derivatives of F(x - alpha g) with centered
4075 // finite differences.
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 *
4082 current_gradient))
4083 .norm_square();
4084 };
4085
4086 for (unsigned int inner_n = 0; inner_n < 10; ++inner_n)
4087 {
4088 const double update_numerator = centered_first_difference(
4089 gradient_weight,
4090 step_size,
4091 gradient_weight_objective_function);
4092 const double update_denominator =
4093 centered_second_difference(
4094 gradient_weight,
4095 step_size,
4096 gradient_weight_objective_function);
4097
4098 // avoid division by zero. Note that we limit the gradient
4099 // weight below
4100 if (std::abs(update_denominator) == 0.0)
4101 break;
4102 gradient_weight =
4103 gradient_weight - update_numerator / update_denominator;
4104
4105 // Put a fairly lenient bound on the largest possible
4106 // gradient (things tend to be locally flat, so the gradient
4107 // itself is usually small)
4108 if (std::abs(gradient_weight) > 10)
4109 {
4110 gradient_weight = -10.0;
4111 break;
4112 }
4113 }
4114
4115 // It only makes sense to take convex combinations with weights
4116 // between zero and one. If the update takes us outside of this
4117 // region then rescale the update to stay within the region and
4118 // try again
4119 Tensor<1, n_vertices_per_cell> tentative_weights =
4120 guess_weights + gradient_weight * current_gradient;
4121
4122 double new_gradient_weight = gradient_weight;
4123 for (unsigned int iteration_count = 0; iteration_count < 40;
4124 ++iteration_count)
4125 {
4126 if (internal::weights_are_ok<structdim>(tentative_weights))
4127 break;
4128
4129 for (unsigned int i = 0; i < n_vertices_per_cell; ++i)
4130 {
4131 if (tentative_weights[i] < 0.0)
4132 {
4133 tentative_weights -=
4134 (tentative_weights[i] / current_gradient[i]) *
4135 current_gradient;
4136 }
4137 if (tentative_weights[i] < 0.0 ||
4138 1.0 < tentative_weights[i])
4139 {
4140 new_gradient_weight /= 2.0;
4141 tentative_weights =
4142 guess_weights +
4143 new_gradient_weight * current_gradient;
4144 }
4145 }
4146 }
4147
4148 // the update might still send us outside the valid region, so
4149 // check again and quit if the update is still not valid
4150 if (!internal::weights_are_ok<structdim>(tentative_weights))
4151 break;
4152
4153 // if we cannot get closer by traveling in the gradient
4154 // direction then quit
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;
4159 else
4160 break;
4161 Assert(internal::weights_are_ok<structdim>(guess_weights),
4163 }
4164 Assert(internal::weights_are_ok<structdim>(guess_weights),
4166 projected_point = get_point_from_weights(guess_weights);
4167 }
4168
4169 // if structdim == 2 and the optimal point is not on the interior then
4170 // we may be able to get a more accurate result by projecting onto the
4171 // lines.
4172 if (structdim == 2)
4173 {
4174 std::array<Point<spacedim>, GeometryInfo<structdim>::lines_per_cell>
4175 line_projections;
4176 for (unsigned int line_n = 0;
4177 line_n < GeometryInfo<structdim>::lines_per_cell;
4178 ++line_n)
4179 {
4180 line_projections[line_n] =
4181 project_to_object(object->line(line_n), trial_point);
4182 }
4183 std::sort(line_projections.begin(),
4184 line_projections.end(),
4185 [&](const Point<spacedim> &a, const Point<spacedim> &b) {
4186 return a.distance(trial_point) <
4187 b.distance(trial_point);
4188 });
4189 if (line_projections[0].distance(trial_point) <
4190 projected_point.distance(trial_point))
4191 projected_point = line_projections[0];
4192 }
4193 }
4194 else
4195 {
4196 Assert(false, ExcNotImplemented());
4197 return projected_point;
4198 }
4199
4200 return projected_point;
4201 }
4202
4203
4204
4205 template <int dim, typename T>
4206 template <class Archive>
4207 void
4208 CellDataTransferBuffer<dim, T>::save(Archive &ar,
4209 const unsigned int /*version*/) const
4210 {
4211 Assert(cell_ids.size() == data.size(),
4212 ExcDimensionMismatch(cell_ids.size(), data.size()));
4213 // archive the cellids in an efficient binary format
4214 const std::size_t n_cells = cell_ids.size();
4215 ar & n_cells;
4216 for (const auto &id : cell_ids)
4217 {
4218 CellId::binary_type binary_cell_id = id.template to_binary<dim>();
4219 ar & binary_cell_id;
4220 }
4221
4222 ar &data;
4223 }
4224
4225
4226
4227 template <int dim, typename T>
4228 template <class Archive>
4229 void
4230 CellDataTransferBuffer<dim, T>::load(Archive &ar,
4231 const unsigned int /*version*/)
4232 {
4233 std::size_t n_cells;
4234 ar & n_cells;
4235 cell_ids.clear();
4236 cell_ids.reserve(n_cells);
4237 for (unsigned int c = 0; c < n_cells; ++c)
4238 {
4240 ar & value;
4241 cell_ids.emplace_back(value);
4242 }
4243 ar &data;
4244 }
4245
4246
4247 namespace internal
4248 {
4249 template <typename DataType,
4250 typename MeshType,
4251 typename MeshCellIteratorType>
4252 inline void
4253 exchange_cell_data(
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 &)>
4258 & unpack,
4259 const std::function<bool(const MeshCellIteratorType &)> & cell_filter,
4260 const std::function<void(
4261 const std::function<void(const MeshCellIteratorType &,
4262 const types::subdomain_id)> &)> &process_cells,
4263 const std::function<std::set<types::subdomain_id>(
4264 const parallel::TriangulationBase<MeshType::dimension,
4265 MeshType::space_dimension> &)>
4266 &compute_ghost_owners)
4267 {
4268# ifndef DEAL_II_WITH_MPI
4269 (void)mesh;
4270 (void)pack;
4271 (void)unpack;
4272 (void)cell_filter;
4273 (void)process_cells;
4274 (void)compute_ghost_owners;
4275 Assert(false, ExcNeedsMPI());
4276# else
4277 constexpr int dim = MeshType::dimension;
4278 constexpr int spacedim = MeshType::space_dimension;
4279 auto tria =
4280 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
4281 &mesh.get_triangulation());
4282 Assert(
4283 tria != nullptr,
4284 ExcMessage(
4285 "The function exchange_cell_data_to_ghosts() only works with parallel triangulations."));
4286
4287 if (const auto tria = dynamic_cast<
4289 &mesh.get_triangulation()))
4290 {
4291 Assert(
4292 tria->with_artificial_cells(),
4293 ExcMessage(
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."));
4300 }
4301
4302 // build list of cells to request for each neighbor
4303 std::set<::types::subdomain_id> ghost_owners =
4304 compute_ghost_owners(*tria);
4305 std::map<::types::subdomain_id,
4306 std::vector<typename CellId::binary_type>>
4307 neighbor_cell_list;
4308
4309 for (const auto ghost_owner : ghost_owners)
4310 neighbor_cell_list[ghost_owner] = {};
4311
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>());
4316 });
4317
4318 Assert(ghost_owners.size() == neighbor_cell_list.size(),
4320
4321
4322 // Before sending & receiving, make sure we protect this section with
4323 // a mutex:
4326 mutex, tria->get_communicator());
4327
4328 const int mpi_tag =
4330 const int mpi_tag_reply =
4332
4333 // send our requests:
4334 std::vector<MPI_Request> requests(ghost_owners.size());
4335 {
4336 unsigned int idx = 0;
4337 for (const auto &it : neighbor_cell_list)
4338 {
4339 // send the data about the relevant cells
4340 const int ierr = MPI_Isend(it.second.data(),
4341 it.second.size() * sizeof(it.second[0]),
4342 MPI_BYTE,
4343 it.first,
4344 mpi_tag,
4345 tria->get_communicator(),
4346 &requests[idx]);
4347 AssertThrowMPI(ierr);
4348 ++idx;
4349 }
4350 }
4351
4352 using DestinationToBufferMap =
4353 std::map<::types::subdomain_id,
4355 DestinationToBufferMap destination_to_data_buffer_map;
4356
4357 // receive requests and reply with the ghost indices
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());
4364
4365 for (unsigned int idx = 0; idx < ghost_owners.size(); ++idx)
4366 {
4367 MPI_Status status;
4368 int ierr = MPI_Probe(MPI_ANY_SOURCE,
4369 mpi_tag,
4370 tria->get_communicator(),
4371 &status);
4372 AssertThrowMPI(ierr);
4373
4374 int len;
4375 ierr = MPI_Get_count(&status, MPI_BYTE, &len);
4376 AssertThrowMPI(ierr);
4377 Assert(len % sizeof(cell_data_to_send[idx][0]) == 0,
4379
4380 const unsigned int n_cells =
4381 len / sizeof(typename CellId::binary_type);
4382 cell_data_to_send[idx].resize(n_cells);
4383
4384 ierr = MPI_Recv(cell_data_to_send[idx].data(),
4385 len,
4386 MPI_BYTE,
4387 status.MPI_SOURCE,
4388 status.MPI_TAG,
4389 tria->get_communicator(),
4390 &status);
4391 AssertThrowMPI(ierr);
4392
4393 // store data for each cell
4394 for (unsigned int c = 0; c < static_cast<unsigned int>(n_cells); ++c)
4395 {
4396 const auto cell =
4397 tria->create_cell_iterator(CellId(cell_data_to_send[idx][c]));
4398
4399 MeshCellIteratorType mesh_it(tria,
4400 cell->level(),
4401 cell->index(),
4402 &mesh);
4403 const std_cxx17::optional<DataType> data = pack(mesh_it);
4404
4405 if (data)
4406 {
4407 typename DestinationToBufferMap::iterator p =
4408 destination_to_data_buffer_map
4409 .insert(std::make_pair(
4410 idx,
4412 .first;
4413
4414 p->second.cell_ids.emplace_back(cell->id());
4415 p->second.data.emplace_back(*data);
4416 }
4417 }
4418
4419 // send reply
4421 destination_to_data_buffer_map[idx];
4422
4423 sendbuffers[idx] =
4424 Utilities::pack(data, /*enable_compression*/ false);
4425 ierr = MPI_Isend(sendbuffers[idx].data(),
4426 sendbuffers[idx].size(),
4427 MPI_BYTE,
4428 status.MPI_SOURCE,
4429 mpi_tag_reply,
4430 tria->get_communicator(),
4431 &reply_requests[idx]);
4432 AssertThrowMPI(ierr);
4433 }
4434
4435 // finally receive the replies
4436 std::vector<char> receive;
4437 for (unsigned int idx = 0; idx < ghost_owners.size(); ++idx)
4438 {
4439 MPI_Status status;
4440 int ierr = MPI_Probe(MPI_ANY_SOURCE,
4441 mpi_tag_reply,
4442 tria->get_communicator(),
4443 &status);
4444 AssertThrowMPI(ierr);
4445
4446 int len;
4447 ierr = MPI_Get_count(&status, MPI_BYTE, &len);
4448 AssertThrowMPI(ierr);
4449
4450 receive.resize(len);
4451
4452 char *ptr = receive.data();
4453 ierr = MPI_Recv(ptr,
4454 len,
4455 MPI_BYTE,
4456 status.MPI_SOURCE,
4457 status.MPI_TAG,
4458 tria->get_communicator(),
4459 &status);
4460 AssertThrowMPI(ierr);
4461
4462 auto cellinfo =
4463 Utilities::unpack<CellDataTransferBuffer<dim, DataType>>(
4464 receive, /*enable_compression*/ false);
4465
4466 DataType *data = cellinfo.data.data();
4467 for (unsigned int c = 0; c < cellinfo.cell_ids.size(); ++c, ++data)
4468 {
4470 tria_cell = tria->create_cell_iterator(cellinfo.cell_ids[c]);
4471
4472 MeshCellIteratorType cell(tria,
4473 tria_cell->level(),
4474 tria_cell->index(),
4475 &mesh);
4476
4477 unpack(cell, *data);
4478 }
4479 }
4480
4481 // make sure that all communication is finished
4482 // when we leave this function.
4483 if (requests.size() > 0)
4484 {
4485 const int ierr =
4486 MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
4487 AssertThrowMPI(ierr);
4488 }
4489 if (reply_requests.size() > 0)
4490 {
4491 const int ierr = MPI_Waitall(reply_requests.size(),
4492 reply_requests.data(),
4493 MPI_STATUSES_IGNORE);
4494 AssertThrowMPI(ierr);
4495 }
4496
4497
4498# endif // DEAL_II_WITH_MPI
4499 }
4500
4501 } // namespace internal
4502
4503 template <typename DataType, typename MeshType>
4504 inline void
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 &)>
4512 &cell_filter)
4513 {
4514# ifndef DEAL_II_WITH_MPI
4515 (void)mesh;
4516 (void)pack;
4517 (void)unpack;
4518 (void)cell_filter;
4519 Assert(false, ExcNeedsMPI());
4520# else
4521 internal::exchange_cell_data<DataType,
4522 MeshType,
4523 typename MeshType::active_cell_iterator>(
4524 mesh,
4525 pack,
4526 unpack,
4527 cell_filter,
4528 [&](const auto &process) {
4529 for (const auto &cell : mesh.active_cell_iterators())
4530 if (cell->is_ghost())
4531 process(cell, cell->subdomain_id());
4532 },
4533 [](const auto &tria) { return tria.ghost_owners(); });
4534# endif
4535 }
4536
4537
4538
4539 template <typename DataType, typename MeshType>
4540 inline void
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 &)>
4548 &cell_filter)
4549 {
4550# ifndef DEAL_II_WITH_MPI
4551 (void)mesh;
4552 (void)pack;
4553 (void)unpack;
4554 (void)cell_filter;
4555 Assert(false, ExcNeedsMPI());
4556# else
4557 internal::exchange_cell_data<DataType,
4558 MeshType,
4559 typename MeshType::level_cell_iterator>(
4560 mesh,
4561 pack,
4562 unpack,
4563 cell_filter,
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());
4570 },
4571 [](const auto &tria) { return tria.level_ghost_owners(); });
4572# endif
4573 }
4574} // namespace GridTools
4575
4576# endif
4577
4579
4580/*---------------------------- grid_tools.h ---------------------------*/
4581/* end of #ifndef dealii_grid_tools_h */
4582#endif
4583/*---------------------------- grid_tools.h ---------------------------*/
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition: array_view.h:697
Definition: cell_id.h:71
std::array< unsigned int, 4 > binary_type
Definition: cell_id.h:80
typename VectorType::value_type value_type
Definition: grid_tools.h:3304
void process_sub_cell(const std::vector< value_type > &ls_values, const std::vector< Point< 2 > > &points, const std::vector< unsigned int > mask, const double iso_level, std::vector< Point< 2 > > &vertices, std::vector< CellData< 1 > > &cells) const
Definition: grid_tools.cc:6797
void process_cell(const typename DoFHandler< dim >::active_cell_iterator &cell, const VectorType &ls_vector, const double iso_level, std::vector< Point< dim > > &vertices, std::vector< CellData< dim - 1 > > &cells) const
Definition: grid_tools.cc:6632
void process(const DoFHandler< dim > &background_dof_handler, const VectorType &ls_vector, const double iso_level, std::vector< Point< dim > > &vertices, std::vector< CellData< dim - 1 > > &cells) const
Definition: grid_tools.cc:6616
const unsigned int n_subdivisions
Definition: grid_tools.h:3385
MarchingCubeAlgorithm(const Mapping< dim, dim > &mapping, const FiniteElement< dim, dim > &fe, const unsigned int n_subdivisions=1, const double tolerance=1e-10)
Definition: grid_tools.cc:6569
static Quadrature< dim > create_quadrature_rule(const unsigned int n_subdivisions)
Definition: grid_tools.cc:6586
Definition: point.h:111
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
Definition: vector.h:110
typename MeshType::active_cell_iterator type
Definition: grid_tools.h:99
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition: config.h:416
#define DEAL_II_DEPRECATED_EARLY
Definition: config.h:170
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition: config.h:454
Point< 3 > center
Point< 3 > vertices[4]
unsigned int level
Definition: grid_out.cc:4590
unsigned int vertex_indices[2]
Definition: grid_tools.cc:1256
const double angle
__global__ void set(Number *val, const Number s, const size_type N)
#define DeclException0(Exception0)
Definition: exceptions.h:470
static ::ExceptionBase & ExcMeshNotOrientable()
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcNeedsMPI()
static ::ExceptionBase & ExcNonExistentSubdomain(int arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcVertexNotUsed(unsigned int arg1)
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1746
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcTriangulationHasBeenRefined()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:515
static ::ExceptionBase & ExcInvalidNumberOfPartitions(int arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcScalingFactorNotPositive(double arg1)
typename ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:438
RTree< std::pair< BoundingBox< spacedim >, unsigned int > > build_global_description_tree(const std::vector< BoundingBox< spacedim > > &local_description, const MPI_Comm &mpi_communicator)
Definition: grid_tools.cc:6325
void copy_boundary_to_manifold_id(Triangulation< dim, spacedim > &tria, const bool reset_boundary_ids=false)
Definition: grid_tools.cc:4808
void copy_material_to_manifold_id(Triangulation< dim, spacedim > &tria, const bool compute_face_ids=false)
Definition: grid_tools.cc:4900
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={})
Definition: grid_tools.cc:4833
std::vector< std::vector< BoundingBox< spacedim > > > exchange_local_bounding_boxes(const std::vector< BoundingBox< spacedim > > &local_bboxes, const MPI_Comm &mpi_communicator)
Definition: grid_tools.cc:6229
std::map< unsigned int, std::set<::types::subdomain_id > > compute_vertices_with_ghost_neighbors(const Triangulation< dim, spacedim > &tria)
Definition: grid_tools.cc:6511
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)
Definition: grid_tools.cc:6388
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)
Definition: grid_tools.cc:4928
DistributedComputePointLocationsInternal< dim, spacedim > distributed_compute_point_locations(const GridTools::Cache< dim, spacedim > &cache, const std::vector< Point< spacedim > > &points, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bboxes, const double tolerance, const bool perform_handshake, const bool enforce_unique_mapping=false)
Definition: grid_tools.cc:5779
void get_face_connectivity_of_cells(const Triangulation< dim, spacedim > &triangulation, DynamicSparsityPattern &connectivity)
Definition: grid_tools.cc:3663
return_type distributed_compute_point_locations(const GridTools::Cache< dim, spacedim > &cache, const std::vector< Point< spacedim > > &local_points, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bboxes, const double tolerance=1e-10)
Definition: grid_tools.cc:5622
void delete_unused_vertices(std::vector< Point< spacedim > > &vertices, std::vector< CellData< dim > > &cells, SubCellData &subcelldata)
Definition: grid_tools.cc:626
std::vector< BoundingBox< MeshType::space_dimension > > compute_mesh_predicate_bounding_box(const MeshType &mesh, const std::function< bool(const typename MeshType::active_cell_iterator &)> &predicate, const unsigned int refinement_level=0, const bool allow_merge=false, const unsigned int max_boxes=numbers::invalid_unsigned_int)
Definition: grid_tools.cc:3016
void scale(const double scaling_factor, Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:2042
std::vector< typename MeshType::active_cell_iterator > get_active_child_cells(const typename MeshType::cell_iterator &cell)
void partition_triangulation_zorder(const unsigned int n_partitions, Triangulation< dim, spacedim > &triangulation, const bool group_siblings=true)
Definition: grid_tools.cc:3989
std::map< unsigned int, Point< spacedim > > get_all_vertices_at_boundary(const Triangulation< dim, spacedim > &tria)
Definition: grid_tools.cc:2185
Triangulation< dim, spacedim >::DistortedCellList fix_up_distorted_child_cells(const typename Triangulation< dim, spacedim >::DistortedCellList &distorted_cells, Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:4766
std::vector< typename Container::cell_iterator > get_cells_at_coarsest_common_level(const std::vector< typename Container::active_cell_iterator > &patch_cells)
return_type compute_point_locations(const Cache< dim, spacedim > &cache, const std::vector< Point< spacedim > > &points, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell_hint=typename Triangulation< dim, spacedim >::active_cell_iterator())
Definition: grid_tools.cc:5425
unsigned int find_closest_vertex(const std::map< unsigned int, Point< spacedim > > &vertices, const Point< spacedim > &p)
Definition: grid_tools.cc:6184
void transform(const Transformation &transformation, Triangulation< dim, spacedim > &triangulation)
unsigned int find_closest_vertex_of_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell, const Point< spacedim > &position, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
Definition: grid_tools.cc:2934
std::pair< typename MeshType< dim, spacedim >::active_cell_iterator, Point< dim > > find_active_cell_around_point(const Mapping< dim, spacedim > &mapping, const MeshType< dim, spacedim > &mesh, const Point< spacedim > &p, const std::vector< bool > &marked_vertices={}, const double tolerance=1.e-10)
std::map< unsigned int, Point< spacedim > > extract_used_vertices(const Triangulation< dim, spacedim > &container, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
Definition: grid_tools.cc:6165
void get_active_neighbors(const typename MeshType::active_cell_iterator &cell, std::vector< typename MeshType::active_cell_iterator > &active_neighbors)
std::list< std::pair< typename MeshType::cell_iterator, typename MeshType::cell_iterator > > get_finest_common_cells(const MeshType &mesh_1, const MeshType &mesh_2)
std::vector< bool > get_locally_owned_vertices(const Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:4246
std::pair< DerivativeForm< 1, dim, spacedim >, Tensor< 1, spacedim > > affine_cell_approximation(const ArrayView< const Point< spacedim > > &vertices)
Definition: grid_tools.cc:288
void regularize_corner_cells(Triangulation< dim, spacedim > &tria, const double limit_angle_fraction=.75)
Definition: grid_tools.cc:5110
void remove_anisotropy(Triangulation< dim, spacedim > &tria, const double max_ratio=1.6180339887, const unsigned int max_iterations=5)
Definition: grid_tools.cc:5081
void consistently_order_cells(std::vector< CellData< dim > > &cells)
Definition: grid_tools.cc:1921
double minimal_cell_diameter(const Triangulation< dim, spacedim > &triangulation, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
Definition: grid_tools.cc:4273
void shift(const Tensor< 1, spacedim > &shift_vector, Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:2022
std::vector< typename MeshType::active_cell_iterator > compute_active_cell_halo_layer(const MeshType &mesh, const std::function< bool(const typename MeshType::active_cell_iterator &)> &predicate)
double cell_measure(const std::vector< Point< dim > > &all_vertices, const unsigned int(&vertex_indices)[GeometryInfo< dim >::vertices_per_cell])
std::vector< typename MeshType::active_cell_iterator > compute_active_cell_layer_within_distance(const MeshType &mesh, const std::function< bool(const typename MeshType::active_cell_iterator &)> &predicate, const double layer_thickness)
void remove_hanging_nodes(Triangulation< dim, spacedim > &tria, const bool isotropic=false, const unsigned int max_iterations=100)
Definition: grid_tools.cc:5048
std::vector< typename MeshType::active_cell_iterator > compute_ghost_cell_halo_layer(const MeshType &mesh)
Point< Iterator::AccessorType::space_dimension > project_to_object(const Iterator &object, const Point< Iterator::AccessorType::space_dimension > &trial_point)
void laplace_transform(const std::map< unsigned int, Point< dim > > &new_points, Triangulation< dim > &tria, const Function< dim, double > *coefficient=nullptr, const bool solve_for_absolute_positions=false)
void partition_multigrid_levels(Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:4094
void delete_duplicated_vertices(std::vector< Point< spacedim > > &all_vertices, std::vector< CellData< dim > > &cells, SubCellData &subcelldata, std::vector< unsigned int > &considered_vertices, const double tol=1e-12)
Definition: grid_tools.cc:731
unsigned int count_cells_with_subdomain_association(const Triangulation< dim, spacedim > &triangulation, const types::subdomain_id subdomain)
Definition: grid_tools.cc:4230
return_type guess_point_owner(const std::vector< std::vector< BoundingBox< spacedim > > > &global_bboxes, const std::vector< Point< spacedim > > &points)
Definition: grid_tools.cc:3163
std::map< types::global_dof_index, std::vector< typename DoFHandler< dim, spacedim >::active_cell_iterator > > get_dof_to_support_patch_map(DoFHandler< dim, spacedim > &dof_handler)
void rotate(const double angle, Triangulation< dim > &triangulation)
std::vector< typename MeshType::cell_iterator > compute_cell_halo_layer_on_level(const MeshType &mesh, const std::function< bool(const typename MeshType::cell_iterator &)> &predicate, const unsigned int level)
void partition_triangulation(const unsigned int n_partitions, Triangulation< dim, spacedim > &triangulation, const SparsityTools::Partitioner partitioner=SparsityTools::Partitioner::metis)
Definition: grid_tools.cc:3762
Vector< double > compute_aspect_ratio_of_cells(const Mapping< dim > &mapping, const Triangulation< dim > &triangulation, const Quadrature< dim > &quadrature)
Definition: grid_tools.cc:312
std::vector< std::set< typename Triangulation< dim, spacedim >::active_cell_iterator > > vertex_to_cell_map(const Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:3261
double volume(const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
Definition: grid_tools.cc:137
std::vector< typename MeshType::active_cell_iterator > compute_ghost_cell_layer_within_distance(const MeshType &mesh, const double layer_thickness)
double compute_maximum_aspect_ratio(const Mapping< dim > &mapping, const Triangulation< dim > &triangulation, const Quadrature< dim > &quadrature)
Definition: grid_tools.cc:378
return_type compute_point_locations_try_all(const Cache< dim, spacedim > &cache, const std::vector< Point< spacedim > > &points, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell_hint=typename Triangulation< dim, spacedim >::active_cell_iterator())
Definition: grid_tools.cc:5454
bool have_same_coarse_mesh(const Triangulation< dim, spacedim > &mesh_1, const Triangulation< dim, spacedim > &mesh_2)
void collect_periodic_faces(const MeshType &mesh, const types::boundary_id b_id1, const types::boundary_id b_id2, const int direction, std::vector< PeriodicFacePair< typename MeshType::cell_iterator > > &matched_pairs, const Tensor< 1, MeshType::space_dimension > &offset=::Tensor< 1, MeshType::space_dimension >(), const FullMatrix< double > &matrix=FullMatrix< double >())
std::vector< typename MeshType::active_cell_iterator > get_patch_around_cell(const typename MeshType::active_cell_iterator &cell)
void get_vertex_connectivity_of_cells(const Triangulation< dim, spacedim > &triangulation, DynamicSparsityPattern &connectivity)
Definition: grid_tools.cc:3697
bool orthogonal_equality(std::bitset< 3 > &orientation, const FaceIterator &face1, const FaceIterator &face2, const int direction, const Tensor< 1, FaceIterator::AccessorType::space_dimension > &offset=Tensor< 1, FaceIterator::AccessorType::space_dimension >(), const FullMatrix< double > &matrix=FullMatrix< double >())
std::vector< types::subdomain_id > get_subdomain_association(const Triangulation< dim, spacedim > &triangulation, const std::vector< CellId > &cell_ids)
Definition: grid_tools.cc:4122
void distort_random(const double factor, Triangulation< dim, spacedim > &triangulation, const bool keep_boundary=true, const unsigned int seed=boost::random::mt19937::default_seed)
Definition: grid_tools.cc:2217
std::vector< std::vector< Tensor< 1, spacedim > > > vertex_to_cell_centers_directions(const Triangulation< dim, spacedim > &mesh, const std::vector< std::set< typename Triangulation< dim, spacedim >::active_cell_iterator > > &vertex_to_cells)
Definition: grid_tools.cc:2696
void build_triangulation_from_patch(const std::vector< typename Container::active_cell_iterator > &patch, Triangulation< Container::dimension, Container::space_dimension > &local_triangulation, std::map< typename Triangulation< Container::dimension, Container::space_dimension >::active_cell_iterator, typename Container::active_cell_iterator > &patch_to_global_tria_map)
double diameter(const Triangulation< dim, spacedim > &tria)
Definition: grid_tools.cc:81
void get_vertex_connectivity_of_cells_on_level(const Triangulation< dim, spacedim > &triangulation, const unsigned int level, DynamicSparsityPattern &connectivity)
Definition: grid_tools.cc:3726
void invert_all_negative_measure_cells(const std::vector< Point< spacedim > > &all_vertices, std::vector< CellData< dim > > &cells)
Definition: grid_tools.cc:863
std::vector< typename MeshType< dim, spacedim >::active_cell_iterator > find_cells_adjacent_to_vertex(const MeshType< dim, spacedim > &container, const unsigned int vertex_index)
Definition: grid_tools.cc:2568
std::map< unsigned int, types::global_vertex_index > compute_local_to_global_vertex_index_map(const parallel::distributed::Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:3310
std::vector< std::pair< typename MeshType< dim, spacedim >::active_cell_iterator, Point< dim > > > find_all_active_cells_around_point(const Mapping< dim, spacedim > &mapping, const MeshType< dim, spacedim > &mesh, const Point< spacedim > &p, const double tolerance, const std::pair< typename MeshType< dim, spacedim >::active_cell_iterator, Point< dim > > &first_cell)
BoundingBox< spacedim > compute_bounding_box(const Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:394
double maximal_cell_diameter(const Triangulation< dim, spacedim > &triangulation, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
Definition: grid_tools.cc:4300
std::pair< unsigned int, double > get_longest_direction(typename Triangulation< dim, spacedim >::active_cell_iterator cell)
Definition: grid_tools.cc:5016
std::tuple< std::vector< Point< spacedim > >, std::vector< CellData< dim > >, SubCellData > get_coarse_mesh_description(const Triangulation< dim, spacedim > &tria)
Definition: grid_tools.cc:534
@ 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)
@ exchange_cell_data_request
grid_tools.h: exchange_cell_ghosts()
Definition: mpi_tags.h:66
@ exchange_cell_data_reply
grid_tools.h: exchange_cell_ghosts()
Definition: mpi_tags.h:69
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)
Definition: utilities.h:1218
T unpack(const std::vector< char > &buffer, const bool allow_compression=true)
Definition: utilities.h:1321
Definition: hp.h:118
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
Definition: tria.cc:12587
const types::subdomain_id artificial_subdomain_id
Definition: types.h:293
static const unsigned int invalid_unsigned_int
Definition: types.h:196
const types::manifold_id flat_manifold_id
Definition: types.h:264
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int global_dof_index
Definition: types.h:76
unsigned int subdomain_id
Definition: types.h:43
unsigned int boundary_id
Definition: types.h:129
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
boost::geometry::index::rtree< LeafType, IndexType, IndexableGetter > RTree
Definition: rtree.h:153
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)
void save(Archive &ar, const unsigned int version) const
void load(Archive &ar, const unsigned int version)
void serialize(Archive &archive, const unsigned int version)
std::vector< CellId > cell_ids
Definition: grid_tools.h:3233
std::bitset< 3 > orientation
Definition: grid_tools.h:2664
FullMatrix< double > matrix
Definition: grid_tools.h:2678
unsigned int face_idx[2]
Definition: grid_tools.h:2657
std::vector< std::tuple< unsigned int, unsigned int, unsigned int > > recv_components
Definition: grid_tools.h:1124
std::vector< std::tuple< std::pair< int, int >, unsigned int, unsigned int, Point< dim >, Point< spacedim >, unsigned int > > send_components
Definition: grid_tools.h:1099