]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Introduce GridTools::get_locally_owned_vertices().
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 9 Nov 2014 21:58:08 +0000 (15:58 -0600)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 12 Nov 2014 13:04:11 +0000 (07:04 -0600)
Also use it in an existing test and in an assertion.

doc/news/changes.h
include/deal.II/distributed/tria.h
include/deal.II/grid/grid_tools.h
source/distributed/tria.cc
source/grid/grid_tools.cc
source/grid/grid_tools.inst.in
tests/mpi/communicate_moved_vertices_01.cc

index 507a9ea3ba24a5c4d6b24977e6b1387d03e8713f..e52bbe33e063282ba721a8ad4ea9a85f8f44fc0a 100644 (file)
@@ -339,6 +339,22 @@ inconvenience this causes.
   (Fahad Alrashed, 2014/11/09)
   </li>
 
+  <li> New: GridTools::get_locally_owned_vertices() allows to query
+  which vertices of a triangulation are owned by the current
+  processor.
+  <br>
+  (Wolfgang Bangerth, 2014/11/09)
+  </li>
+
+  <li> New: parallel::distributed::Triangulation::communicate_locally_moved_vertices()
+  allows to
+  update vertex positions that have been moved just locally on distributed
+  meshes. GridTools::distort_random now works for distributed meshes and 
+  hanging nodes in 3D as well.
+  <br>
+  (Fahad Alrashed, 2014/11/09)
+  </li>
+
   <li> New: TableHandler objects can be cleared - i.e. reset to a 
   zero-sized state.
   <br>
index 4fb06a477e0dc1f75cbdd8d8236f8b1a62078ee7..fc13758a915c48bcb9c523f74e221519cdc2cccc 100644 (file)
@@ -512,6 +512,11 @@ namespace parallel
        * every vertex, call this function, and before refining or coarsening
        * the mesh apply the opposite offset and call this function again.
        *
+       * @param vertex_locally_moved A bitmap indicating which vertices have
+       *   been moved. The size of this array must be equal to
+       *   Triangulation::n_vertices() and must be a subset of those vertices
+       *   flagged by GridTools::get_locally_owned_vertices().
+       *
        * @see This function is used, for example, in GridTools::distort_random().
        */
       void
index eb00e598fcacda7358e9cd3b5975076471fd5b6b..460c9886180de5640e8b91702ee1d3de9dcdc8fc 100644 (file)
@@ -767,6 +767,38 @@ namespace GridTools
   count_cells_with_subdomain_association (const Triangulation<dim, spacedim> &triangulation,
                                           const types::subdomain_id         subdomain);
 
