--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 - 2019 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.
+//
+// ---------------------------------------------------------------------
+#ifndef dealii_particles_data_out_h
+#define dealii_particles_data_out_h
+
+#include <deal.II/base/array_view.h>
+#include <deal.II/base/data_out_base.h>
+#include <deal.II/base/mpi.h>
+#include <deal.II/base/smartpointer.h>
+#include <deal.II/base/subscriptor.h>
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/fe/mapping.h>
+
+#include <deal.II/particles/particle.h>
+#include <deal.II/particles/particle_iterator.h>
+#include <deal.II/particles/property_pool.h>
+
+#include <boost/range/iterator_range.hpp>
+#include <boost/serialization/map.hpp>
+
+DEAL_II_NAMESPACE_OPEN
+
+#ifdef DEAL_II_WITH_P4EST
+
+namespace Particles
+{
+ /**
+ * This class manages the DataOut of a Particle Handler
+ * From a particle handler, it generates patches which can then be used to
+ * write traditional output files. This class currently only supports witing
+ * the particle position and their ID and does not allow to write the
+ * properties attached to the particles
+ *
+ * @ingroup Particle
+ *
+ * @author Bruno Blais, Luca Heltai 2019
+ */
+
+ template <int dim, int spacedim>
+ class DataOut : public dealii::DataOutInterface<0, spacedim>
+ {
+ public:
+ DataOut() = default;
+
+ ~DataOut() = default;
+
+
+ /**
+ * Build the particles for a given partcle handler
+ *
+ * @param [in] particle_handler A particle handler for which the patches will be build
+ * A dim=0 patch is build for each particle. The position of the particle is
+ * used to build the node position and the ID of the particle is added as a
+ * single data element
+ *
+ *
+ * @author Bruno Blais, Luca Heltai 2019
+ */
+ void
+ build_patches(const Particles::ParticleHandler<dim, spacedim> &particles);
+
+ protected:
+ /**
+ * Returns the patches built by the data_out class which was previously
+ * built using a particle handler
+ */
+ virtual const std::vector<DataOutBase::Patch<0, spacedim>> &
+ get_patches() const override;
+ // {
+ // return patches;
+ // }
+
+ /**
+ * Returns the name of the data sets associated with the patches. In the
+ * current implementation the particles only contain the ID
+ */
+ virtual std::vector<std::string>
+ get_dataset_names() const override;
+ // {
+ // 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;
+ };
+
+} // namespace Particles
+
+#endif // DEAL_II_WITH_P4EST
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
&data_range);
};
- /**
- * 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 ------------------
particle_handler.cc
generators.cc
property_pool.cc
+ data_out.cc
)
SET(_inst
particle_iterator.inst.in
particle_handler.inst.in
generators.inst.in
+ data_out.inst.in
)
FILE(GLOB _header
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 - 2019 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.
+//
+// ---------------------------------------------------------------------
+
+#include <deal.II/particles/data_out.h>
+#include <deal.II/particles/particle_handler.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+#ifdef DEAL_II_WITH_P4EST
+
+namespace Particles
+{
+ template <int dim, int spacedim>
+ void
+ DataOut<dim, spacedim>::build_patches(
+ const Particles::ParticleHandler<dim, spacedim> &particles)
+ {
+ dataset_names.clear();
+ dataset_names.push_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();
+ }
+ }
+
+ template <int dim, int spacedim>
+ const std::vector<DataOutBase::Patch<0, spacedim>> &
+ DataOut<dim, spacedim>::get_patches() const
+ {
+ return patches;
+ }
+
+ template <int dim, int spacedim>
+ std::vector<std::string>
+ DataOut<dim, spacedim>::get_dataset_names() const
+ {
+ return dataset_names;
+ }
+
+} // namespace Particles
+
+#endif // DEAL_II_WITH_P4EST
+
+DEAL_II_NAMESPACE_CLOSE
+
+DEAL_II_NAMESPACE_OPEN
+
+#ifdef DEAL_II_WITH_P4EST
+# include "data_out.inst"
+#endif
+
+DEAL_II_NAMESPACE_CLOSE
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 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.
+//
+// ---------------------------------------------------------------------
+
+
+for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
+ {
+#if deal_II_dimension <= deal_II_space_dimension
+ namespace Particles
+ \{
+ template class DataOut<deal_II_dimension, deal_II_space_dimension>;
+ \}
+#endif
+ }
// Create a particle handler then generate an output
// Tests the Visualization class of the particle handler
+#include <deal.II/base/data_out_base.h>
+
#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/data_out.h>
#include <deal.II/particles/particle_handler.h>
#include "../tests.h"
particle_handler.insert_particle(particle1, cell1);
particle_handler.insert_particle(particle2, cell2);
- Particles::ParticleOutput<dim, spacedim> particle_viz;
+ Particles::DataOut<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;
+ // 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;
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.
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.
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