Reference documentation for deal.II version 9.1.0-pre
grid_tools.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2018 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 
22 # include <deal.II/base/bounding_box.h>
23 # include <deal.II/base/geometry_info.h>
24 
25 # include <deal.II/dofs/dof_handler.h>
26 
27 # include <deal.II/fe/mapping.h>
28 # include <deal.II/fe/mapping_q1.h>
29 
30 # include <deal.II/grid/manifold.h>
31 # include <deal.II/grid/tria.h>
32 # include <deal.II/grid/tria_accessor.h>
33 # include <deal.II/grid/tria_iterator.h>
34 
35 # include <deal.II/hp/dof_handler.h>
36 
37 # include <deal.II/lac/sparsity_tools.h>
38 
39 # include <boost/archive/binary_iarchive.hpp>
40 # include <boost/archive/binary_oarchive.hpp>
41 # include <boost/optional.hpp>
42 # include <boost/serialization/array.hpp>
43 # include <boost/serialization/vector.hpp>
44 
45 # ifdef DEAL_II_WITH_ZLIB
46 # include <boost/iostreams/device/back_inserter.hpp>
47 # include <boost/iostreams/filter/gzip.hpp>
48 # include <boost/iostreams/filtering_stream.hpp>
49 # include <boost/iostreams/stream.hpp>
50 # endif
51 
52 
53 # include <bitset>
54 # include <list>
55 # include <set>
56 
57 DEAL_II_NAMESPACE_OPEN
58 
59 namespace parallel
60 {
61  namespace distributed
62  {
63  template <int, int>
64  class Triangulation;
65  }
66 } // namespace parallel
67 
68 namespace hp
69 {
70  template <int, int>
71  class MappingCollection;
72 }
73 
74 class SparsityPattern;
75 
76 namespace internal
77 {
78  template <int dim, int spacedim, class MeshType>
79  class ActiveCellIterator
80  {
81  public:
82 # ifndef _MSC_VER
83  using type = typename MeshType::active_cell_iterator;
84 # else
86 # endif
87  };
88 
89 # ifdef _MSC_VER
90  template <int dim, int spacedim>
91  class ActiveCellIterator<dim, spacedim, ::DoFHandler<dim, spacedim>>
92  {
93  public:
94  using type = TriaActiveIterator<
96  };
97 
98  template <int dim, int spacedim>
99  class ActiveCellIterator<dim, spacedim, ::hp::DoFHandler<dim, spacedim>>
100  {
101  public:
102  using type = TriaActiveIterator<
104  };
105 # endif
106 } // namespace internal
107 
116 namespace GridTools
117 {
118  template <int dim, int spacedim>
119  class Cache;
120 
125 
132  template <int dim, int spacedim>
133  double
135 
162  template <int dim, int spacedim>
163  double
165  const Mapping<dim, spacedim> & mapping =
167 
178  template <int dim, int spacedim>
179  double
181  const Mapping<dim, spacedim> & mapping =
183 
194  template <int dim, int spacedim>
195  double
197  const Mapping<dim, spacedim> & mapping =
199 
209  template <int dim>
210  double
211  cell_measure(
212  const std::vector<Point<dim>> &all_vertices,
213  const unsigned int (&vertex_indices)[GeometryInfo<dim>::vertices_per_cell]);
214 
220  template <int dim, typename T>
221  double
222  cell_measure(const T &, ...);
223 
237  template <int dim, int spacedim>
240 
260  template <typename Iterator>
263  const Iterator & object,
265 
271 
288  template <int dim, int spacedim>
289  void
290  delete_unused_vertices(std::vector<Point<spacedim>> &vertices,
291  std::vector<CellData<dim>> & cells,
292  SubCellData & subcelldata);
293 
310  template <int dim, int spacedim>
311  void
312  delete_duplicated_vertices(std::vector<Point<spacedim>> &all_vertices,
313  std::vector<CellData<dim>> & cells,
314  SubCellData & subcelldata,
315  std::vector<unsigned int> & considered_vertices,
316  const double tol = 1e-12);
317 
323 
347  template <int dim, typename Transformation, int spacedim>
348  void
349  transform(const Transformation & transformation,
350  Triangulation<dim, spacedim> &triangulation);
351 
357  template <int dim, int spacedim>
358  void
359  shift(const Tensor<1, spacedim> & shift_vector,
360  Triangulation<dim, spacedim> &triangulation);
361 
362 
370  void
371  rotate(const double angle, Triangulation<2> &triangulation);
372 
385  template <int dim>
386  void
387  rotate(const double angle,
388  const unsigned int axis,
389  Triangulation<dim, 3> &triangulation);
390 
448  template <int dim>
449  void
450  laplace_transform(const std::map<unsigned int, Point<dim>> &new_points,
451  Triangulation<dim> & tria,
452  const Function<dim, double> *coefficient = nullptr,
453  const bool solve_for_absolute_positions = false);
454 
460  template <int dim, int spacedim>
461  std::map<unsigned int, Point<spacedim>>
463 
471  template <int dim, int spacedim>
472  void
473  scale(const double scaling_factor,
474  Triangulation<dim, spacedim> &triangulation);
475 
486  template <int dim, int spacedim>
487  void
488  distort_random(const double factor,
489  Triangulation<dim, spacedim> &triangulation,
490  const bool keep_boundary = true);
491 
527  template <int dim, int spacedim>
528  void
530  const bool isotropic = false,
531  const unsigned int max_iterations = 100);
532 
559  template <int dim, int spacedim>
560  void
562  const double max_ratio = 1.6180339887,
563  const unsigned int max_iterations = 5);
564 
656  template <int dim, int spacedim>
657  void
659  const double limit_angle_fraction = .75);
660 
666 
708  template <int dim, int spacedim>
709 # ifndef DOXYGEN
710  std::tuple<
711  std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
712  std::vector<std::vector<Point<dim>>>,
713  std::vector<std::vector<unsigned int>>>
714 # else
715  return_type
716 # endif
718  const Cache<dim, spacedim> & cache,
719  const std::vector<Point<spacedim>> &points,
721  &cell_hint =
723 
788  template <int dim, int spacedim>
789 # ifndef DOXYGEN
790  std::tuple<
791  std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>,
792  std::vector<std::vector<Point<dim>>>,
793  std::vector<std::vector<unsigned int>>,
794  std::vector<std::vector<Point<spacedim>>>,
795  std::vector<std::vector<unsigned int>>>
796 # else
797  return_type
798 # endif
800  const GridTools::Cache<dim, spacedim> & cache,
801  const std::vector<Point<spacedim>> & local_points,
802  const std::vector<std::vector<BoundingBox<spacedim>>> &global_bboxes);
803 
842  template <int dim, int spacedim>
843  std::map<unsigned int, Point<spacedim>>
845  const Mapping<dim, spacedim> & mapping =
847 
859  template <int spacedim>
860  unsigned int
861  find_closest_vertex(const std::map<unsigned int, Point<spacedim>> &vertices,
862  const Point<spacedim> & p);
863 
890  template <int dim, template <int, int> class MeshType, int spacedim>
891  unsigned int
893  const MeshType<dim, spacedim> &mesh,
894  const Point<spacedim> & p,
895  const std::vector<bool> & marked_vertices = std::vector<bool>());
896 
922  template <int dim, template <int, int> class MeshType, int spacedim>
923  unsigned int
925  const Mapping<dim, spacedim> & mapping,
926  const MeshType<dim, spacedim> &mesh,
927  const Point<spacedim> & p,
928  const std::vector<bool> & marked_vertices = std::vector<bool>());
929 
930 
955  template <int dim, template <int, int> class MeshType, int spacedim>
956 # ifndef _MSC_VER
957  std::vector<typename MeshType<dim, spacedim>::active_cell_iterator>
958 # else
959  std::vector<
960  typename ::internal::
961  ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type>
962 # endif
963  find_cells_adjacent_to_vertex(const MeshType<dim, spacedim> &container,
964  const unsigned int vertex_index);
965 
966 
991  template <int dim, template <int, int> class MeshType, int spacedim>
992 # ifndef _MSC_VER
993  typename MeshType<dim, spacedim>::active_cell_iterator
994 # else
995  typename ::internal::
996  ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type
997 # endif
999  const MeshType<dim, spacedim> &mesh,
1000  const Point<spacedim> & p,
1001  const std::vector<bool> & marked_vertices = std::vector<bool>());
1002 
1089  template <int dim, template <int, int> class MeshType, int spacedim>
1090 # ifndef _MSC_VER
1091  std::pair<typename MeshType<dim, spacedim>::active_cell_iterator, Point<dim>>
1092 # else
1093  std::pair<typename ::internal::
1094  ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1095  Point<dim>>
1096 # endif
1098  const Mapping<dim, spacedim> & mapping,
1099  const MeshType<dim, spacedim> &mesh,
1100  const Point<spacedim> & p,
1101  const std::vector<bool> & marked_vertices = std::vector<bool>());
1102 
1112  template <int dim, template <int, int> class MeshType, int spacedim>
1113 # ifndef _MSC_VER
1114  std::pair<typename MeshType<dim, spacedim>::active_cell_iterator, Point<dim>>
1115 # else
1116  std::pair<typename ::internal::
1117  ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1118  Point<dim>>
1119 # endif
1121  const Mapping<dim, spacedim> & mapping,
1122  const MeshType<dim, spacedim> &mesh,
1123  const Point<spacedim> & p,
1124  const std::vector<
1125  std::set<typename MeshType<dim, spacedim>::active_cell_iterator>>
1126  & vertex_to_cell_map,
1127  const std::vector<std::vector<Tensor<1, spacedim>>> &vertex_to_cell_centers,
1128  const typename MeshType<dim, spacedim>::active_cell_iterator &cell_hint =
1129  typename MeshType<dim, spacedim>::active_cell_iterator(),
1130  const std::vector<bool> &marked_vertices = std::vector<bool>());
1131 
1153  template <int dim, int spacedim>
1154  std::pair<typename hp::DoFHandler<dim, spacedim>::active_cell_iterator,
1155  Point<dim>>
1157  const hp::MappingCollection<dim, spacedim> &mapping,
1158  const hp::DoFHandler<dim, spacedim> & mesh,
1159  const Point<spacedim> & p);
1160 
1167  template <int dim, int spacedim>
1168  std::pair<typename Triangulation<dim, spacedim>::active_cell_iterator,
1169  Point<dim>>
1171  const Cache<dim, spacedim> &cache,
1172  const Point<spacedim> & p,
1175  const std::vector<bool> &marked_vertices = std::vector<bool>());
1176 
1194  template <int dim, template <int, int> class MeshType, int spacedim>
1195 # ifndef _MSC_VER
1196  std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
1197  Point<dim>>>
1198 # else
1199  std::vector<std::pair<
1200  typename ::internal::
1201  ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
1202  Point<dim>>>
1203 # endif
1205  const Mapping<dim, spacedim> & mapping,
1206  const MeshType<dim, spacedim> &mesh,
1207  const Point<spacedim> & p,
1208  const double tolerance = 1e-12,
1209  const std::vector<bool> & marked_vertices = std::vector<bool>());
1210 
1232  template <class MeshType>
1233  std::vector<typename MeshType::active_cell_iterator>
1234  get_active_child_cells(const typename MeshType::cell_iterator &cell);
1235 
1246  template <class MeshType>
1247  void
1249  const typename MeshType::active_cell_iterator & cell,
1250  std::vector<typename MeshType::active_cell_iterator> &active_neighbors);
1251 
1304  template <class MeshType>
1305  std::vector<typename MeshType::active_cell_iterator>
1307  const MeshType &mesh,
1308  const std::function<bool(const typename MeshType::active_cell_iterator &)>
1309  &predicate);
1310 
1311 
1319  template <class MeshType>
1320  std::vector<typename MeshType::cell_iterator>
1322  const MeshType &mesh,
1323  const std::function<bool(const typename MeshType::cell_iterator &)>
1324  & predicate,
1325  const unsigned int level);
1326 
1327 
1343  template <class MeshType>
1344  std::vector<typename MeshType::active_cell_iterator>
1345  compute_ghost_cell_halo_layer(const MeshType &mesh);
1346 
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 &)>
1403  & predicate,
1404  const double layer_thickness);
1405 
1431  template <class MeshType>
1432  std::vector<typename MeshType::active_cell_iterator>
1433  compute_ghost_cell_layer_within_distance(const MeshType &mesh,
1434  const double layer_thickness);
1435 
1451  template <class MeshType>
1452  std::pair<Point<MeshType::space_dimension>, Point<MeshType::space_dimension>>
1454  const MeshType &mesh,
1455  const std::function<bool(const typename MeshType::active_cell_iterator &)>
1456  &predicate);
1457 
1510  template <class MeshType>
1511  std::vector<BoundingBox<MeshType::space_dimension>>
1513  const MeshType &mesh,
1514  const std::function<bool(const typename MeshType::active_cell_iterator &)>
1515  & predicate,
1516  const unsigned int refinement_level = 0,
1517  const bool allow_merge = false,
1518  const unsigned int max_boxes = numbers::invalid_unsigned_int);
1519 
1549  template <int spacedim>
1550 # ifndef DOXYGEN
1551  std::tuple<std::vector<std::vector<unsigned int>>,
1552  std::map<unsigned int, unsigned int>,
1553  std::map<unsigned int, std::vector<unsigned int>>>
1554 # else
1555  return_type
1556 # endif
1558  const std::vector<std::vector<BoundingBox<spacedim>>> &global_bboxes,
1559  const std::vector<Point<spacedim>> & points);
1560 
1561 
1570  template <int dim, int spacedim>
1571  std::vector<
1572  std::set<typename Triangulation<dim, spacedim>::active_cell_iterator>>
1573  vertex_to_cell_map(const Triangulation<dim, spacedim> &triangulation);
1574 
1589  template <int dim, int spacedim>
1590  std::vector<std::vector<Tensor<1, spacedim>>>
1592  const Triangulation<dim, spacedim> &mesh,
1593  const std::vector<
1595  &vertex_to_cells);
1596 
1597 
1604  template <int dim, int spacedim>
1605  unsigned int
1608  const Point<spacedim> &position);
1609 
1621  template <int dim, int spacedim>
1622  std::map<unsigned int, types::global_vertex_index>
1625 
1639  template <int dim, int spacedim>
1640  std::pair<unsigned int, double>
1643 
1649 
1658  template <int dim, int spacedim>
1659  void
1661  const Triangulation<dim, spacedim> &triangulation,
1662  DynamicSparsityPattern & connectivity);
1663 
1672  template <int dim, int spacedim>
1673  void
1675  const Triangulation<dim, spacedim> &triangulation,
1676  DynamicSparsityPattern & connectivity);
1677 
1686  template <int dim, int spacedim>
1687  void
1689  const Triangulation<dim, spacedim> &triangulation,
1690  const unsigned int level,
1691  DynamicSparsityPattern & connectivity);
1692 
1713  template <int dim, int spacedim>
1714  void
1715  partition_triangulation(const unsigned int n_partitions,
1716  Triangulation<dim, spacedim> & triangulation,
1717  const SparsityTools::Partitioner partitioner =
1719 
1730  template <int dim, int spacedim>
1731  void
1732  partition_triangulation(const unsigned int n_partitions,
1733  const std::vector<unsigned int> &cell_weights,
1734  Triangulation<dim, spacedim> & triangulation,
1735  const SparsityTools::Partitioner partitioner =
1737 
1783  template <int dim, int spacedim>
1784  void
1785  partition_triangulation(const unsigned int n_partitions,
1786  const SparsityPattern & cell_connection_graph,
1787  Triangulation<dim, spacedim> &triangulation,
1788  const SparsityTools::Partitioner partitioner =
1790 
1801  template <int dim, int spacedim>
1802  void
1803  partition_triangulation(const unsigned int n_partitions,
1804  const std::vector<unsigned int> &cell_weights,
1805  const SparsityPattern & cell_connection_graph,
1806  Triangulation<dim, spacedim> &triangulation,
1807  const SparsityTools::Partitioner partitioner =
1809 
1817  template <int dim, int spacedim>
1818  void
1819  partition_triangulation_zorder(const unsigned int n_partitions,
1820  Triangulation<dim, spacedim> &triangulation);
1821 
1833  template <int dim, int spacedim>
1834  void
1836 
1847  template <int dim, int spacedim>
1848  void
1849  get_subdomain_association(const Triangulation<dim, spacedim> &triangulation,
1850  std::vector<types::subdomain_id> & subdomain);
1851 
1866  template <int dim, int spacedim>
1867  unsigned int
1869  const Triangulation<dim, spacedim> &triangulation,
1870  const types::subdomain_id subdomain);
1871 
1872 
1902  template <int dim, int spacedim>
1903  std::vector<bool>
1905 
1911 
1940  template <typename MeshType>
1941  std::list<std::pair<typename MeshType::cell_iterator,
1942  typename MeshType::cell_iterator>>
1943  get_finest_common_cells(const MeshType &mesh_1, const MeshType &mesh_2);
1944 
1954  template <int dim, int spacedim>
1955  bool
1957  const Triangulation<dim, spacedim> &mesh_2);
1958 
1968  template <typename MeshType>
1969  bool
1970  have_same_coarse_mesh(const MeshType &mesh_1, const MeshType &mesh_2);
1971 
1977 
1993  template <int dim, int spacedim>
1997  & distorted_cells,
1998  Triangulation<dim, spacedim> &triangulation);
1999 
2000 
2001 
2008 
2009 
2049  template <class MeshType>
2050  std::vector<typename MeshType::active_cell_iterator>
2051  get_patch_around_cell(const typename MeshType::active_cell_iterator &cell);
2052 
2053 
2077  template <class Container>
2078  std::vector<typename Container::cell_iterator>
2080  const std::vector<typename Container::active_cell_iterator> &patch_cells);
2081 
2150  template <class Container>
2151  void
2153  const std::vector<typename Container::active_cell_iterator> &patch,
2155  &local_triangulation,
2156  std::map<
2157  typename Triangulation<Container::dimension,
2158  Container::space_dimension>::active_cell_iterator,
2159  typename Container::active_cell_iterator> &patch_to_global_tria_map);
2160 
2197  template <class DoFHandlerType>
2198  std::map<types::global_dof_index,
2199  std::vector<typename DoFHandlerType::active_cell_iterator>>
2200  get_dof_to_support_patch_map(DoFHandlerType &dof_handler);
2201 
2202 
2209 
2215  template <typename CellIterator>
2216  struct PeriodicFacePair
2217  {
2221  CellIterator cell[2];
2222 
2227  unsigned int face_idx[2];
2228 
2234  std::bitset<3> orientation;
2235 
2249  };
2250 
2251 
2317  template <typename FaceIterator>
2318  bool orthogonal_equality(
2319  std::bitset<3> & orientation,
2320  const FaceIterator & face1,
2321  const FaceIterator & face2,
2322  const int direction,
2325  const FullMatrix<double> &matrix = FullMatrix<double>());
2326 
2327 
2331  template <typename FaceIterator>
2332  bool
2334  const FaceIterator & face1,
2335  const FaceIterator & face2,
2336  const int direction,
2339  const FullMatrix<double> &matrix = FullMatrix<double>());
2340 
2341 
2400  template <typename MeshType>
2401  void
2403  const MeshType & mesh,
2404  const types::boundary_id b_id1,
2405  const types::boundary_id b_id2,
2406  const int direction,
2408  & matched_pairs,
2409  const Tensor<1, MeshType::space_dimension> &offset =
2411  const FullMatrix<double> &matrix = FullMatrix<double>());
2412 
2413 
2436  template <typename MeshType>
2437  void
2439  const MeshType & mesh,
2440  const types::boundary_id b_id,
2441  const int direction,
2443  & matched_pairs,
2444  const ::Tensor<1, MeshType::space_dimension> &offset =
2446  const FullMatrix<double> &matrix = FullMatrix<double>());
2447 
2453 
2476  template <int dim, int spacedim>
2477  void
2479  const bool reset_boundary_ids = false);
2480 
2504  template <int dim, int spacedim>
2505  void
2507  const std::vector<types::boundary_id> &src_boundary_ids,
2508  const std::vector<types::manifold_id> &dst_manifold_ids,
2510  const std::vector<types::boundary_id> &reset_boundary_ids = {});
2511 
2543  template <int dim, int spacedim>
2544  void
2546  const bool compute_face_ids = false);
2547 
2548 
2633  template <typename DataType, typename MeshType>
2634  void
2636  const MeshType & mesh,
2637  const std::function<boost::optional<DataType>(
2638  const typename MeshType::active_cell_iterator &)> &pack,
2639  const std::function<void(const typename MeshType::active_cell_iterator &,
2640  const DataType &)> & unpack);
2641 
2642  /* Exchange with all processors of the MPI communicator @p mpi_communicator the vector of bounding
2643  * boxes @p local_bboxes.
2644  *
2645  * This function is meant to exchange bounding boxes describing the locally
2646  * owned cells in a distributed triangulation obtained with the function
2647  * GridTools::compute_mesh_predicate_bounding_box .
2648  *
2649  * The output vector's size is the number of processes of the MPI
2650  * communicator:
2651  * its i-th entry contains the vector @p local_bboxes of the i-th process.
2652  */
2653  template <int spacedim>
2654  std::vector<std::vector<BoundingBox<spacedim>>>
2655  exchange_local_bounding_boxes(
2656  const std::vector<BoundingBox<spacedim>> &local_bboxes,
2657  MPI_Comm mpi_communicator);
2658 
2671  template <int dim, typename T>
2673  {
2677  std::vector<CellId> cell_ids;
2678 
2682  std::vector<T> data;
2683 
2691  template <class Archive>
2692  void
2693  save(Archive &ar, const unsigned int version) const;
2694 
2699  template <class Archive>
2700  void
2701  load(Archive &ar, const unsigned int version);
2702 
2703  BOOST_SERIALIZATION_SPLIT_MEMBER()
2704  };
2705 
2710 
2715  int,
2716  << "The number of partitions you gave is " << arg1
2717  << ", but must be greater than zero.");
2722  int,
2723  << "The subdomain id " << arg1
2724  << " has no cells associated with it.");
2729 
2734  double,
2735  << "The scaling factor must be positive, but it is " << arg1
2736  << ".");
2740  template <int N>
2742  Point<N>,
2743  << "The point <" << arg1
2744  << "> could not be found inside any of the "
2745  << "coarse grid cells.");
2749  template <int N>
2751  Point<N>,
2752  << "The point <" << arg1
2753  << "> could not be found inside any of the "
2754  << "subcells of a coarse grid cell.");
2755 
2760  unsigned int,
2761  << "The given vertex with index " << arg1
2762  << " is not used in the given triangulation.");
2763 
2764 
2767 } /*namespace GridTools*/
2768 
2769 
2770 
2771 /* ----------------- Template function --------------- */
2772 
2773 # ifndef DOXYGEN
2774 
2775 namespace GridTools
2776 {
2777  template <int dim, typename T>
2778  double
2779  cell_measure(const T &, ...)
2780  {
2781  Assert(false, ExcNotImplemented());
2782  return std::numeric_limits<double>::quiet_NaN();
2783  }
2784 
2785  template <int dim, typename Predicate, int spacedim>
2786  void
2787  transform(const Predicate & predicate,
2788  Triangulation<dim, spacedim> &triangulation)
2789  {
2790  std::vector<bool> treated_vertices(triangulation.n_vertices(), false);
2791 
2792  // loop over all active cells, and
2793  // transform those vertices that
2794  // have not yet been touched. note
2795  // that we get to all vertices in
2796  // the triangulation by only
2797  // visiting the active cells.
2799  cell = triangulation.begin_active(),
2800  endc = triangulation.end();
2801  for (; cell != endc; ++cell)
2802  for (unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell; ++v)
2803  if (treated_vertices[cell->vertex_index(v)] == false)
2804  {
2805  // transform this vertex
2806  cell->vertex(v) = predicate(cell->vertex(v));
2807  // and mark it as treated
2808  treated_vertices[cell->vertex_index(v)] = true;
2809  };
2810 
2811 
2812  // now fix any vertices on hanging nodes so that we don't create any holes
2813  if (dim == 2)
2814  {
2816  cell = triangulation.begin_active(),
2817  endc = triangulation.end();
2818  for (; cell != endc; ++cell)
2819  for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell;
2820  ++face)
2821  if (cell->face(face)->has_children() &&
2822  !cell->face(face)->at_boundary())
2823  {
2824  // this line has children
2825  cell->face(face)->child(0)->vertex(1) =
2826  (cell->face(face)->vertex(0) + cell->face(face)->vertex(1)) /
2827  2;
2828  }
2829  }
2830  else if (dim == 3)
2831  {
2833  cell = triangulation.begin_active(),
2834  endc = triangulation.end();
2835  for (; cell != endc; ++cell)
2836  for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell;
2837  ++face)
2838  if (cell->face(face)->has_children() &&
2839  !cell->face(face)->at_boundary())
2840  {
2841  // this face has hanging nodes
2842  cell->face(face)->child(0)->vertex(1) =
2843  (cell->face(face)->vertex(0) + cell->face(face)->vertex(1)) /
2844  2.0;
2845  cell->face(face)->child(0)->vertex(2) =
2846  (cell->face(face)->vertex(0) + cell->face(face)->vertex(2)) /
2847  2.0;
2848  cell->face(face)->child(1)->vertex(3) =
2849  (cell->face(face)->vertex(1) + cell->face(face)->vertex(3)) /
2850  2.0;
2851  cell->face(face)->child(2)->vertex(3) =
2852  (cell->face(face)->vertex(2) + cell->face(face)->vertex(3)) /
2853  2.0;
2854 
2855  // center of the face
2856  cell->face(face)->child(0)->vertex(3) =
2857  (cell->face(face)->vertex(0) + cell->face(face)->vertex(1) +
2858  cell->face(face)->vertex(2) + cell->face(face)->vertex(3)) /
2859  4.0;
2860  }
2861  }
2862 
2863  // Make sure FEValues notices that the mesh has changed
2864  triangulation.signals.mesh_movement();
2865  }
2866 
2867 
2868 
2869  template <class MeshType>
2870  std::vector<typename MeshType::active_cell_iterator>
2871  get_active_child_cells(const typename MeshType::cell_iterator &cell)
2872  {
2873  std::vector<typename MeshType::active_cell_iterator> child_cells;
2874 
2875  if (cell->has_children())
2876  {
2877  for (unsigned int child = 0; child < cell->n_children(); ++child)
2878  if (cell->child(child)->has_children())
2879  {
2880  const std::vector<typename MeshType::active_cell_iterator>
2881  children = get_active_child_cells<MeshType>(cell->child(child));
2882  child_cells.insert(child_cells.end(),
2883  children.begin(),
2884  children.end());
2885  }
2886  else
2887  child_cells.push_back(cell->child(child));
2888  }
2889 
2890  return child_cells;
2891  }
2892 
2893 
2894 
2895  template <class MeshType>
2896  void
2898  const typename MeshType::active_cell_iterator & cell,
2899  std::vector<typename MeshType::active_cell_iterator> &active_neighbors)
2900  {
2901  active_neighbors.clear();
2902  for (unsigned int n = 0;
2903  n < GeometryInfo<MeshType::dimension>::faces_per_cell;
2904  ++n)
2905  if (!cell->at_boundary(n))
2906  {
2907  if (MeshType::dimension == 1)
2908  {
2909  // check children of neighbor. note
2910  // that in 1d children of the neighbor
2911  // may be further refined. In 1d the
2912  // case is simple since we know what
2913  // children bound to the present cell
2914  typename MeshType::cell_iterator neighbor_child =
2915  cell->neighbor(n);
2916  if (!neighbor_child->active())
2917  {
2918  while (neighbor_child->has_children())
2919  neighbor_child = neighbor_child->child(n == 0 ? 1 : 0);
2920 
2921  Assert(neighbor_child->neighbor(n == 0 ? 1 : 0) == cell,
2922  ExcInternalError());
2923  }
2924  active_neighbors.push_back(neighbor_child);
2925  }
2926  else
2927  {
2928  if (cell->face(n)->has_children())
2929  // this neighbor has children. find
2930  // out which border to the present
2931  // cell
2932  for (unsigned int c = 0;
2933  c < cell->face(n)->number_of_children();
2934  ++c)
2935  active_neighbors.push_back(
2936  cell->neighbor_child_on_subface(n, c));
2937  else
2938  {
2939  // the neighbor must be active
2940  // himself
2941  Assert(cell->neighbor(n)->active(), ExcInternalError());
2942  active_neighbors.push_back(cell->neighbor(n));
2943  }
2944  }
2945  }
2946  }
2947 
2948 
2949 
2950  namespace internal
2951  {
2952  namespace ProjectToObject
2953  {
2966  struct CrossDerivative
2967  {
2968  const unsigned int direction_0;
2969  const unsigned int direction_1;
2970 
2971  CrossDerivative(const unsigned int d0, const unsigned int d1);
2972  };
2973 
2974  inline CrossDerivative::CrossDerivative(const unsigned int d0,
2975  const unsigned int d1)
2976  : direction_0(d0)
2977  , direction_1(d1)
2978  {}
2979 
2980 
2981 
2986  template <typename F>
2987  inline auto
2988  centered_first_difference(const double center,
2989  const double step,
2990  const F &f) -> decltype(f(center) - f(center))
2991  {
2992  return (f(center + step) - f(center - step)) / (2.0 * step);
2993  }
2994 
2995 
2996 
3001  template <typename F>
3002  inline auto
3003  centered_second_difference(const double center,
3004  const double step,
3005  const F &f) -> decltype(f(center) - f(center))
3006  {
3007  return (f(center + step) - 2.0 * f(center) + f(center - step)) /
3008  (step * step);
3009  }
3010 
3011 
3012 
3022  template <int structdim, typename F>
3023  inline auto
3024  cross_stencil(
3025  const CrossDerivative cross_derivative,
3027  const double step,
3028  const F &f) -> decltype(f(center) - f(center))
3029  {
3031  simplex_vector[cross_derivative.direction_0] = 0.5 * step;
3032  simplex_vector[cross_derivative.direction_1] = -0.5 * step;
3033  return (-4.0 * f(center) - 1.0 * f(center + simplex_vector) -
3034  1.0 / 3.0 * f(center - simplex_vector) +
3035  16.0 / 3.0 * f(center + 0.5 * simplex_vector)) /
3036  step;
3037  }
3038 
3039 
3040 
3047  template <int spacedim, int structdim, typename F>
3048  inline double
3049  gradient_entry(
3050  const unsigned int row_n,
3051  const unsigned int dependent_direction,
3052  const Point<spacedim> &p0,
3054  const double step,
3055  const F & f)
3056  {
3058  dependent_direction <
3060  ExcMessage("This function assumes that the last weight is a "
3061  "dependent variable (and hence we cannot take its "
3062  "derivative directly)."));
3063  Assert(row_n != dependent_direction,
3064  ExcMessage(
3065  "We cannot differentiate with respect to the variable "
3066  "that is assumed to be dependent."));
3067 
3068  const Point<spacedim> manifold_point = f(center);
3069  const Tensor<1, spacedim> stencil_value = cross_stencil<structdim>(
3070  {row_n, dependent_direction}, center, step, f);
3071  double entry = 0.0;
3072  for (unsigned int dim_n = 0; dim_n < spacedim; ++dim_n)
3073  entry +=
3074  -2.0 * (p0[dim_n] - manifold_point[dim_n]) * stencil_value[dim_n];
3075  return entry;
3076  }
3077 
3083  template <typename Iterator, int spacedim, int structdim>
3085  project_to_d_linear_object(const Iterator & object,
3086  const Point<spacedim> &trial_point)
3087  {
3088  // let's look at this for simplicity for a quad (structdim==2) in a
3089  // space with spacedim>2 (notate trial_point by y): all points on the
3090  // surface are given by
3091  // x(\xi) = sum_i v_i phi_x(\xi)
3092  // where v_i are the vertices of the quad, and \xi=(\xi_1,\xi_2) are the
3093  // reference coordinates of the quad. so what we are trying to do is
3094  // find a point x on the surface that is closest to the point y. there
3095  // are different ways to solve this problem, but in the end it's a
3096  // nonlinear problem and we have to find reference coordinates \xi so
3097  // that J(\xi) = 1/2 || x(\xi)-y ||^2 is minimal. x(\xi) is a function
3098  // that is structdim-linear in \xi, so J(\xi) is a polynomial of degree
3099  // 2*structdim that we'd like to minimize. unless structdim==1, we'll
3100  // have to use a Newton method to find the answer. This leads to the
3101  // following formulation of Newton steps:
3102  //
3103  // Given \xi_k, find \delta\xi_k so that
3104  // H_k \delta\xi_k = - F_k
3105  // where H_k is an approximation to the second derivatives of J at
3106  // \xi_k, and F_k is the first derivative of J. We'll iterate this a
3107  // number of times until the right hand side is small enough. As a
3108  // stopping criterion, we terminate if ||\delta\xi||<eps.
3109  //
3110  // As for the Hessian, the best choice would be
3111  // H_k = J''(\xi_k)
3112  // but we'll opt for the simpler Gauss-Newton form
3113  // H_k = A^T A
3114  // i.e.
3115  // (H_k)_{nm} = \sum_{i,j} v_i*v_j *
3116  // \partial_n phi_i *
3117  // \partial_m phi_j
3118  // we start at xi=(0.5, 0.5).
3119  Point<structdim> xi;
3120  for (unsigned int d = 0; d < structdim; ++d)
3121  xi[d] = 0.5;
3122 
3123  Point<spacedim> x_k;
3124  for (unsigned int i = 0; i < GeometryInfo<structdim>::vertices_per_cell;
3125  ++i)
3126  x_k += object->vertex(i) *
3128 
3129  do
3130  {
3132  for (unsigned int i = 0;
3133  i < GeometryInfo<structdim>::vertices_per_cell;
3134  ++i)
3135  F_k +=
3136  (x_k - trial_point) * object->vertex(i) *
3138  i);
3139 
3141  for (unsigned int i = 0;
3142  i < GeometryInfo<structdim>::vertices_per_cell;
3143  ++i)
3144  for (unsigned int j = 0;
3145  j < GeometryInfo<structdim>::vertices_per_cell;
3146  ++j)
3147  {
3150  xi, i),
3152  xi, j));
3153  H_k += (object->vertex(i) * object->vertex(j)) * tmp;
3154  }
3155 
3156  const Tensor<1, structdim> delta_xi = -invert(H_k) * F_k;
3157  xi += delta_xi;
3158 
3159  x_k = Point<spacedim>();
3160  for (unsigned int i = 0;
3161  i < GeometryInfo<structdim>::vertices_per_cell;
3162  ++i)
3163  x_k += object->vertex(i) *
3165 
3166  if (delta_xi.norm() < 1e-7)
3167  break;
3168  }
3169  while (true);
3170 
3171  return x_k;
3172  }
3173  } // namespace ProjectToObject
3174  } // namespace internal
3175 
3176 
3177 
3178  namespace internal
3179  {
3180  // We hit an internal compiler error in ICC 15 if we define this as a lambda
3181  // inside the project_to_object function below.
3182  template <int structdim>
3183  inline bool
3184  weights_are_ok(
3186  {
3187  // clang has trouble figuring out structdim here, so define it
3188  // again:
3189  static const std::size_t n_vertices_per_cell =
3191  n_independent_components;
3192  std::array<double, n_vertices_per_cell> copied_weights;
3193  for (unsigned int i = 0; i < n_vertices_per_cell; ++i)
3194  {
3195  copied_weights[i] = v[i];
3196  if (v[i] < 0.0 || v[i] > 1.0)
3197  return false;
3198  }
3199 
3200  // check the sum: try to avoid some roundoff errors by summing in order
3201  std::sort(copied_weights.begin(), copied_weights.end());
3202  const double sum =
3203  std::accumulate(copied_weights.begin(), copied_weights.end(), 0.0);
3204  return std::abs(sum - 1.0) < 1e-10; // same tolerance used in manifold.cc
3205  }
3206  } // namespace internal
3207 
3208  template <typename Iterator>
3211  const Iterator & object,
3213  {
3214  const int spacedim = Iterator::AccessorType::space_dimension;
3215  const int structdim = Iterator::AccessorType::structure_dimension;
3216 
3217  Point<spacedim> projected_point = trial_point;
3218 
3219  if (structdim >= spacedim)
3220  return projected_point;
3221  else if (structdim == 1 || structdim == 2)
3222  {
3223  using namespace internal::ProjectToObject;
3224  // Try to use the special flat algorithm for quads (this is better
3225  // than the general algorithm in 3D). This does not take into account
3226  // whether projected_point is outside the quad, but we optimize along
3227  // lines below anyway:
3228  const int dim = Iterator::AccessorType::dimension;
3229  const Manifold<dim, spacedim> &manifold = object->get_manifold();
3230  if (structdim == 2 && dynamic_cast<const FlatManifold<dim, spacedim> *>(
3231  &manifold) != nullptr)
3232  {
3233  projected_point =
3234  project_to_d_linear_object<Iterator, spacedim, structdim>(
3235  object, trial_point);
3236  }
3237  else
3238  {
3239  // We want to find a point on the convex hull (defined by the
3240  // vertices of the object and the manifold description) that is
3241  // relatively close to the trial point. This has a few issues:
3242  //
3243  // 1. For a general convex hull we are not guaranteed that a unique
3244  // minimum exists.
3245  // 2. The independent variables in the optimization process are the
3246  // weights given to Manifold::get_new_point, which must sum to 1,
3247  // so we cannot use standard finite differences to approximate a
3248  // gradient.
3249  //
3250  // There is not much we can do about 1., but for 2. we can derive
3251  // finite difference stencils that work on a structdim-dimensional
3252  // simplex and rewrite the optimization problem to use those
3253  // instead. Consider the structdim 2 case and let
3254  //
3255  // F(c0, c1, c2, c3) = Manifold::get_new_point(vertices, {c0, c1,
3256  // c2, c3})
3257  //
3258  // where {c0, c1, c2, c3} are the weights for the four vertices on
3259  // the quadrilateral. We seek to minimize the Euclidean distance
3260  // between F(...) and trial_point. We can solve for c3 in terms of
3261  // the other weights and get, for one coordinate direction
3262  //
3263  // d/dc0 ((x0 - F(c0, c1, c2, 1 - c0 - c1 - c2))^2)
3264  // = -2(x0 - F(...)) (d/dc0 F(...) - d/dc3 F(...))
3265  //
3266  // where we substitute back in for c3 after taking the
3267  // derivative. We can compute a stencil for the cross derivative
3268  // d/dc0 - d/dc3: this is exactly what cross_stencil approximates
3269  // (and gradient_entry computes the sum over the independent
3270  // variables). Below, we somewhat arbitrarily pick the last
3271  // component as the dependent one.
3272  //
3273  // Since we can now calculate derivatives of the objective
3274  // function we can use gradient descent to minimize it.
3275  //
3276  // Of course, this is much simpler in the structdim = 1 case (we
3277  // could rewrite the projection as a 1D optimization problem), but
3278  // to reduce the potential for bugs we use the same code in both
3279  // cases.
3280  const double step_size = object->diameter() / 64.0;
3281 
3282  constexpr unsigned int n_vertices_per_cell =
3284 
3285  std::array<Point<spacedim>, n_vertices_per_cell> vertices;
3286  for (unsigned int vertex_n = 0; vertex_n < n_vertices_per_cell;
3287  ++vertex_n)
3288  vertices[vertex_n] = object->vertex(vertex_n);
3289 
3290  auto get_point_from_weights =
3291  [&](const Tensor<1, n_vertices_per_cell> &weights)
3292  -> Point<spacedim> {
3293  return object->get_manifold().get_new_point(
3294  make_array_view(vertices.begin(), vertices.end()),
3295  make_array_view(&weights[0],
3296  &weights[n_vertices_per_cell - 1] + 1));
3297  };
3298 
3299  // pick the initial weights as (normalized) inverse distances from
3300  // the trial point:
3301  Tensor<1, n_vertices_per_cell> guess_weights;
3302  double guess_weights_sum = 0.0;
3303  for (unsigned int vertex_n = 0; vertex_n < n_vertices_per_cell;
3304  ++vertex_n)
3305  {
3306  const double distance =
3307  vertices[vertex_n].distance(trial_point);
3308  if (distance == 0.0)
3309  {
3310  guess_weights = 0.0;
3311  guess_weights[vertex_n] = 1.0;
3312  guess_weights_sum = 1.0;
3313  break;
3314  }
3315  else
3316  {
3317  guess_weights[vertex_n] = 1.0 / distance;
3318  guess_weights_sum += guess_weights[vertex_n];
3319  }
3320  }
3321  guess_weights /= guess_weights_sum;
3322  Assert(internal::weights_are_ok<structdim>(guess_weights),
3323  ExcInternalError());
3324 
3325  // The optimization algorithm consists of two parts:
3326  //
3327  // 1. An outer loop where we apply the gradient descent algorithm.
3328  // 2. An inner loop where we do a line search to find the optimal
3329  // length of the step one should take in the gradient direction.
3330  //
3331  for (unsigned int outer_n = 0; outer_n < 40; ++outer_n)
3332  {
3333  const unsigned int dependent_direction =
3334  n_vertices_per_cell - 1;
3335  Tensor<1, n_vertices_per_cell> current_gradient;
3336  for (unsigned int row_n = 0; row_n < n_vertices_per_cell;
3337  ++row_n)
3338  {
3339  if (row_n != dependent_direction)
3340  {
3341  current_gradient[row_n] =
3342  gradient_entry<spacedim, structdim>(
3343  row_n,
3344  dependent_direction,
3345  trial_point,
3346  guess_weights,
3347  step_size,
3348  get_point_from_weights);
3349 
3350  current_gradient[dependent_direction] -=
3351  current_gradient[row_n];
3352  }
3353  }
3354 
3355  // We need to travel in the -gradient direction, as noted
3356  // above, but we may not want to take a full step in that
3357  // direction; instead, guess that we will go -0.5*gradient and
3358  // do quasi-Newton iteration to pick the best multiplier. The
3359  // goal is to find a scalar alpha such that
3360  //
3361  // F(x - alpha g)
3362  //
3363  // is minimized, where g is the gradient and F is the
3364  // objective function. To find the optimal value we find roots
3365  // of the derivative of the objective function with respect to
3366  // alpha by Newton iteration, where we approximate the first
3367  // and second derivatives of F(x - alpha g) with centered
3368  // finite differences.
3369  double gradient_weight = -0.5;
3370  auto gradient_weight_objective_function =
3371  [&](const double gradient_weight_guess) -> double {
3372  return (trial_point -
3373  get_point_from_weights(guess_weights +
3374  gradient_weight_guess *
3375  current_gradient))
3376  .norm_square();
3377  };
3378 
3379  for (unsigned int inner_n = 0; inner_n < 10; ++inner_n)
3380  {
3381  const double update_numerator = centered_first_difference(
3382  gradient_weight,
3383  step_size,
3384  gradient_weight_objective_function);
3385  const double update_denominator =
3386  centered_second_difference(
3387  gradient_weight,
3388  step_size,
3389  gradient_weight_objective_function);
3390 
3391  // avoid division by zero. Note that we limit the gradient
3392  // weight below
3393  if (std::abs(update_denominator) == 0.0)
3394  break;
3395  gradient_weight =
3396  gradient_weight - update_numerator / update_denominator;
3397 
3398  // Put a fairly lenient bound on the largest possible
3399  // gradient (things tend to be locally flat, so the gradient
3400  // itself is usually small)
3401  if (std::abs(gradient_weight) > 10)
3402  {
3403  gradient_weight = -10.0;
3404  break;
3405  }
3406  }
3407 
3408  // It only makes sense to take convex combinations with weights
3409  // between zero and one. If the update takes us outside of this
3410  // region then rescale the update to stay within the region and
3411  // try again
3412  Tensor<1, n_vertices_per_cell> tentative_weights =
3413  guess_weights + gradient_weight * current_gradient;
3414 
3415  double new_gradient_weight = gradient_weight;
3416  for (unsigned int iteration_count = 0; iteration_count < 40;
3417  ++iteration_count)
3418  {
3419  if (internal::weights_are_ok<structdim>(tentative_weights))
3420  break;
3421 
3422  for (unsigned int i = 0; i < n_vertices_per_cell; ++i)
3423  {
3424  if (tentative_weights[i] < 0.0)
3425  {
3426  tentative_weights -=
3427  (tentative_weights[i] / current_gradient[i]) *
3428  current_gradient;
3429  }
3430  if (tentative_weights[i] < 0.0 ||
3431  1.0 < tentative_weights[i])
3432  {
3433  new_gradient_weight /= 2.0;
3434  tentative_weights =
3435  guess_weights +
3436  new_gradient_weight * current_gradient;
3437  }
3438  }
3439  }
3440 
3441  // the update might still send us outside the valid region, so
3442  // check again and quit if the update is still not valid
3443  if (!internal::weights_are_ok<structdim>(tentative_weights))
3444  break;
3445 
3446  // if we cannot get closer by traveling in the gradient
3447  // direction then quit
3448  if (get_point_from_weights(tentative_weights)
3449  .distance(trial_point) <
3450  get_point_from_weights(guess_weights).distance(trial_point))
3451  guess_weights = tentative_weights;
3452  else
3453  break;
3454  Assert(internal::weights_are_ok<structdim>(guess_weights),
3455  ExcInternalError());
3456  }
3457  Assert(internal::weights_are_ok<structdim>(guess_weights),
3458  ExcInternalError());
3459  projected_point = get_point_from_weights(guess_weights);
3460  }
3461 
3462  // if structdim == 2 and the optimal point is not on the interior then
3463  // we may be able to get a more accurate result by projecting onto the
3464  // lines.
3465  if (structdim == 2)
3466  {
3467  std::array<Point<spacedim>, GeometryInfo<structdim>::lines_per_cell>
3468  line_projections;
3469  for (unsigned int line_n = 0;
3470  line_n < GeometryInfo<structdim>::lines_per_cell;
3471  ++line_n)
3472  {
3473  line_projections[line_n] =
3474  project_to_object(object->line(line_n), trial_point);
3475  }
3476  std::sort(line_projections.begin(),
3477  line_projections.end(),
3478  [&](const Point<spacedim> &a, const Point<spacedim> &b) {
3479  return a.distance(trial_point) <
3480  b.distance(trial_point);
3481  });
3482  if (line_projections[0].distance(trial_point) <
3483  projected_point.distance(trial_point))
3484  projected_point = line_projections[0];
3485  }
3486  }
3487  else
3488  {
3489  Assert(false, ExcNotImplemented());
3490  return projected_point;
3491  }
3492 
3493  return projected_point;
3494  }
3495 
3496 
3497 
3498  template <int dim, typename T>
3499  template <class Archive>
3500  void
3502  const unsigned int /*version*/) const
3503  {
3504  Assert(cell_ids.size() == data.size(),
3505  ExcDimensionMismatch(cell_ids.size(), data.size()));
3506  // archive the cellids in an efficient binary format
3507  const size_t n_cells = cell_ids.size();
3508  ar & n_cells;
3509  for (auto &it : cell_ids)
3510  {
3511  CellId::binary_type binary_cell_id = it.template to_binary<dim>();
3512  ar & binary_cell_id;
3513  }
3514 
3515  ar &data;
3516  }
3517 
3518 
3519 
3520  template <int dim, typename T>
3521  template <class Archive>
3522  void
3524  const unsigned int /*version*/)
3525  {
3526  size_t n_cells;
3527  ar & n_cells;
3528  cell_ids.clear();
3529  cell_ids.reserve(n_cells);
3530  for (unsigned int c = 0; c < n_cells; ++c)
3531  {
3532  CellId::binary_type value;
3533  ar & value;
3534  cell_ids.emplace_back(value);
3535  }
3536  ar &data;
3537  }
3538 
3539 
3540 
3541  template <typename DataType, typename MeshType>
3542  void
3544  const MeshType & mesh,
3545  const std::function<boost::optional<DataType>(
3546  const typename MeshType::active_cell_iterator &)> &pack,
3547  const std::function<void(const typename MeshType::active_cell_iterator &,
3548  const DataType &)> & unpack)
3549  {
3550 # ifndef DEAL_II_WITH_MPI
3551  (void)mesh;
3552  (void)pack;
3553  (void)unpack;
3554  Assert(false,
3555  ExcMessage(
3556  "GridTools::exchange_cell_data_to_ghosts() requires MPI."));
3557 # else
3558  constexpr int dim = MeshType::dimension;
3559  constexpr int spacedim = MeshType::space_dimension;
3560  auto tria = static_cast<const parallel::Triangulation<dim, spacedim> *>(
3561  &mesh.get_triangulation());
3562  Assert(
3563  tria != nullptr,
3564  ExcMessage(
3565  "The function exchange_cell_data_to_ghosts() only works with parallel triangulations."));
3566 
3567  // map neighbor_id -> data_buffer where we accumulate the data to send
3568  using DestinationToBufferMap =
3569  std::map<::types::subdomain_id,
3571  DestinationToBufferMap destination_to_data_buffer_map;
3572 
3573  std::map<unsigned int, std::set<::types::subdomain_id>>
3574  vertices_with_ghost_neighbors =
3575  tria->compute_vertices_with_ghost_neighbors();
3576 
3577  for (auto cell : tria->active_cell_iterators())
3578  if (cell->is_locally_owned())
3579  {
3580  std::set<::types::subdomain_id> send_to;
3581  for (unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell;
3582  ++v)
3583  {
3584  const std::map<unsigned int,
3585  std::set<::types::subdomain_id>>::
3586  const_iterator neighbor_subdomains_of_vertex =
3587  vertices_with_ghost_neighbors.find(cell->vertex_index(v));
3588 
3589  if (neighbor_subdomains_of_vertex ==
3590  vertices_with_ghost_neighbors.end())
3591  continue;
3592 
3593  Assert(neighbor_subdomains_of_vertex->second.size() != 0,
3594  ExcInternalError());
3595 
3596  send_to.insert(neighbor_subdomains_of_vertex->second.begin(),
3597  neighbor_subdomains_of_vertex->second.end());
3598  }
3599 
3600  if (send_to.size() > 0)
3601  {
3602  // this cell's data needs to be sent to someone
3603  typename MeshType::active_cell_iterator mesh_it(tria,
3604  cell->level(),
3605  cell->index(),
3606  &mesh);
3607 
3608  const boost::optional<DataType> data = pack(mesh_it);
3609 
3610  if (data)
3611  {
3612  const CellId cellid = cell->id();
3613 
3614  for (auto it : send_to)
3615  {
3616  const ::types::subdomain_id subdomain = it;
3617 
3618  // find the data buffer for proc "subdomain" if it exists
3619  // or create an empty one otherwise
3620  typename DestinationToBufferMap::iterator p =
3621  destination_to_data_buffer_map
3622  .insert(std::make_pair(
3623  subdomain, CellDataTransferBuffer<dim, DataType>()))
3624  .first;
3625 
3626  p->second.cell_ids.emplace_back(cellid);
3627  p->second.data.emplace_back(data.get());
3628  }
3629  }
3630  }
3631  }
3632 
3633 
3634  // 2. send our messages
3635  std::set<::types::subdomain_id> ghost_owners = tria->ghost_owners();
3636  const unsigned int n_ghost_owners = ghost_owners.size();
3637  std::vector<std::vector<char>> sendbuffers(n_ghost_owners);
3638  std::vector<MPI_Request> requests(n_ghost_owners);
3639 
3640  unsigned int idx = 0;
3641  for (auto it = ghost_owners.begin(); it != ghost_owners.end(); ++it, ++idx)
3642  {
3643  CellDataTransferBuffer<dim, DataType> &data =
3644  destination_to_data_buffer_map[*it];
3645 
3646  // pack all the data into the buffer for this recipient and send it.
3647  // keep data around till we can make sure that the packet has been
3648  // received
3649  sendbuffers[idx] = Utilities::pack(data);
3650  const int ierr = MPI_Isend(sendbuffers[idx].data(),
3651  sendbuffers[idx].size(),
3652  MPI_BYTE,
3653  *it,
3654  786,
3655  tria->get_communicator(),
3656  &requests[idx]);
3657  AssertThrowMPI(ierr);
3658  }
3659 
3660  // 3. receive messages
3661  std::vector<char> receive;
3662  for (unsigned int idx = 0; idx < n_ghost_owners; ++idx)
3663  {
3664  MPI_Status status;
3665  int len;
3666  int ierr =
3667  MPI_Probe(MPI_ANY_SOURCE, 786, tria->get_communicator(), &status);
3668  AssertThrowMPI(ierr);
3669  ierr = MPI_Get_count(&status, MPI_BYTE, &len);
3670  AssertThrowMPI(ierr);
3671 
3672  receive.resize(len);
3673 
3674  char *ptr = receive.data();
3675  ierr = MPI_Recv(ptr,
3676  len,
3677  MPI_BYTE,
3678  status.MPI_SOURCE,
3679  status.MPI_TAG,
3680  tria->get_communicator(),
3681  &status);
3682  AssertThrowMPI(ierr);
3683 
3684  auto cellinfo =
3685  Utilities::unpack<CellDataTransferBuffer<dim, DataType>>(receive);
3686 
3687  DataType *data = cellinfo.data.data();
3688  for (unsigned int c = 0; c < cellinfo.cell_ids.size(); ++c, ++data)
3689  {
3691  tria_cell = cellinfo.cell_ids[c].to_cell(*tria);
3692 
3693  const typename MeshType::active_cell_iterator cell(
3694  tria, tria_cell->level(), tria_cell->index(), &mesh);
3695 
3696  unpack(cell, *data);
3697  }
3698  }
3699 
3700  // make sure that all communication is finished
3701  // when we leave this function.
3702  if (requests.size())
3703  {
3704  const int ierr =
3705  MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
3706  AssertThrowMPI(ierr);
3707  }
3708 # endif // DEAL_II_WITH_MPI
3709  }
3710 } // namespace GridTools
3711 
3712 # endif
3713 
3714 DEAL_II_NAMESPACE_CLOSE
3715 
3716 /*---------------------------- grid_tools.h ---------------------------*/
3717 /* end of #ifndef dealii_grid_tools_h */
3718 #endif
3719 /*---------------------------- grid_tools.h ---------------------------*/
void remove_hanging_nodes(Triangulation< dim, spacedim > &tria, const bool isotropic=false, const unsigned int max_iterations=100)
Definition: grid_tools.cc:3608
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:3475
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)
static::ExceptionBase & ExcScalingFactorNotPositive(double arg1)
static const unsigned int invalid_unsigned_int
Definition: types.h:173
std::map< unsigned int, Point< spacedim > > get_all_vertices_at_boundary(const Triangulation< dim, spacedim > &tria)
Definition: grid_tools.cc:875
void copy_boundary_to_manifold_id(Triangulation< dim, spacedim > &tria, const bool reset_boundary_ids=false)
Definition: grid_tools.cc:3450
return_type guess_point_owner(const std::vector< std::vector< BoundingBox< spacedim >>> &global_bboxes, const std::vector< Point< spacedim >> &points)
Definition: grid_tools.cc:1867
double diameter(const Triangulation< dim, spacedim > &tria)
Definition: grid_tools.cc:75
void distort_random(const double factor, Triangulation< dim, spacedim > &triangulation, const bool keep_boundary=true)
Definition: grid_tools.cc:908
cell_iterator end() const
Definition: tria.cc:12004
static::ExceptionBase & ExcPointNotFoundInCoarseGrid(Point< N > arg1)
std::map< unsigned int, Point< spacedim > > extract_used_vertices(const Triangulation< dim, spacedim > &container, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim, spacedim >::mapping)
Definition: grid_tools.cc:4740
bool have_same_coarse_mesh(const Triangulation< dim, spacedim > &mesh_1, const Triangulation< dim, spacedim > &mesh_2)
MeshType< dim, spacedim >::active_cell_iterator find_active_cell_around_point(const MeshType< dim, spacedim > &mesh, const Point< spacedim > &p, const std::vector< bool > &marked_vertices=std::vector< bool >())
std::pair< unsigned int, double > get_longest_direction(typename Triangulation< dim, spacedim >::active_cell_iterator cell)
Definition: grid_tools.cc:3576
BoundingBox< spacedim > compute_bounding_box(const Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:420
unsigned int find_closest_vertex_of_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell, const Point< spacedim > &position)
Definition: grid_tools.cc:1637
void load(Archive &ar, const unsigned int version)
std::vector< typename MeshType::active_cell_iterator > compute_ghost_cell_layer_within_distance(const MeshType &mesh, const double layer_thickness)
IteratorRange< active_cell_iterator > active_cell_iterators() const
Definition: tria.cc:12060
void regularize_corner_cells(Triangulation< dim, spacedim > &tria, const double limit_angle_fraction=.75)
Definition: grid_tools.cc:3678
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 >())
void scale(const double scaling_factor, Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:725
double volume(const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping=(StaticMappingQ1< dim, spacedim >::mapping))
Definition: grid_tools.cc:133
std::vector< std::set< typename Triangulation< dim, spacedim >::active_cell_iterator > > vertex_to_cell_map(const Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:1919
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=1e-12, const std::vector< bool > &marked_vertices=std::vector< bool >())
SymmetricTensor< 4, dim, Number > outer_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2)
static::ExceptionBase & ExcNonExistentSubdomain(int arg1)
std::vector< typename MeshType::active_cell_iterator > get_active_child_cells(const typename MeshType::cell_iterator &cell)
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)
Definition: grid_tools.cc:4506
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:1970
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
static Tensor< 1, dim > d_linear_shape_function_gradient(const Point< dim > &xi, const unsigned int i)
numbers::NumberTraits< Number >::real_type norm() const
Definition: tensor.h:1301
std::array< unsigned int, 4 > binary_type
Definition: cell_id.h:73
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)
double maximal_cell_diameter(const Triangulation< dim, spacedim > &triangulation, const Mapping< dim, spacedim > &mapping=(StaticMappingQ1< dim, spacedim >::mapping))
Definition: grid_tools.cc:2916
std::bitset< 3 > orientation
Definition: grid_tools.h:2234
void partition_multigrid_levels(Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:2763
double minimal_cell_diameter(const Triangulation< dim, spacedim > &triangulation, const Mapping< dim, spacedim > &mapping=(StaticMappingQ1< dim, spacedim >::mapping))
Definition: grid_tools.cc:2888
void partition_triangulation_zorder(const unsigned int n_partitions, Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:2657
unsigned long long int global_dof_index
Definition: types.h:72
void partition_triangulation(const unsigned int n_partitions, Triangulation< dim, spacedim > &triangulation, const SparsityTools::Partitioner partitioner=SparsityTools::Partitioner::metis)
Definition: grid_tools.cc:2454
void delete_unused_vertices(std::vector< Point< spacedim >> &vertices, std::vector< CellData< dim >> &cells, SubCellData &subcelldata)
Definition: grid_tools.cc:434
SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
void get_vertex_connectivity_of_cells(const Triangulation< dim, spacedim > &triangulation, DynamicSparsityPattern &connectivity)
Definition: grid_tools.cc:2383
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)
static::ExceptionBase & ExcTriangulationHasBeenRefined()
Definition: tria.h:81
static::ExceptionBase & ExcInvalidNumberOfPartitions(int arg1)
static::ExceptionBase & ExcMessage(std::string arg1)
unsigned int subdomain_id
Definition: types.h:43
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:408
void get_vertex_connectivity_of_cells_on_level(const Triangulation< dim, spacedim > &triangulation, const unsigned int level, DynamicSparsityPattern &connectivity)
Definition: grid_tools.cc:2418
#define Assert(cond, exc)
Definition: exceptions.h:1227
Signals signals
Definition: tria.h:2287
static::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::vector< CellId > cell_ids
Definition: grid_tools.h:2677
Abstract base class for mapping classes.
Definition: dof_tools.h:57
void get_active_neighbors(const typename MeshType::active_cell_iterator &cell, std::vector< typename MeshType::active_cell_iterator > &active_neighbors)
#define DeclException0(Exception0)
Definition: exceptions.h:385
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:1720
void remove_anisotropy(Triangulation< dim, spacedim > &tria, const double max_ratio=1.6180339887, const unsigned int max_iterations=5)
Definition: grid_tools.cc:3645
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)
void copy_material_to_manifold_id(Triangulation< dim, spacedim > &tria, const bool compute_face_ids=false)
Definition: grid_tools.cc:3549
void rotate(const double angle, Triangulation< 2 > &triangulation)
Definition: grid_tools.cc:707
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:3411
size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
Definition: utilities.h:1046
Definition: hp.h:102
std::vector< typename MeshType::active_cell_iterator > compute_ghost_cell_halo_layer(const MeshType &mesh)
Definition: cell_id.h:63
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:1436
Point< Iterator::AccessorType::space_dimension > project_to_object(const Iterator &object, const Point< Iterator::AccessorType::space_dimension > &trial_point)
void transform(const InputIterator &begin_in, const InputIterator &end_in, OutputIterator out, Predicate &predicate, const unsigned int grainsize)
Definition: parallel.h:173
double cell_measure(const std::vector< Point< dim >> &all_vertices, const unsigned int(&vertex_indices)[GeometryInfo< dim >::vertices_per_cell])
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1443
std::map< types::global_dof_index, std::vector< typename DoFHandlerType::active_cell_iterator > > get_dof_to_support_patch_map(DoFHandlerType &dof_handler)
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
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:1298
void get_face_connectivity_of_cells(const Triangulation< dim, spacedim > &triangulation, DynamicSparsityPattern &connectivity)
Definition: grid_tools.cc:2327
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 >())
void save(Archive &ar, const unsigned int version) const
FullMatrix< double > matrix
Definition: grid_tools.h:2248
static::ExceptionBase & ExcNotImplemented()
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:3999
std::vector< bool > get_locally_owned_vertices(const Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:2826
void exchange_cell_data_to_ghosts(const MeshType &mesh, const std::function< boost::optional< DataType >(const typename MeshType::active_cell_iterator &)> &pack, const std::function< void(const typename MeshType::active_cell_iterator &, const DataType &)> &unpack)
std::vector< typename Container::cell_iterator > get_cells_at_coarsest_common_level(const std::vector< typename Container::active_cell_iterator > &patch_cells)
static::ExceptionBase & ExcVertexNotUsed(unsigned int arg1)
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: tria.cc:11935
std::list< std::pair< typename MeshType::cell_iterator, typename MeshType::cell_iterator > > get_finest_common_cells(const MeshType &mesh_1, const MeshType &mesh_2)
unsigned int n_vertices() const
unsigned int boundary_id
Definition: types.h:111
static::ExceptionBase & ExcPointNotFound(Point< N > arg1)
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)
unsigned int count_cells_with_subdomain_association(const Triangulation< dim, spacedim > &triangulation, const types::subdomain_id subdomain)
Definition: grid_tools.cc:2807
unsigned int find_closest_vertex(const std::map< unsigned int, Point< spacedim >> &vertices, const Point< spacedim > &p)
Definition: grid_tools.cc:4759
std::vector< typename MeshType::active_cell_iterator > get_patch_around_cell(const typename MeshType::active_cell_iterator &cell)
static::ExceptionBase & ExcInternalError()
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:533