+
+  /**
+   * For a triangulation, return a mask that represents which of its vertices
+   * are "owned" by the current process in the same way as we talk about
+   * locally owned cells or degrees of freedom (see @ref GlossLocallyOwnedCell
+   * and @ref GlossLocallyOwnedDof). For the purpose of this function,
+   * we define a locally owned vertex as follows: a vertex is owned by
+   * that processor with the smallest subdomain id (which equals the MPI
+   * rank of that processor) among all owners of cells adjacent to this vertex.
+   * In other words, vertices that are in the interior of a partition
+   * of the triangulation are owned by the owner of this partition; for
+   * vertices that lie on the boundary between two or more partitions,
+   * the owner is the processor with the least subdomain_id among all
+   * adjacent subdomains.
+   *
+   * For sequential triangulations (as opposed to, for example,
+   * parallel::distributed::Triangulation), every user vertex is of course
+   * owned by the current processor, i.e., the function returns
+   * Triangulation::get_used_vertices(). For parallel triangulations, the
+   * returned mask is a subset of what Triangulation::get_used_vertices()
+   * returns.
+   *
+   * @param triangulation The triangulation of which the function evaluates
+   *   which vertices are locally owned.
+   * @return The subset of vertices, as described above. The length of the
+   *   returned array equals Triangulation.n_vertices() and may, consequently,
+   *   be larger than Triangulation::n_used_vertices().
+   */
+  template <int dim, int spacedim>
+  std::vector<bool>
+  get_locally_owned_vertices (const Triangulation<dim,spacedim> &triangulation);
+
   /*@}*/
   /**
    *  @name Comparing different meshes
index 643e3e3a64a6aa0753b1d1d5e1788b95b4e61c72..8e0e70df137a72217ffbb3b0bf1d55d05f3d8992 100644 (file)
@@ -3496,6 +3496,19 @@ namespace parallel
       Assert (vertex_locally_moved.size() == this->n_vertices(),
               ExcDimensionMismatch(vertex_locally_moved.size(),
                                    this->n_vertices()));
+#ifdef DEBUG
+      {
+        const std::vector<bool> locally_owned_vertices
+          = GridTools::get_locally_owned_vertices (*this);
+        for (unsigned int i=0; i<locally_owned_vertices.size(); ++i)
+          Assert ((vertex_locally_moved[i] == false)
+                  ||
+                  (locally_owned_vertices[i] == true),
+                  ExcMessage ("The vertex_locally_moved argument must not "
+                              "contain vertices that are not locally owned"));
+      }
+#endif
+
 
       std::map<unsigned int, std::set<dealii::types::subdomain_id> >
       vertices_with_ghost_neighbors;
index b736aadf55cae0edaf487ef7900f934dc718fe6f..2bf29901e13242f2a5156143ccd8d919159ba596 100644 (file)
@@ -37,6 +37,8 @@
 
 #include <cmath>
 #include <numeric>
+#include <list>
+#include <set>
 
 
 DEAL_II_NAMESPACE_OPEN
@@ -1468,7 +1470,6 @@ next_cell:
 
   template <int dim, int spacedim>
   unsigned int
-
   count_cells_with_subdomain_association (const Triangulation<dim, spacedim> &triangulation,
                                           const types::subdomain_id       subdomain)
   {
@@ -1484,6 +1485,35 @@ next_cell:
 
 
 
+  template <int dim, int spacedim>
+  std::vector<bool>
+  get_locally_owned_vertices (const Triangulation<dim,spacedim> &triangulation)
+  {
+    // start with all vertices
+    std::vector<bool> locally_owned_vertices = triangulation.get_used_vertices();
+
+    // if the triangulation is distributed, eliminate those that
+    // are owned by other processors -- either because the vertex is
+    // on an artificial cell, or because it is on a ghost cell with
+    // a smaller subdomain
+    if (const parallel::distributed::Triangulation<dim,spacedim> *tr
+        = dynamic_cast<const parallel::distributed::Triangulation<dim,spacedim> *>
+        (&triangulation))
+      for (typename Triangulation<dim,spacedim>::active_cell_iterator
+           cell = triangulation.begin_active();
+           cell != triangulation.end(); ++cell)
+        if (cell->is_artificial()
+            ||
+            (cell->is_ghost() &&
+             (cell->subdomain_id() < tr->locally_owned_subdomain())))
+          for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
+            locally_owned_vertices[cell->vertex_index(v)] = false;
+
+    return locally_owned_vertices;
+  }
+
+
+
   template <typename Container>
   std::list<std::pair<typename Container::cell_iterator,
       typename Container::cell_iterator> >
index a11c8cbf729e2b89d29e762068d215fb06cad897..b949d29f2d19b3a86d882bbfc30a4dc4005c35ee 100644 (file)
@@ -219,6 +219,10 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS
       const Triangulation<deal_II_dimension, deal_II_space_dimension> &,
       const types::subdomain_id);
 
+    template
+    std::vector<bool>
+    get_locally_owned_vertices (const Triangulation<deal_II_dimension, deal_II_space_dimension> &);
+
     template
       double
       minimal_cell_diameter (const Triangulation<deal_II_dimension, deal_II_space_dimension> &triangulation);
index e500951aeaa6a21a55b28e4bc48f4381d2ffb9b6..0c7081fa6a3e74fec58e78eb21bbcf7bd2666ebe 100644 (file)
 #include <fstream>
 
 
-template <int dim, int spacedim>
-std::vector<bool>
-get_locally_owned_vertices (const Triangulation<dim,spacedim> &triangulation)
-{
-  std::vector<bool> locally_owned_vertices = triangulation.get_used_vertices();
-
-  if (const parallel::distributed::Triangulation<dim,spacedim> *tr
-      = dynamic_cast<const parallel::distributed::Triangulation<dim,spacedim> *>
-      (&triangulation))
-    for (typename Triangulation<dim,spacedim>::active_cell_iterator
-        cell = triangulation.begin_active();
-       cell != triangulation.end(); ++cell)
-    if (cell->is_artificial()
-       ||
-       (cell->is_ghost() && (cell->subdomain_id() < tr->locally_owned_subdomain())))
-      for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
-       locally_owned_vertices[cell->vertex_index(v)] = false;
-
-  return locally_owned_vertices;
-}
-
-
 template<int dim>
 void test()
 {
@@ -67,7 +45,7 @@ void test()
   tr.refine_global(2);
 
   const std::vector<bool> locally_owned_vertices
-    = get_locally_owned_vertices (tr);
+    = GridTools::get_locally_owned_vertices (tr);
   
   if (myid == 0)
     deallog << "#vertices = "

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.