]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Hand back used quadratures during conversion of intersections to points 16391/head
authorJohannes Heinz <johannes.heinz@tuwien.ac.at>
Fri, 29 Dec 2023 09:47:44 +0000 (10:47 +0100)
committerJohannes Heinz <johannes.heinz@tuwien.ac.at>
Fri, 29 Dec 2023 20:50:57 +0000 (21:50 +0100)
doc/news/changes/minor/20231229Heinz [new file with mode: 0644]
include/deal.II/grid/grid_tools.h
source/grid/grid_tools.cc
source/grid/grid_tools.inst.in
tests/remote_point_evaluation/convert_intersections_to_rpe_data_01.cc
tests/remote_point_evaluation/convert_intersections_to_rpe_data_02.cc
tests/remote_point_evaluation/convert_intersections_to_rpe_data_03.cc
tests/remote_point_evaluation/convert_intersections_to_rpe_data_04.cc
tests/remote_point_evaluation/setup_rpe_with_intersections_02.cc

diff --git a/doc/news/changes/minor/20231229Heinz b/doc/news/changes/minor/20231229Heinz
new file mode 100644 (file)
index 0000000..5ad8dfe
--- /dev/null
@@ -0,0 +1,6 @@
+Improved: convert_to_distributed_compute_point_locations_internal() takes a pointer to a vector
+of quadratures as an additional argument. In case this pointer is not a null pointer, the vector is
+populated by the mapped quadratures that are internally used during the conversion.
+<br>
+(Johannes Heinz, 2023/12/29)
+
index 80e6c4104e9d02fce9c176c4ac42f129b68bfe44..81d806d99481b5c19504c48025e6cdb26a90908c 100644 (file)
@@ -877,6 +877,11 @@ namespace GridTools
        * from class members. This can be done without searching for points again
        * since all information is locally known.
        *
+       * @p mapped_quadratures_recv_comp is a pointer to an empty vector of
+       * mapped quadratures. By default it is a `nullptr` and the parameter is
+       * ignored. Otherwise, the vector is filled with the mapped quadrature
+       * rules (in real coordinates) corresponding to recv_components.
+       *
        * The parameter @p consistent_numbering_of_sender_and_receiver can be used to ensure
        * points on sender and receiver side are numbered consistently.
        * This parameter is optional if DistributedComputePointLocationsInternal
@@ -892,6 +897,8 @@ namespace GridTools
         const unsigned int                  n_points_1D,
         const Triangulation<dim, spacedim> &tria,
         const Mapping<dim, spacedim>       &mapping,
+        std::vector<Quadrature<spacedim>>  *mapped_quadratures_recv_comp =
+          nullptr,
         const bool consistent_numbering_of_sender_and_receiver = false) const;
 
     private:
index 43b778a793180454f6d028ef2b4ea5950a5c3653..40da559d17a8da6ab3cb969204807d78357a6729 100644 (file)
@@ -4105,11 +4105,18 @@ namespace GridTools
         const unsigned int                  n_points_1D,
         const Triangulation<dim, spacedim> &tria,
         const Mapping<dim, spacedim>       &mapping,
+        std::vector<Quadrature<spacedim>>  *mapped_quadratures_recv_comp,
         const bool consistent_numbering_of_sender_and_receiver) const
     {
       using CellIterator =
         typename Triangulation<dim, spacedim>::active_cell_iterator;
 
+      if (mapped_quadratures_recv_comp != nullptr)
+        {
+          AssertDimension(mapped_quadratures_recv_comp->size(), 0);
+          mapped_quadratures_recv_comp->reserve(recv_components.size());
+        }
+
       GridTools::internal::DistributedComputePointLocationsInternal<dim,
                                                                     spacedim>
         result;
@@ -4141,6 +4148,10 @@ namespace GridTools
                 result.recv_components.size(), // number of point
                 numbers::invalid_unsigned_int);
             }
+
+          // append quadrature
+          if (mapped_quadratures_recv_comp != nullptr)
+            mapped_quadratures_recv_comp->push_back(quad);
         }
 
       // since empty quadratures might be present we have to set the number
