]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Changed GridTools::distributed_compute_ptloc to use global description of the mesh... 6251/head
authorGiovanni Alzetta <giovannialzetta@hotmail.it>
Mon, 16 Apr 2018 11:43:34 +0000 (13:43 +0200)
committerGiovanni Alzetta <giovannialzetta@hotmail.it>
Mon, 16 Apr 2018 20:40:29 +0000 (22:40 +0200)
doc/news/changes/minor/20180416GiovanniAlzetta [new file with mode: 0644]
include/deal.II/grid/grid_tools.h
source/grid/grid_tools.cc
source/grid/grid_tools.inst.in
tests/grid/distributed_compute_point_locations_01.cc
tests/grid/distributed_compute_point_locations_02.cc

diff --git a/doc/news/changes/minor/20180416GiovanniAlzetta b/doc/news/changes/minor/20180416GiovanniAlzetta
new file mode 100644 (file)
index 0000000..bf2f880
--- /dev/null
@@ -0,0 +1,5 @@
+Changed: GridTools::distributed_compute_point_locations now takes as input the
+global description of the manifold using bounding boxes.
+<br>
+(Giovanni Alzetta, 2018/04/16)
+
index 3c9a3081baf7f49c410133e48a39269a9da3b0f8..2ee6f31ed03e9036c252ff88b994dcc3b96883e3 100644 (file)
@@ -683,9 +683,14 @@ namespace GridTools
    * @param[in] local_points the array of points owned by the current process. Every
    *  process can have a different array of points which can be empty and not
    *  contained within the locally owned part of the triangulation
-   * @param[in] local_bbox the description of the locally owned part of the mesh made
-   *  with bounding boxes. It can be obtained from
-   *  GridTools::compute_mesh_predicate_bounding_box
+   * @param[in] global_bboxes a vector of vectors of bounding boxes; it describes
+   *  the locally owned part of the mesh for each process.
+   *  The bounding boxes describing which part of the mesh is locally owned by
+   *  process with rank rk are contained in global_bboxes[rk].
+   *  The local description can be obtained from
+   *  GridTools::compute_mesh_predicate_bounding_box ; then the global one can
+   *  be obtained using either GridTools::exchange_local_bounding_boxes
+   *  or Utilities::MPI::all_gather
    * @return A tuple containing the quadrature information
    *
    * The elements of the output tuple are:
@@ -739,9 +744,9 @@ namespace GridTools
   return_type
 #endif
       distributed_compute_point_locations
