almost_infinite_length);
// also note if a vertex is at the boundary
- std::vector<bool> at_boundary (triangulation.n_vertices(), false);
+ std::vector<bool> at_boundary (keep_boundary ? triangulation.n_vertices() : 0, false);
+ // for parallel::shared::Triangulation we need to work on all vertices,
+ // not just the ones related to loacally owned cells;
+ const bool is_parallel_shared
+ = (dynamic_cast<parallel::shared::Triangulation<dim,spacedim>*> (&triangulation) != nullptr);
for (typename Triangulation<dim,spacedim>::active_cell_iterator
cell=triangulation.begin_active(); cell!=triangulation.end(); ++cell)
- if (cell->is_locally_owned())
+ if (is_parallel_shared || cell->is_locally_owned())
{
if (dim>1)
{
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// check GridTools::distort_random for parallel::shared::Triangulation
+
+
+#include "../tests.h"
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/distributed/shared_tria.h>
+#include <deal.II/numerics/data_out.h>
+
+
+
+
+template <int dim>
+void test1 (const bool keep_boundary)
+{
+ const unsigned int my_id = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
+ parallel::shared::Triangulation<dim> tria (MPI_COMM_WORLD);
+ GridGenerator::hyper_cube(tria);
+ tria.refine_global(2);
+ GridTools::distort_random (0.1, tria, keep_boundary);
+
+ deallog << "dim=" << dim << ", keep_boundary=" << keep_boundary << std::endl;
+ std::string filename;
+ if (keep_boundary)
+ filename = "keep_true-";
+ else
+ filename = "keep_false-";
+ filename += Utilities::int_to_string(dim);
+
+ std::stringstream ss;
+ GridOut().write_gnuplot (tria, ss);
+ const std::string local_grid = ss.str();
+ deallog << checksum(local_grid.begin(), local_grid.end()) << std::endl;
+
+ deallog << "OK" << std::endl;
+}
+
+
+
+int main (int argc, char *argv[])
+{
+ Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, testing_max_num_threads());
+
+ MPILogInitAll log_all;
+
+ test1<2> (true);
+ test1<2> (false);
+ test1<3> (true);
+ test1<3> (false);
+
+ return 0;
+}
--- /dev/null
+
+DEAL:0::dim=2, keep_boundary=1
+DEAL:0::561638746
+DEAL:0::OK
+DEAL:0::dim=2, keep_boundary=0
+DEAL:0::1946764437
+DEAL:0::OK
+DEAL:0::dim=3, keep_boundary=1
+DEAL:0::3893493253
+DEAL:0::OK
+DEAL:0::dim=3, keep_boundary=0
+DEAL:0::2979128258
+DEAL:0::OK
+
+DEAL:1::dim=2, keep_boundary=1
+DEAL:1::561638746
+DEAL:1::OK
+DEAL:1::dim=2, keep_boundary=0
+DEAL:1::1946764437
+DEAL:1::OK
+DEAL:1::dim=3, keep_boundary=1
+DEAL:1::3893493253
+DEAL:1::OK
+DEAL:1::dim=3, keep_boundary=0
+DEAL:1::2979128258
+DEAL:1::OK
+
+
+DEAL:2::dim=2, keep_boundary=1
+DEAL:2::561638746
+DEAL:2::OK
+DEAL:2::dim=2, keep_boundary=0
+DEAL:2::1946764437
+DEAL:2::OK
+DEAL:2::dim=3, keep_boundary=1
+DEAL:2::3893493253
+DEAL:2::OK
+DEAL:2::dim=3, keep_boundary=0
+DEAL:2::2979128258
+DEAL:2::OK
+