index ea2f87577ae315a6c7be62d5e7235d455463738f..23c56beef66257c98916c26315f36ec5d7a5c586 100644 (file)
@@ -193,6 +193,7 @@ for (deal_II_struct_dimension : DIMENSIONS; deal_II_dimension : DIMENSIONS;
           const unsigned int,
           const Triangulation<deal_II_dimension, deal_II_space_dimension> &,
           const Mapping<deal_II_dimension, deal_II_space_dimension> &,
+          std::vector<Quadrature<deal_II_space_dimension>> *,
           const bool) const;
 
       template internal::DistributedComputeIntersectionLocationsInternal<
index c1e31e45e635553037c9637836154af90969583e..13cb52cf8df779adabd694c7d7f3fdfc6193dc09 100644 (file)
@@ -99,7 +99,7 @@ do_test(const unsigned int n_quad_points)
 
   auto rpe_data = intersection_location
                     .convert_to_distributed_compute_point_locations_internal(
-                      n_quad_points, tria, mapping, true);
+                      n_quad_points, tria, mapping, nullptr, true);
 
   deallog << "Recv Components " << std::endl;
   for (const auto &rc : rpe_data.recv_components)
index 88993c09712d8a787b316771dc2ce67f77c3aefa..518beb66caba404e45e3c317d5fcc13e69203dbe 100644 (file)
@@ -99,7 +99,7 @@ do_test(const unsigned int n_quad_points)
 
   auto rpe_data = intersection_location
                     .convert_to_distributed_compute_point_locations_internal(
-                      n_quad_points, tria, mapping, true);
+                      n_quad_points, tria, mapping, nullptr, true);
 
   deallog << "Recv Components " << std::endl;
   for (const auto &rc : rpe_data.recv_components)
index d55e666cd3d3fe3324346ab8d917d3af1e863676..c070dc142b405b2ab2a3232472a0dd0d3e267772 100644 (file)
@@ -96,7 +96,7 @@ do_test(const unsigned int n_quad_points)
 
   auto rpe_data = intersection_location
                     .convert_to_distributed_compute_point_locations_internal(
-                      n_quad_points, tria, mapping, true);
+                      n_quad_points, tria, mapping, nullptr, true);
 
   deallog << "Recv Components " << std::endl;
   for (const auto &rc : rpe_data.recv_components)
index 6bd08f880d3da47306ed257892ca4cccb7b458d7..f210ebf5cd73deb18c0262825cf972f8ad12354b 100644 (file)
@@ -96,7 +96,7 @@ do_test(const unsigned int n_quad_points)
 
   auto rpe_data = intersection_location
                     .convert_to_distributed_compute_point_locations_internal(
-                      n_quad_points, tria, mapping, true);
+                      n_quad_points, tria, mapping, nullptr, true);
 
   deallog << "Recv Components " << std::endl;
   for (const auto &rc : rpe_data.recv_components)
index 2abfb15fb84e9d8e65073fe88eb19704b870f78c..02caa53cd5219bbf879e83a7333c71d3e987ad6a 100644 (file)
@@ -91,7 +91,7 @@ do_test(const unsigned int n_quad_points)
           intersection_requests.emplace_back(vertices);
         }
     }
-  std::cout << "CALLING FUNCTION" << std::endl;
+
   auto intersection_location =
     GridTools::internal::distributed_compute_intersection_locations<structdim>(
       cache,
@@ -100,15 +100,11 @@ do_test(const unsigned int n_quad_points)
       std::vector<bool>(),
       1.0e-9);
 
-  std::cout << "FINISHED" << std::endl;
-
   auto rpe_intersection_data =
     intersection_location
       .convert_to_distributed_compute_point_locations_internal(n_quad_points,
                                                                tria,
                                                                mapping);
-  std::cout << "WTF" << std::endl;
-
   // setup rpe without additional search
   dealii::Utilities::MPI::RemotePointEvaluation<dim> rpe_intersections;
   rpe_intersections.reinit(rpe_intersection_data, tria, mapping);

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.