const Mapping<dim, spacedim> & mapping,
const unsigned int n_properties = 0);
+ /**
+ * Copy the state of particle handler @p particle_handler into the
+ * current object. This will copy
+ * all particles and properties and leave this object
+ * as an identical copy of @p particle_handler. Existing
+ * particles in this object are deleted. Be aware that this
+ * does not copy functions that are connected to the signals of
+ * @p particle_handler, nor does it connect the current object's member
+ * functions to triangulation signals, which must be done by the caller
+ * if necessary, that is if the @p particle_handler had
+ * connected functions.
+ *
+ * This function is expensive as it has to duplicate all data
+ * in @p particle_handler, and insert it into this object,
+ * which may be a significant amount of data. However, it can
+ * be useful to save the state of a particle
+ * collection at a certain point in time and reset this
+ * state later under certain conditions, for example if
+ * a timestep has to be undone and repeated.
+ */
+ void
+ copy_from(const ParticleHandler<dim, spacedim> &particle_handler);
+
/**
* Clear all particle related data.
*/
+ template <int dim, int spacedim>
+ void
+ ParticleHandler<dim, spacedim>::copy_from(
+ const ParticleHandler<dim, spacedim> &particle_handler)
+ {
+ // clear and initialize this object before copying particles
+ clear();
+ const unsigned int n_properties =
+ particle_handler.property_pool->n_properties_per_slot();
+ initialize(*particle_handler.triangulation,
+ *particle_handler.mapping,
+ n_properties);
+
+ // copy static members
+ global_number_of_particles = particle_handler.global_number_of_particles;
+ global_max_particles_per_cell =
+ particle_handler.global_max_particles_per_cell;
+ next_free_particle_index = particle_handler.next_free_particle_index;
+ particles = particle_handler.particles;
+ ghost_particles = particle_handler.ghost_particles;
+ handle = particle_handler.handle;
+
+ // copy dynamic properties
+ auto from_particle = particle_handler.begin();
+ for (auto &particle : *this)
+ {
+ particle.set_property_pool(*property_pool);
+ particle.set_properties(from_particle->get_properties());
+ ++from_particle;
+ }
+
+ auto from_ghost = particle_handler.begin_ghost();
+ for (auto ghost = begin_ghost(); ghost != end_ghost();
+ ++ghost, ++from_ghost)
+ {
+ ghost->set_property_pool(*property_pool);
+ ghost->set_properties(from_ghost->get_properties());
+ }
+ }
+
+
+
template <int dim, int spacedim>
void
ParticleHandler<dim, spacedim>::clear()
--- /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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// like particle_handler_08, but check that the properties of particles are
+// correctly copied in the ParticleHandler::copy_from function.
+
+#include <deal.II/base/quadrature_lib.h>
+
+#include <deal.II/fe/mapping_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/particles/generators.h>
+#include <deal.II/particles/particle_handler.h>
+
+#include "../tests.h"
+
+template <int dim, int spacedim>
+void
+test()
+{
+ {
+ Triangulation<dim, spacedim> tr;
+
+ GridGenerator::hyper_cube(tr);
+ MappingQ<dim, spacedim> mapping(1);
+ const unsigned int n_properties = spacedim;
+ Particles::ParticleHandler<dim, spacedim> particle_handler(tr,
+ mapping,
+ n_properties);
+
+ std::vector<Point<dim>> particle_reference_locations =
+ QGauss<dim>(3).get_points();
+
+ Particles::Generators::regular_reference_locations(
+ tr, particle_reference_locations, particle_handler, mapping);
+
+ for (auto &particle : particle_handler)
+ {
+ particle.get_properties()[spacedim - 1] = particle.get_location()[0];
+ deallog << "Before copying particle id " << particle.get_id()
+ << " has first property " << particle.get_properties()[0]
+ << " and last property "
+ << particle.get_properties()[spacedim - 1] << " and position "
+ << particle.get_location() << std::endl;
+ }
+
+ Particles::ParticleHandler<dim, spacedim> particle_handler_copy;
+ particle_handler_copy.copy_from(particle_handler);
+
+ for (const auto &particle : particle_handler_copy)
+ deallog << "After copying particle id " << particle.get_id()
+ << " has first property " << particle.get_properties()[0]
+ << " and last property "
+ << particle.get_properties()[spacedim - 1] << " and position "
+ << particle.get_location() << 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
+
+DEAL:0:2d/2d::Before copying particle id 0 has first property 0.00000 and last property 0.112702 and position 0.112702 0.112702
+DEAL:0:2d/2d::Before copying particle id 1 has first property 0.00000 and last property 0.500000 and position 0.500000 0.112702
+DEAL:0:2d/2d::Before copying particle id 2 has first property 0.00000 and last property 0.887298 and position 0.887298 0.112702
+DEAL:0:2d/2d::Before copying particle id 3 has first property 0.00000 and last property 0.112702 and position 0.112702 0.500000
+DEAL:0:2d/2d::Before copying particle id 4 has first property 0.00000 and last property 0.500000 and position 0.500000 0.500000
+DEAL:0:2d/2d::Before copying particle id 5 has first property 0.00000 and last property 0.887298 and position 0.887298 0.500000
+DEAL:0:2d/2d::Before copying particle id 6 has first property 0.00000 and last property 0.112702 and position 0.112702 0.887298
+DEAL:0:2d/2d::Before copying particle id 7 has first property 0.00000 and last property 0.500000 and position 0.500000 0.887298
+DEAL:0:2d/2d::Before copying particle id 8 has first property 0.00000 and last property 0.887298 and position 0.887298 0.887298
+DEAL:0:2d/2d::After copying particle id 0 has first property 0.00000 and last property 0.112702 and position 0.112702 0.112702
+DEAL:0:2d/2d::After copying particle id 1 has first property 0.00000 and last property 0.500000 and position 0.500000 0.112702
+DEAL:0:2d/2d::After copying particle id 2 has first property 0.00000 and last property 0.887298 and position 0.887298 0.112702
+DEAL:0:2d/2d::After copying particle id 3 has first property 0.00000 and last property 0.112702 and position 0.112702 0.500000
+DEAL:0:2d/2d::After copying particle id 4 has first property 0.00000 and last property 0.500000 and position 0.500000 0.500000
+DEAL:0:2d/2d::After copying particle id 5 has first property 0.00000 and last property 0.887298 and position 0.887298 0.500000
+DEAL:0:2d/2d::After copying particle id 6 has first property 0.00000 and last property 0.112702 and position 0.112702 0.887298
+DEAL:0:2d/2d::After copying particle id 7 has first property 0.00000 and last property 0.500000 and position 0.500000 0.887298
+DEAL:0:2d/2d::After copying particle id 8 has first property 0.00000 and last property 0.887298 and position 0.887298 0.887298
+DEAL:0:2d/2d::OK
+DEAL:0:2d/3d::Before copying particle id 0 has first property 0.00000 and last property 0.112702 and position 0.112702 0.112702 0.00000
+DEAL:0:2d/3d::Before copying particle id 1 has first property 0.00000 and last property 0.500000 and position 0.500000 0.112702 0.00000
+DEAL:0:2d/3d::Before copying particle id 2 has first property 0.00000 and last property 0.887298 and position 0.887298 0.112702 0.00000
+DEAL:0:2d/3d::Before copying particle id 3 has first property 0.00000 and last property 0.112702 and position 0.112702 0.500000 0.00000
+DEAL:0:2d/3d::Before copying particle id 4 has first property 0.00000 and last property 0.500000 and position 0.500000 0.500000 0.00000
+DEAL:0:2d/3d::Before copying particle id 5 has first property 0.00000 and last property 0.887298 and position 0.887298 0.500000 0.00000
+DEAL:0:2d/3d::Before copying particle id 6 has first property 0.00000 and last property 0.112702 and position 0.112702 0.887298 0.00000
+DEAL:0:2d/3d::Before copying particle id 7 has first property 0.00000 and last property 0.500000 and position 0.500000 0.887298 0.00000
+DEAL:0:2d/3d::Before copying particle id 8 has first property 0.00000 and last property 0.887298 and position 0.887298 0.887298 0.00000
+DEAL:0:2d/3d::After copying particle id 0 has first property 0.00000 and last property 0.112702 and position 0.112702 0.112702 0.00000
+DEAL:0:2d/3d::After copying particle id 1 has first property 0.00000 and last property 0.500000 and position 0.500000 0.112702 0.00000
+DEAL:0:2d/3d::After copying particle id 2 has first property 0.00000 and last property 0.887298 and position 0.887298 0.112702 0.00000
+DEAL:0:2d/3d::After copying particle id 3 has first property 0.00000 and last property 0.112702 and position 0.112702 0.500000 0.00000
+DEAL:0:2d/3d::After copying particle id 4 has first property 0.00000 and last property 0.500000 and position 0.500000 0.500000 0.00000
+DEAL:0:2d/3d::After copying particle id 5 has first property 0.00000 and last property 0.887298 and position 0.887298 0.500000 0.00000
+DEAL:0:2d/3d::After copying particle id 6 has first property 0.00000 and last property 0.112702 and position 0.112702 0.887298 0.00000
+DEAL:0:2d/3d::After copying particle id 7 has first property 0.00000 and last property 0.500000 and position 0.500000 0.887298 0.00000
+DEAL:0:2d/3d::After copying particle id 8 has first property 0.00000 and last property 0.887298 and position 0.887298 0.887298 0.00000
+DEAL:0:2d/3d::OK
+DEAL:0:3d/3d::Before copying particle id 0 has first property 0.00000 and last property 0.112702 and position 0.112702 0.112702 0.112702
+DEAL:0:3d/3d::Before copying particle id 1 has first property 0.00000 and last property 0.500000 and position 0.500000 0.112702 0.112702
+DEAL:0:3d/3d::Before copying particle id 2 has first property 0.00000 and last property 0.887298 and position 0.887298 0.112702 0.112702
+DEAL:0:3d/3d::Before copying particle id 3 has first property 0.00000 and last property 0.112702 and position 0.112702 0.500000 0.112702
+DEAL:0:3d/3d::Before copying particle id 4 has first property 0.00000 and last property 0.500000 and position 0.500000 0.500000 0.112702
+DEAL:0:3d/3d::Before copying particle id 5 has first property 0.00000 and last property 0.887298 and position 0.887298 0.500000 0.112702
+DEAL:0:3d/3d::Before copying particle id 6 has first property 0.00000 and last property 0.112702 and position 0.112702 0.887298 0.112702
+DEAL:0:3d/3d::Before copying particle id 7 has first property 0.00000 and last property 0.500000 and position 0.500000 0.887298 0.112702
+DEAL:0:3d/3d::Before copying particle id 8 has first property 0.00000 and last property 0.887298 and position 0.887298 0.887298 0.112702
+DEAL:0:3d/3d::Before copying particle id 9 has first property 0.00000 and last property 0.112702 and position 0.112702 0.112702 0.500000
+DEAL:0:3d/3d::Before copying particle id 10 has first property 0.00000 and last property 0.500000 and position 0.500000 0.112702 0.500000
+DEAL:0:3d/3d::Before copying particle id 11 has first property 0.00000 and last property 0.887298 and position 0.887298 0.112702 0.500000
+DEAL:0:3d/3d::Before copying particle id 12 has first property 0.00000 and last property 0.112702 and position 0.112702 0.500000 0.500000
+DEAL:0:3d/3d::Before copying particle id 13 has first property 0.00000 and last property 0.500000 and position 0.500000 0.500000 0.500000
+DEAL:0:3d/3d::Before copying particle id 14 has first property 0.00000 and last property 0.887298 and position 0.887298 0.500000 0.500000
+DEAL:0:3d/3d::Before copying particle id 15 has first property 0.00000 and last property 0.112702 and position 0.112702 0.887298 0.500000
+DEAL:0:3d/3d::Before copying particle id 16 has first property 0.00000 and last property 0.500000 and position 0.500000 0.887298 0.500000
+DEAL:0:3d/3d::Before copying particle id 17 has first property 0.00000 and last property 0.887298 and position 0.887298 0.887298 0.500000
+DEAL:0:3d/3d::Before copying particle id 18 has first property 0.00000 and last property 0.112702 and position 0.112702 0.112702 0.887298
+DEAL:0:3d/3d::Before copying particle id 19 has first property 0.00000 and last property 0.500000 and position 0.500000 0.112702 0.887298
+DEAL:0:3d/3d::Before copying particle id 20 has first property 0.00000 and last property 0.887298 and position 0.887298 0.112702 0.887298
+DEAL:0:3d/3d::Before copying particle id 21 has first property 0.00000 and last property 0.112702 and position 0.112702 0.500000 0.887298
+DEAL:0:3d/3d::Before copying particle id 22 has first property 0.00000 and last property 0.500000 and position 0.500000 0.500000 0.887298
+DEAL:0:3d/3d::Before copying particle id 23 has first property 0.00000 and last property 0.887298 and position 0.887298 0.500000 0.887298
+DEAL:0:3d/3d::Before copying particle id 24 has first property 0.00000 and last property 0.112702 and position 0.112702 0.887298 0.887298
+DEAL:0:3d/3d::Before copying particle id 25 has first property 0.00000 and last property 0.500000 and position 0.500000 0.887298 0.887298
+DEAL:0:3d/3d::Before copying particle id 26 has first property 0.00000 and last property 0.887298 and position 0.887298 0.887298 0.887298
+DEAL:0:3d/3d::After copying particle id 0 has first property 0.00000 and last property 0.112702 and position 0.112702 0.112702 0.112702
+DEAL:0:3d/3d::After copying particle id 1 has first property 0.00000 and last property 0.500000 and position 0.500000 0.112702 0.112702
+DEAL:0:3d/3d::After copying particle id 2 has first property 0.00000 and last property 0.887298 and position 0.887298 0.112702 0.112702
+DEAL:0:3d/3d::After copying particle id 3 has first property 0.00000 and last property 0.112702 and position 0.112702 0.500000 0.112702
+DEAL:0:3d/3d::After copying particle id 4 has first property 0.00000 and last property 0.500000 and position 0.500000 0.500000 0.112702
+DEAL:0:3d/3d::After copying particle id 5 has first property 0.00000 and last property 0.887298 and position 0.887298 0.500000 0.112702
+DEAL:0:3d/3d::After copying particle id 6 has first property 0.00000 and last property 0.112702 and position 0.112702 0.887298 0.112702
+DEAL:0:3d/3d::After copying particle id 7 has first property 0.00000 and last property 0.500000 and position 0.500000 0.887298 0.112702
+DEAL:0:3d/3d::After copying particle id 8 has first property 0.00000 and last property 0.887298 and position 0.887298 0.887298 0.112702
+DEAL:0:3d/3d::After copying particle id 9 has first property 0.00000 and last property 0.112702 and position 0.112702 0.112702 0.500000
+DEAL:0:3d/3d::After copying particle id 10 has first property 0.00000 and last property 0.500000 and position 0.500000 0.112702 0.500000
+DEAL:0:3d/3d::After copying particle id 11 has first property 0.00000 and last property 0.887298 and position 0.887298 0.112702 0.500000
+DEAL:0:3d/3d::After copying particle id 12 has first property 0.00000 and last property 0.112702 and position 0.112702 0.500000 0.500000
+DEAL:0:3d/3d::After copying particle id 13 has first property 0.00000 and last property 0.500000 and position 0.500000 0.500000 0.500000
+DEAL:0:3d/3d::After copying particle id 14 has first property 0.00000 and last property 0.887298 and position 0.887298 0.500000 0.500000
+DEAL:0:3d/3d::After copying particle id 15 has first property 0.00000 and last property 0.112702 and position 0.112702 0.887298 0.500000
+DEAL:0:3d/3d::After copying particle id 16 has first property 0.00000 and last property 0.500000 and position 0.500000 0.887298 0.500000
+DEAL:0:3d/3d::After copying particle id 17 has first property 0.00000 and last property 0.887298 and position 0.887298 0.887298 0.500000
+DEAL:0:3d/3d::After copying particle id 18 has first property 0.00000 and last property 0.112702 and position 0.112702 0.112702 0.887298
+DEAL:0:3d/3d::After copying particle id 19 has first property 0.00000 and last property 0.500000 and position 0.500000 0.112702 0.887298
+DEAL:0:3d/3d::After copying particle id 20 has first property 0.00000 and last property 0.887298 and position 0.887298 0.112702 0.887298
+DEAL:0:3d/3d::After copying particle id 21 has first property 0.00000 and last property 0.112702 and position 0.112702 0.500000 0.887298
+DEAL:0:3d/3d::After copying particle id 22 has first property 0.00000 and last property 0.500000 and position 0.500000 0.500000 0.887298
+DEAL:0:3d/3d::After copying particle id 23 has first property 0.00000 and last property 0.887298 and position 0.887298 0.500000 0.887298
+DEAL:0:3d/3d::After copying particle id 24 has first property 0.00000 and last property 0.112702 and position 0.112702 0.887298 0.887298
+DEAL:0:3d/3d::After copying particle id 25 has first property 0.00000 and last property 0.500000 and position 0.500000 0.887298 0.887298
+DEAL:0:3d/3d::After copying particle id 26 has first property 0.00000 and last property 0.887298 and position 0.887298 0.887298 0.887298
+DEAL:0:3d/3d::OK