&data_range);
};
- /* ---------------------- inline and template functions ------------------ */
+ /**
+ * This class manages the DataOut of a Particle Handler
+ * It currently only supports witing the particle position
+ * and their ID
+ *
+ * @ingroup Particle
+ *
+ * @author : Bruno Blais, Luca Heltai 2019
+ */
+
+ template <int dim, int spacedim>
+ class ParticleOutput : public dealii::DataOutInterface<0, spacedim>
+ {
+ public:
+ ParticleOutput() = default;
+
+ void
+ build_patches(const Particles::ParticleHandler<dim, spacedim> &particles)
+ {
+ dataset_names.reserve(1);
+ dataset_names.emplace_back("id");
+ patches.resize(particles.n_locally_owned_particles());
+
+ auto particle = particles.begin();
+ for (unsigned int i = 0; particle != particles.end(); ++particle, ++i)
+ {
+ patches[i].vertices[0] = particle->get_location();
+ patches[i].patch_index = i;
+ patches[i].n_subdivisions = 1;
+ patches[i].data.reinit(dataset_names.size(), 1);
+
+ patches[i].data(0, 0) = particle->get_id();
+ }
+ }
+
+
+
+ ~ParticleOutput() = default;
+
+ protected:
+ /**
+ * Implementation of the corresponding function of the base class.
+ */
+ virtual const std::vector<DataOutBase::Patch<0, spacedim>> &
+ get_patches() const
+ {
+ return patches;
+ }
+
+ /**
+ * Implementation of the corresponding function of the base class.
+ */
+ virtual std::vector<std::string>
+ get_dataset_names() const
+ {
+ return dataset_names;
+ }
+
+ private:
+ std::vector<DataOutBase::Patch<0, spacedim>> patches;
+
+ /**
+ * A list of field names for all data components stored in patches.
+ */
+ std::vector<std::string> dataset_names;
+ };
+
+
+ /* ---------------------- inline and template functions ------------------
+ */
+
+
template <int dim, int spacedim>
template <class Archive>
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 - 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.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Create a particle handler then generate an output
+// Tests the Visualization class of the particle handler
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/fe/mapping_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_tools.h>
+
+#include <deal.II/particles/particle_handler.h>
+
+#include "../tests.h"
+
+template <int dim, int spacedim>
+void
+test()
+{
+ {
+ parallel::distributed::Triangulation<dim, spacedim> tr(MPI_COMM_WORLD);
+
+ GridGenerator::hyper_cube(tr);
+ tr.refine_global(2);
+ MappingQ<dim, spacedim> mapping(1);
+
+ // both processes create a particle handler, but only the first creates
+ // particles
+ Particles::ParticleHandler<dim, spacedim> particle_handler(tr, mapping);
+ std::vector<Point<spacedim>> position(2);
+ std::vector<Point<dim>> reference_position(2);
+
+ for (unsigned int i = 0; i < dim; ++i)
+ {
+ position[0](i) = 0.125;
+ position[1](i) = 0.525;
+ }
+
+ Particles::Particle<dim, spacedim> particle1(position[0],
+ reference_position[0],
+ 0);
+ Particles::Particle<dim, spacedim> particle2(position[1],
+ reference_position[1],
+ 1);
+
+ typename Triangulation<dim, spacedim>::active_cell_iterator cell1(&tr,
+ 2,
+ 0);
+ typename Triangulation<dim, spacedim>::active_cell_iterator cell2(&tr,
+ 2,
+ 0);
+
+ particle_handler.insert_particle(particle1, cell1);
+ particle_handler.insert_particle(particle2, cell2);
+
+ Particles::ParticleOutput<dim, spacedim> particle_viz;
+ particle_viz.build_patches(particle_handler);
+ particle_viz.write_gnuplot(deallog.get_file_stream());
+
+ for (const auto &particle : particle_handler)
+ deallog << "Particle id " << particle.get_id() << " is in cell "
+ << particle.get_surrounding_cell(tr) << std::endl;
+ }
+
+ deallog << "OK" << std::endl;
+}
+
+
+
+int
+main(int argc, char *argv[])
+{
+ Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+
+ MPILogInitAll all;
+
+ deallog.push("2d/2d");
+ test<2, 2>();
+ deallog.pop();
+ deallog.push("2d/3d");
+ test<2, 3>();
+ deallog.pop();
+ deallog.push("3d/3d");
+ test<3, 3>();
+ deallog.pop();
+}
--- /dev/null
+
+# This file was generated by the deal.II library.
+
+
+#
+# For a description of the GNUPLOT format see the GNUPLOT manual.
+#
+# <x> <y> <id>
+0.125000 0.125000 0.00000
+
+0.525000 0.525000 1.00000
+
+DEAL:0:2d/2d::Particle id 0 is in cell 2.0
+DEAL:0:2d/2d::Particle id 1 is in cell 2.0
+DEAL:0:2d/2d::OK
+# This file was generated by the deal.II library.
+
+
+#
+# For a description of the GNUPLOT format see the GNUPLOT manual.
+#
+# <x> <y> <z> <id>
+0.125000 0.125000 0.00000 0.00000
+
+0.525000 0.525000 0.00000 1.00000
+
+DEAL:0:2d/3d::Particle id 0 is in cell 2.0
+DEAL:0:2d/3d::Particle id 1 is in cell 2.0
+DEAL:0:2d/3d::OK
+# This file was generated by the deal.II library.
+
+
+#
+# For a description of the GNUPLOT format see the GNUPLOT manual.
+#
+# <x> <y> <z> <id>
+0.125000 0.125000 0.125000 0.00000
+
+0.525000 0.525000 0.525000 1.00000
+
+DEAL:0:3d/3d::Particle id 0 is in cell 2.0
+DEAL:0:3d/3d::Particle id 1 is in cell 2.0
+DEAL:0:3d/3d::OK