From de6084af18472dba42159c024d90aafb19562d63 Mon Sep 17 00:00:00 2001 From: Wolfgang Bangerth Date: Sun, 9 Nov 2014 15:58:08 -0600 Subject: [PATCH] Introduce GridTools::get_locally_owned_vertices(). Also use it in an existing test and in an assertion. --- doc/news/changes.h | 16 +++++++++++ include/deal.II/distributed/tria.h | 5 ++++ include/deal.II/grid/grid_tools.h | 32 ++++++++++++++++++++++ source/distributed/tria.cc | 13 +++++++++ source/grid/grid_tools.cc | 32 +++++++++++++++++++++- source/grid/grid_tools.inst.in | 4 +++ tests/mpi/communicate_moved_vertices_01.cc | 24 +--------------- 7 files changed, 102 insertions(+), 24 deletions(-) diff --git a/doc/news/changes.h b/doc/news/changes.h index 507a9ea3ba..e52bbe33e0 100644 --- a/doc/news/changes.h +++ b/doc/news/changes.h @@ -339,6 +339,22 @@ inconvenience this causes. (Fahad Alrashed, 2014/11/09) +
  • New: GridTools::get_locally_owned_vertices() allows to query + which vertices of a triangulation are owned by the current + processor. +
    + (Wolfgang Bangerth, 2014/11/09) +
  • + +
  • 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. +
    + (Fahad Alrashed, 2014/11/09) +
  • +
  • New: TableHandler objects can be cleared - i.e. reset to a zero-sized state.
    diff --git a/include/deal.II/distributed/tria.h b/include/deal.II/distributed/tria.h index 4fb06a477e..fc13758a91 100644 --- a/include/deal.II/distributed/tria.h +++ b/include/deal.II/distributed/tria.h @@ -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 diff --git a/include/deal.II/grid/grid_tools.h b/include/deal.II/grid/grid_tools.h index eb00e598fc..460c988618 100644 --- a/include/deal.II/grid/grid_tools.h +++ b/include/deal.II/grid/grid_tools.h @@ -767,6 +767,38 @@ namespace GridTools count_cells_with_subdomain_association (const Triangulation &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 + std::vector + get_locally_owned_vertices (const Triangulation &triangulation); + /*@}*/ /** * @name Comparing different meshes diff --git a/source/distributed/tria.cc b/source/distributed/tria.cc index 643e3e3a64..8e0e70df13 100644 --- a/source/distributed/tria.cc +++ b/source/distributed/tria.cc @@ -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 locally_owned_vertices + = GridTools::get_locally_owned_vertices (*this); + for (unsigned int i=0; i > vertices_with_ghost_neighbors; diff --git a/source/grid/grid_tools.cc b/source/grid/grid_tools.cc index b736aadf55..2bf29901e1 100644 --- a/source/grid/grid_tools.cc +++ b/source/grid/grid_tools.cc @@ -37,6 +37,8 @@ #include #include +#include +#include DEAL_II_NAMESPACE_OPEN @@ -1468,7 +1470,6 @@ next_cell: template unsigned int - count_cells_with_subdomain_association (const Triangulation &triangulation, const types::subdomain_id subdomain) { @@ -1484,6 +1485,35 @@ next_cell: + template + std::vector + get_locally_owned_vertices (const Triangulation &triangulation) + { + // start with all vertices + std::vector 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 *tr + = dynamic_cast *> + (&triangulation)) + for (typename Triangulation::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::vertices_per_cell; ++v) + locally_owned_vertices[cell->vertex_index(v)] = false; + + return locally_owned_vertices; + } + + + template std::list > diff --git a/source/grid/grid_tools.inst.in b/source/grid/grid_tools.inst.in index a11c8cbf72..b949d29f2d 100644 --- a/source/grid/grid_tools.inst.in +++ b/source/grid/grid_tools.inst.in @@ -219,6 +219,10 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS const Triangulation &, const types::subdomain_id); + template + std::vector + get_locally_owned_vertices (const Triangulation &); + template double minimal_cell_diameter (const Triangulation &triangulation); diff --git a/tests/mpi/communicate_moved_vertices_01.cc b/tests/mpi/communicate_moved_vertices_01.cc index e500951aea..0c7081fa6a 100644 --- a/tests/mpi/communicate_moved_vertices_01.cc +++ b/tests/mpi/communicate_moved_vertices_01.cc @@ -34,28 +34,6 @@ #include -template -std::vector -get_locally_owned_vertices (const Triangulation &triangulation) -{ - std::vector locally_owned_vertices = triangulation.get_used_vertices(); - - if (const parallel::distributed::Triangulation *tr - = dynamic_cast *> - (&triangulation)) - for (typename Triangulation::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::vertices_per_cell; ++v) - locally_owned_vertices[cell->vertex_index(v)] = false; - - return locally_owned_vertices; -} - - template void test() { @@ -67,7 +45,7 @@ void test() tr.refine_global(2); const std::vector locally_owned_vertices - = get_locally_owned_vertices (tr); + = GridTools::get_locally_owned_vertices (tr); if (myid == 0) deallog << "#vertices = " -- 2.39.5