Also use it in an existing test and in an assertion.
(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>
* 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
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
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;
#include <cmath>
#include <numeric>
+#include <list>
+#include <set>
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim>
unsigned int
-
count_cells_with_subdomain_association (const Triangulation<dim, spacedim> &triangulation,
const types::subdomain_id subdomain)
{
+ 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> >
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);
#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()
{
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 = "