--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2003 - 2017 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.
+//
+// ---------------------------------------------------------------------
+
+
+// Test that GridTools::extract_used_vertices works with MappingQEulerian
+// in parallel. This is a variant of _07
+
+#include "../tests.h"
+
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/lac/full_matrix.h>
+#include <deal.II/lac/identity_matrix.h>
+#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/base/thread_management.h>
+#include <deal.II/base/function.h>
+#include <deal.II/lac/vector.h>
+#include <deal.II/lac/vector_memory.h>
+#include <deal.II/lac/filtered_matrix.h>
+#include <deal.II/lac/precondition.h>
+#include <deal.II/lac/solver_cg.h>
+#include <deal.II/lac/sparse_matrix.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_reordering.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/distributed/tria.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/dofs/dof_tools.h>
+#include <deal.II/lac/constraint_matrix.h>
+#include <deal.II/lac/trilinos_vector.h>
+#include <deal.II/fe/mapping_q1.h>
+#include <deal.II/fe/mapping_q_eulerian.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/numerics/vector_tools.h>
+#include <deal.II/numerics/data_out.h>
+#include <deal.II/base/multithread_info.h>
+#include <iostream>
+#include <vector>
+
+
+using namespace dealii;
+
+
+template <int dim>
+class Displacement : public Function<dim>
+{
+public:
+ Displacement() :
+ Function<dim>(dim)
+ {}
+
+ double value (const Point<dim> &p,
+ const unsigned int component) const
+ {
+ return p[component]*p[0]/2;
+ }
+
+ void vector_value (const Point<dim> &p,
+ Vector<double> &v) const
+ {
+ for (unsigned int i=0; i<dim; ++i)
+ v(i) = p[i]*p[0]/2;
+ }
+};
+
+
+template <int dim>
+void test ()
+{
+ deallog << "dim=" << dim << std::endl;
+ unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+ unsigned int numproc = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD);
+
+ if (myid==0)
+ deallog << "numproc=" << numproc << std::endl;
+
+ parallel::distributed::Triangulation<dim>
+ triangulation (MPI_COMM_WORLD,
+ typename Triangulation<dim>::MeshSmoothing
+ (Triangulation<dim>::smoothing_on_refinement |
+ Triangulation<dim>::smoothing_on_coarsening));
+ GridGenerator::hyper_cube (triangulation, -1, 1);
+ triangulation.refine_global(2);
+
+ FESystem<dim> fe(FE_Q<dim>(2), dim);
+ DoFHandler<dim> dof_handler(triangulation);
+ dof_handler.distribute_dofs(fe);
+
+ IndexSet locally_owned_dofs = dof_handler.locally_owned_dofs ();
+ IndexSet locally_relevant_dofs;
+ DoFTools::extract_locally_relevant_dofs (dof_handler,locally_relevant_dofs);
+
+ TrilinosWrappers::MPI::Vector x(locally_owned_dofs, MPI_COMM_WORLD);
+ TrilinosWrappers::MPI::Vector x_relevant(locally_owned_dofs,locally_relevant_dofs, MPI_COMM_WORLD);
+
+ VectorTools::interpolate (dof_handler,
+ Displacement<dim>(),
+ x);
+ x_relevant = x;
+
+ MappingQEulerian<dim, TrilinosWrappers::MPI::Vector> euler(2, dof_handler, x_relevant);
+
+ // now the actual test
+ std::map<unsigned int, Point<dim> > active_vertices = GridTools::extract_used_vertices(triangulation, euler);
+
+ for (auto point = active_vertices.begin(); point != active_vertices.end(); ++point)
+ deallog << point->second << " ";
+
+ deallog << std::endl << "Ok." << std::endl;
+}
+
+
+
+int main (int argc, char *argv[])
+{
+ Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, 1);
+ MPILogInitAll log;
+
+ test<2> ();
+}
+
+