-      (const GridTools::Cache<dim,spacedim>                &cache,
-       const std::vector<Point<spacedim> >                 &local_points,
-       const std::vector< BoundingBox<spacedim> >          &local_bbox);
+      (const GridTools::Cache<dim,spacedim>                        &cache,
+       const std::vector<Point<spacedim> >                         &local_points,
+       const std::vector< std::vector< BoundingBox<spacedim> > >   &global_bboxes);
 
   /**
    * Return a map of index:Point<spacedim>, containing the used vertices of the
index be39c27e8ba784cd8cfe3ce9183e45dc3332eb47..e63e6de417f2dfe4890d6c22b5cb3d4bbcd7a776 100644 (file)
@@ -4213,14 +4213,14 @@ next_cell:
       std::vector< std::vector< unsigned int > >
       >
       distributed_compute_point_locations
-      (const GridTools::Cache<dim,spacedim>                &cache,
-       const std::vector<Point<spacedim> >                 &local_points,
-       const std::vector< BoundingBox<spacedim> >          &local_bbox)
+      (const GridTools::Cache<dim,spacedim>                         &cache,
+       const std::vector<Point<spacedim> >                          &local_points,
+       const std::vector< std::vector< BoundingBox<spacedim> > >    &global_bboxes)
   {
 #ifndef DEAL_II_WITH_MPI
     (void)cache;
     (void)local_points;
-    (void)local_bbox;
+    (void)global_bboxes;
     Assert(false, ExcMessage("GridTools::distributed_compute_point_locations() requires MPI."));
     std::tuple<
     std::vector< typename Triangulation<dim, spacedim>::active_cell_iterator >,
@@ -4257,10 +4257,6 @@ next_cell:
         internal::distributed_cptloc::cell_hash<dim,spacedim> >
         temporary_unmap;
 
-    // Obtaining the global mesh description through an all to all communication
-    std::vector< std::vector< BoundingBox<spacedim> > > global_bounding_boxes;
-    global_bounding_boxes = Utilities::MPI::all_gather(mpi_communicator,local_bbox);
-
     // Step 1 (part 1): Using the bounding boxes to guess the owner of each points
     // in local_points
     unsigned int my_rank = Utilities::MPI::this_mpi_process(mpi_communicator);
@@ -4269,7 +4265,7 @@ next_cell:
     std::tuple< std::vector< std::vector< unsigned int > >, std::map< unsigned int, unsigned int >,
         std::map< unsigned int, std::vector< unsigned int > > > guessed_points;
     guessed_points =
-      GridTools::guess_point_owner(global_bounding_boxes, local_points);
+      GridTools::guess_point_owner(global_bboxes, local_points);
 
     // Preparing to call compute point locations on points which are/might be
     // local
index 5268a15fd7ad888007e463f1c347eff4dc80865f..3170d290c477afddba64a07ba119ea323cab5473 100644 (file)
@@ -88,9 +88,9 @@ for (deal_II_dimension : DIMENSIONS ; deal_II_space_dimension : SPACE_DIMENSIONS
         std::vector< std::vector< unsigned int > >
         >
     distributed_compute_point_locations
-    (const Cache< deal_II_dimension, deal_II_space_dimension >   &,
-     const std::vector< Point< deal_II_space_dimension > >       &,
-     const std::vector< BoundingBox< deal_II_space_dimension > > &);
+    (const Cache< deal_II_dimension, deal_II_space_dimension >                  &,
+     const std::vector< Point< deal_II_space_dimension > >                      &,
+     const std::vector< std::vector< BoundingBox< deal_II_space_dimension > > > &);
     \}
 
 #endif
index cf65d62a46f1753b103aaad9f1942a8fc37c4b51..ecf0cb195cd2f657cd18be3c321406d7b9bdc2d5 100644 (file)
@@ -55,11 +55,14 @@ void test_compute_pt_loc(unsigned int n_points)
   std::vector< BoundingBox<dim> > local_bbox = GridTools::compute_mesh_predicate_bounding_box
                                                (cache.get_triangulation(),
                                                 std::function<bool (const typename Triangulation<dim>::active_cell_iterator &)>(locally_owned_cell_predicate),
-                                                1, true, 4); // These options should be passed
-  // Using the distributed version of compute point location
+                                                1, true, 4);
+
+  // Obtaining the global mesh description through an all to all communication
+  std::vector< std::vector< BoundingBox<dim> > > global_bboxes;
+  global_bboxes = Utilities::MPI::all_gather(mpi_communicator,local_bbox);
 
-  // Using the distributed version
-  auto output_tuple = distributed_compute_point_locations(cache,points,local_bbox);
+  // Using the distributed version of compute point location
+  auto output_tuple = distributed_compute_point_locations(cache,points,global_bboxes);
   // Testing in serial against the serial version
   auto cell_qpoint_map = GridTools::compute_point_locations(cache,points);
 
index 08208c6eb3402dd2ea676e3e6a0247c6ae32ca72..1ead16c1f06f5fac803471b38a73e7089faeb95b 100644 (file)
@@ -121,9 +121,13 @@ void test_compute_pt_loc(unsigned int ref_cube, unsigned int ref_sphere)
                                                (cache.get_triangulation(), locally_owned_cell_predicate,
                                                 1, true, 4);
 
+  // Obtaining the global mesh description through an all to all communication
+  std::vector< std::vector< BoundingBox<dim> > > global_bboxes;
+  global_bboxes = Utilities::MPI::all_gather(mpi_communicator,local_bbox);
+
   // Using the distributed version of compute point location
   auto output_tuple = distributed_compute_point_locations
-                      (cache,loc_owned_points,local_bbox);
+                      (cache,loc_owned_points,global_bboxes);
   deallog << "Comparing results" << std::endl;
   const auto &output_cells = std::get<0>(output_tuple);
   const auto &output_qpoints = std::get<1>(output_tuple);

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.