From: Rene Gassmoeller Date: Thu, 31 Jan 2019 00:42:20 +0000 (-0800) Subject: Fix invalid memory access in particle properties after transfer X-Git-Tag: v9.1.0-rc1~386^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=7d99efb93bb43ea9927ec1efeae9f0362ae2ad0a;p=dealii.git Fix invalid memory access in particle properties after transfer --- diff --git a/include/deal.II/particles/particle.h b/include/deal.II/particles/particle.h index b4d3de0581..df6ebd1741 100644 --- a/include/deal.II/particles/particle.h +++ b/include/deal.II/particles/particle.h @@ -256,7 +256,12 @@ namespace Particles set_properties(const ArrayView &new_properties); /** - * Get write-access to properties of this particle. + * Get write-access to properties of this particle. If the + * particle has no properties yet, but has access to a + * PropertyPool object it will allocate properties to + * allow writing into them. If it has no properties and + * has no access to a PropertyPool this function will + * throw an exception. * * @return An ArrayView of the properties of this particle. */ @@ -264,7 +269,8 @@ namespace Particles get_properties(); /** - * Get read-access to properties of this particle. + * Get read-access to properties of this particle. If the particle + * has no properties this function throws an exception. * * @return An ArrayView of the properties of this particle. */ diff --git a/source/particles/particle.cc b/source/particles/particle.cc index b5dea11fa8..03f90811b3 100644 --- a/source/particles/particle.cc +++ b/source/particles/particle.cc @@ -284,6 +284,8 @@ namespace Particles Particle::set_properties( const ArrayView &new_properties) { + Assert(property_pool != nullptr, ExcInternalError()); + if (properties == PropertyPool::invalid_handle) properties = property_pool->allocate_properties_array(); @@ -313,7 +315,7 @@ namespace Particles const ArrayView Particle::get_properties() const { - Assert(property_pool != nullptr, ExcInternalError()); + Assert(has_properties(), ExcInternalError()); return property_pool->get_properties(properties); } @@ -326,6 +328,9 @@ namespace Particles { Assert(property_pool != nullptr, ExcInternalError()); + if (properties == PropertyPool::invalid_handle) + properties = property_pool->allocate_properties_array(); + return property_pool->get_properties(properties); } diff --git a/source/particles/particle_handler.cc b/source/particles/particle_handler.cc index 84f61ff4da..af3b69a031 100644 --- a/source/particles/particle_handler.cc +++ b/source/particles/particle_handler.cc @@ -1200,6 +1200,12 @@ namespace Particles data_range.end(), /*allow_compression=*/true); + // Update the reference to the current property pool for all particles. + // This was not stored, because they might be transported across process + // domains. + for (auto &particle : loaded_particles_on_cell) + particle.set_property_pool(*property_pool); + switch (status) { case parallel::distributed::Triangulation::CELL_PERSIST: diff --git a/tests/particles/particle_handler_08.cc b/tests/particles/particle_handler_08.cc new file mode 100644 index 0000000000..d5cd3be2b6 --- /dev/null +++ b/tests/particles/particle_handler_08.cc @@ -0,0 +1,183 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 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.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + + + +// like particle_handler_05, but check that the properties of particles are +// correctly transferred during grid +// refinement/coarsening. Note that this in particular tests for a bug, where +// references to the property pool where lost during particle transfer. + +#include + +#include + +#include +#include + +#include + +#include "../tests.h" + +template +void +create_regular_particle_distribution( + Particles::ParticleHandler & particle_handler, + const parallel::distributed::Triangulation &tr, + const unsigned int particles_per_direction = 3) +{ + for (unsigned int i = 0; i < particles_per_direction; ++i) + for (unsigned int j = 0; j < particles_per_direction; ++j) + { + Point position; + Point reference_position; + unsigned int id = i * particles_per_direction + j; + + position[0] = static_cast(i) / + static_cast(particles_per_direction - 1); + position[1] = static_cast(j) / + static_cast(particles_per_direction - 1); + + if (dim > 2) + for (unsigned int k = 0; k < particles_per_direction; ++k) + { + position[2] = static_cast(j) / + static_cast(particles_per_direction - 1); + id = i * particles_per_direction * particles_per_direction + + j * particles_per_direction + k; + Particles::Particle particle(position, + reference_position, + id); + + typename parallel::distributed::Triangulation:: + active_cell_iterator cell = + GridTools::find_active_cell_around_point( + tr, particle.get_location()); + + Particles::ParticleIterator pit = + particle_handler.insert_particle(particle, cell); + + for (unsigned int i = 0; i < spacedim; ++i) + pit->get_properties()[i] = pit->get_location()[i]; + } + else + { + Particles::Particle particle(position, + reference_position, + id); + + typename parallel::distributed::Triangulation:: + active_cell_iterator cell = + GridTools::find_active_cell_around_point( + tr, particle.get_location()); + + Particles::ParticleIterator pit = + particle_handler.insert_particle(particle, cell); + + for (unsigned int i = 0; i < spacedim; ++i) + pit->get_properties()[i] = pit->get_location()[i]; + } + } +} + +template +void +test() +{ + { + parallel::distributed::Triangulation tr(MPI_COMM_WORLD); + + GridGenerator::hyper_cube(tr); + MappingQ mapping(1); + + const unsigned int n_properties = spacedim; + Particles::ParticleHandler particle_handler(tr, + mapping, + n_properties); + + // TODO: Move this into the Particle handler class. Unfortunately, there are + // some interactions with the SolutionTransfer class that prevent us from + // doing this at the moment. When doing this, check that transferring a + // solution and particles during the same refinement is possible (in + // particular that the order of serialization/deserialization is preserved). + tr.signals.pre_distributed_refinement.connect(std::bind( + &Particles::ParticleHandler::register_store_callback_function, + &particle_handler)); + + tr.signals.post_distributed_refinement.connect(std::bind( + &Particles::ParticleHandler::register_load_callback_function, + &particle_handler, + false)); + + create_regular_particle_distribution(particle_handler, tr); + + for (auto particle = particle_handler.begin(); + particle != particle_handler.end(); + ++particle) + deallog << "Before refinement particle id " << particle->get_id() + << " has first property " << particle->get_properties()[0] + << " and second property " << particle->get_properties()[1] + << std::endl; + + // Check that all particles are moved to children + tr.refine_global(1); + + for (auto particle = particle_handler.begin(); + particle != particle_handler.end(); + ++particle) + deallog << "After refinement particle id " << particle->get_id() + << " has first property " << particle->get_properties()[0] + << " and second property " << particle->get_properties()[1] + << std::endl; + + // Reverse the refinement and check again + for (auto cell = tr.begin_active(); cell != tr.end(); ++cell) + cell->set_coarsen_flag(); + + tr.execute_coarsening_and_refinement(); + + for (auto particle = particle_handler.begin(); + particle != particle_handler.end(); + ++particle) + deallog << "After coarsening particle id " << particle->get_id() + << " has first property " << particle->get_properties()[0] + << " and second property " << particle->get_properties()[1] + << 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(); +} diff --git a/tests/particles/particle_handler_08.with_p4est=true.output b/tests/particles/particle_handler_08.with_p4est=true.output new file mode 100644 index 0000000000..f54a16d67e --- /dev/null +++ b/tests/particles/particle_handler_08.with_p4est=true.output @@ -0,0 +1,139 @@ + +DEAL:0:2d/2d::Before refinement particle id 0 has first property 0.00000 and second property 0.00000 +DEAL:0:2d/2d::Before refinement particle id 1 has first property 0.00000 and second property 0.500000 +DEAL:0:2d/2d::Before refinement particle id 2 has first property 0.00000 and second property 1.00000 +DEAL:0:2d/2d::Before refinement particle id 3 has first property 0.500000 and second property 0.00000 +DEAL:0:2d/2d::Before refinement particle id 4 has first property 0.500000 and second property 0.500000 +DEAL:0:2d/2d::Before refinement particle id 5 has first property 0.500000 and second property 1.00000 +DEAL:0:2d/2d::Before refinement particle id 6 has first property 1.00000 and second property 0.00000 +DEAL:0:2d/2d::Before refinement particle id 7 has first property 1.00000 and second property 0.500000 +DEAL:0:2d/2d::Before refinement particle id 8 has first property 1.00000 and second property 1.00000 +DEAL:0:2d/2d::After refinement particle id 0 has first property 0.00000 and second property 0.00000 +DEAL:0:2d/2d::After refinement particle id 1 has first property 0.00000 and second property 0.500000 +DEAL:0:2d/2d::After refinement particle id 3 has first property 0.500000 and second property 0.00000 +DEAL:0:2d/2d::After refinement particle id 4 has first property 0.500000 and second property 0.500000 +DEAL:0:2d/2d::After refinement particle id 6 has first property 1.00000 and second property 0.00000 +DEAL:0:2d/2d::After refinement particle id 7 has first property 1.00000 and second property 0.500000 +DEAL:0:2d/2d::After refinement particle id 2 has first property 0.00000 and second property 1.00000 +DEAL:0:2d/2d::After refinement particle id 5 has first property 0.500000 and second property 1.00000 +DEAL:0:2d/2d::After refinement particle id 8 has first property 1.00000 and second property 1.00000 +DEAL:0:2d/2d::After coarsening particle id 0 has first property 0.00000 and second property 0.00000 +DEAL:0:2d/2d::After coarsening particle id 1 has first property 0.00000 and second property 0.500000 +DEAL:0:2d/2d::After coarsening particle id 3 has first property 0.500000 and second property 0.00000 +DEAL:0:2d/2d::After coarsening particle id 4 has first property 0.500000 and second property 0.500000 +DEAL:0:2d/2d::After coarsening particle id 6 has first property 1.00000 and second property 0.00000 +DEAL:0:2d/2d::After coarsening particle id 7 has first property 1.00000 and second property 0.500000 +DEAL:0:2d/2d::After coarsening particle id 2 has first property 0.00000 and second property 1.00000 +DEAL:0:2d/2d::After coarsening particle id 5 has first property 0.500000 and second property 1.00000 +DEAL:0:2d/2d::After coarsening particle id 8 has first property 1.00000 and second property 1.00000 +DEAL:0:2d/2d::OK +DEAL:0:2d/3d::Before refinement particle id 0 has first property 0.00000 and second property 0.00000 +DEAL:0:2d/3d::Before refinement particle id 1 has first property 0.00000 and second property 0.500000 +DEAL:0:2d/3d::Before refinement particle id 2 has first property 0.00000 and second property 1.00000 +DEAL:0:2d/3d::Before refinement particle id 3 has first property 0.500000 and second property 0.00000 +DEAL:0:2d/3d::Before refinement particle id 4 has first property 0.500000 and second property 0.500000 +DEAL:0:2d/3d::Before refinement particle id 5 has first property 0.500000 and second property 1.00000 +DEAL:0:2d/3d::Before refinement particle id 6 has first property 1.00000 and second property 0.00000 +DEAL:0:2d/3d::Before refinement particle id 7 has first property 1.00000 and second property 0.500000 +DEAL:0:2d/3d::Before refinement particle id 8 has first property 1.00000 and second property 1.00000 +DEAL:0:2d/3d::After refinement particle id 0 has first property 0.00000 and second property 0.00000 +DEAL:0:2d/3d::After refinement particle id 1 has first property 0.00000 and second property 0.500000 +DEAL:0:2d/3d::After refinement particle id 3 has first property 0.500000 and second property 0.00000 +DEAL:0:2d/3d::After refinement particle id 4 has first property 0.500000 and second property 0.500000 +DEAL:0:2d/3d::After refinement particle id 6 has first property 1.00000 and second property 0.00000 +DEAL:0:2d/3d::After refinement particle id 7 has first property 1.00000 and second property 0.500000 +DEAL:0:2d/3d::After refinement particle id 2 has first property 0.00000 and second property 1.00000 +DEAL:0:2d/3d::After refinement particle id 5 has first property 0.500000 and second property 1.00000 +DEAL:0:2d/3d::After refinement particle id 8 has first property 1.00000 and second property 1.00000 +DEAL:0:2d/3d::After coarsening particle id 0 has first property 0.00000 and second property 0.00000 +DEAL:0:2d/3d::After coarsening particle id 1 has first property 0.00000 and second property 0.500000 +DEAL:0:2d/3d::After coarsening particle id 3 has first property 0.500000 and second property 0.00000 +DEAL:0:2d/3d::After coarsening particle id 4 has first property 0.500000 and second property 0.500000 +DEAL:0:2d/3d::After coarsening particle id 6 has first property 1.00000 and second property 0.00000 +DEAL:0:2d/3d::After coarsening particle id 7 has first property 1.00000 and second property 0.500000 +DEAL:0:2d/3d::After coarsening particle id 2 has first property 0.00000 and second property 1.00000 +DEAL:0:2d/3d::After coarsening particle id 5 has first property 0.500000 and second property 1.00000 +DEAL:0:2d/3d::After coarsening particle id 8 has first property 1.00000 and second property 1.00000 +DEAL:0:2d/3d::OK +DEAL:0:3d/3d::Before refinement particle id 0 has first property 0.00000 and second property 0.00000 +DEAL:0:3d/3d::Before refinement particle id 1 has first property 0.00000 and second property 0.00000 +DEAL:0:3d/3d::Before refinement particle id 2 has first property 0.00000 and second property 0.00000 +DEAL:0:3d/3d::Before refinement particle id 3 has first property 0.00000 and second property 0.500000 +DEAL:0:3d/3d::Before refinement particle id 4 has first property 0.00000 and second property 0.500000 +DEAL:0:3d/3d::Before refinement particle id 5 has first property 0.00000 and second property 0.500000 +DEAL:0:3d/3d::Before refinement particle id 6 has first property 0.00000 and second property 1.00000 +DEAL:0:3d/3d::Before refinement particle id 7 has first property 0.00000 and second property 1.00000 +DEAL:0:3d/3d::Before refinement particle id 8 has first property 0.00000 and second property 1.00000 +DEAL:0:3d/3d::Before refinement particle id 9 has first property 0.500000 and second property 0.00000 +DEAL:0:3d/3d::Before refinement particle id 10 has first property 0.500000 and second property 0.00000 +DEAL:0:3d/3d::Before refinement particle id 11 has first property 0.500000 and second property 0.00000 +DEAL:0:3d/3d::Before refinement particle id 12 has first property 0.500000 and second property 0.500000 +DEAL:0:3d/3d::Before refinement particle id 13 has first property 0.500000 and second property 0.500000 +DEAL:0:3d/3d::Before refinement particle id 14 has first property 0.500000 and second property 0.500000 +DEAL:0:3d/3d::Before refinement particle id 15 has first property 0.500000 and second property 1.00000 +DEAL:0:3d/3d::Before refinement particle id 16 has first property 0.500000 and second property 1.00000 +DEAL:0:3d/3d::Before refinement particle id 17 has first property 0.500000 and second property 1.00000 +DEAL:0:3d/3d::Before refinement particle id 18 has first property 1.00000 and second property 0.00000 +DEAL:0:3d/3d::Before refinement particle id 19 has first property 1.00000 and second property 0.00000 +DEAL:0:3d/3d::Before refinement particle id 20 has first property 1.00000 and second property 0.00000 +DEAL:0:3d/3d::Before refinement particle id 21 has first property 1.00000 and second property 0.500000 +DEAL:0:3d/3d::Before refinement particle id 22 has first property 1.00000 and second property 0.500000 +DEAL:0:3d/3d::Before refinement particle id 23 has first property 1.00000 and second property 0.500000 +DEAL:0:3d/3d::Before refinement particle id 24 has first property 1.00000 and second property 1.00000 +DEAL:0:3d/3d::Before refinement particle id 25 has first property 1.00000 and second property 1.00000 +DEAL:0:3d/3d::Before refinement particle id 26 has first property 1.00000 and second property 1.00000 +DEAL:0:3d/3d::After refinement particle id 0 has first property 0.00000 and second property 0.00000 +DEAL:0:3d/3d::After refinement particle id 1 has first property 0.00000 and second property 0.00000 +DEAL:0:3d/3d::After refinement particle id 2 has first property 0.00000 and second property 0.00000 +DEAL:0:3d/3d::After refinement particle id 3 has first property 0.00000 and second property 0.500000 +DEAL:0:3d/3d::After refinement particle id 4 has first property 0.00000 and second property 0.500000 +DEAL:0:3d/3d::After refinement particle id 5 has first property 0.00000 and second property 0.500000 +DEAL:0:3d/3d::After refinement particle id 9 has first property 0.500000 and second property 0.00000 +DEAL:0:3d/3d::After refinement particle id 10 has first property 0.500000 and second property 0.00000 +DEAL:0:3d/3d::After refinement particle id 11 has first property 0.500000 and second property 0.00000 +DEAL:0:3d/3d::After refinement particle id 12 has first property 0.500000 and second property 0.500000 +DEAL:0:3d/3d::After refinement particle id 13 has first property 0.500000 and second property 0.500000 +DEAL:0:3d/3d::After refinement particle id 14 has first property 0.500000 and second property 0.500000 +DEAL:0:3d/3d::After refinement particle id 18 has first property 1.00000 and second property 0.00000 +DEAL:0:3d/3d::After refinement particle id 19 has first property 1.00000 and second property 0.00000 +DEAL:0:3d/3d::After refinement particle id 20 has first property 1.00000 and second property 0.00000 +DEAL:0:3d/3d::After refinement particle id 21 has first property 1.00000 and second property 0.500000 +DEAL:0:3d/3d::After refinement particle id 22 has first property 1.00000 and second property 0.500000 +DEAL:0:3d/3d::After refinement particle id 23 has first property 1.00000 and second property 0.500000 +DEAL:0:3d/3d::After refinement particle id 6 has first property 0.00000 and second property 1.00000 +DEAL:0:3d/3d::After refinement particle id 7 has first property 0.00000 and second property 1.00000 +DEAL:0:3d/3d::After refinement particle id 8 has first property 0.00000 and second property 1.00000 +DEAL:0:3d/3d::After refinement particle id 15 has first property 0.500000 and second property 1.00000 +DEAL:0:3d/3d::After refinement particle id 16 has first property 0.500000 and second property 1.00000 +DEAL:0:3d/3d::After refinement particle id 17 has first property 0.500000 and second property 1.00000 +DEAL:0:3d/3d::After refinement particle id 24 has first property 1.00000 and second property 1.00000 +DEAL:0:3d/3d::After refinement particle id 25 has first property 1.00000 and second property 1.00000 +DEAL:0:3d/3d::After refinement particle id 26 has first property 1.00000 and second property 1.00000 +DEAL:0:3d/3d::After coarsening particle id 0 has first property 0.00000 and second property 0.00000 +DEAL:0:3d/3d::After coarsening particle id 1 has first property 0.00000 and second property 0.00000 +DEAL:0:3d/3d::After coarsening particle id 2 has first property 0.00000 and second property 0.00000 +DEAL:0:3d/3d::After coarsening particle id 3 has first property 0.00000 and second property 0.500000 +DEAL:0:3d/3d::After coarsening particle id 4 has first property 0.00000 and second property 0.500000 +DEAL:0:3d/3d::After coarsening particle id 5 has first property 0.00000 and second property 0.500000 +DEAL:0:3d/3d::After coarsening particle id 9 has first property 0.500000 and second property 0.00000 +DEAL:0:3d/3d::After coarsening particle id 10 has first property 0.500000 and second property 0.00000 +DEAL:0:3d/3d::After coarsening particle id 11 has first property 0.500000 and second property 0.00000 +DEAL:0:3d/3d::After coarsening particle id 12 has first property 0.500000 and second property 0.500000 +DEAL:0:3d/3d::After coarsening particle id 13 has first property 0.500000 and second property 0.500000 +DEAL:0:3d/3d::After coarsening particle id 14 has first property 0.500000 and second property 0.500000 +DEAL:0:3d/3d::After coarsening particle id 18 has first property 1.00000 and second property 0.00000 +DEAL:0:3d/3d::After coarsening particle id 19 has first property 1.00000 and second property 0.00000 +DEAL:0:3d/3d::After coarsening particle id 20 has first property 1.00000 and second property 0.00000 +DEAL:0:3d/3d::After coarsening particle id 21 has first property 1.00000 and second property 0.500000 +DEAL:0:3d/3d::After coarsening particle id 22 has first property 1.00000 and second property 0.500000 +DEAL:0:3d/3d::After coarsening particle id 23 has first property 1.00000 and second property 0.500000 +DEAL:0:3d/3d::After coarsening particle id 6 has first property 0.00000 and second property 1.00000 +DEAL:0:3d/3d::After coarsening particle id 7 has first property 0.00000 and second property 1.00000 +DEAL:0:3d/3d::After coarsening particle id 8 has first property 0.00000 and second property 1.00000 +DEAL:0:3d/3d::After coarsening particle id 15 has first property 0.500000 and second property 1.00000 +DEAL:0:3d/3d::After coarsening particle id 16 has first property 0.500000 and second property 1.00000 +DEAL:0:3d/3d::After coarsening particle id 17 has first property 0.500000 and second property 1.00000 +DEAL:0:3d/3d::After coarsening particle id 24 has first property 1.00000 and second property 1.00000 +DEAL:0:3d/3d::After coarsening particle id 25 has first property 1.00000 and second property 1.00000 +DEAL:0:3d/3d::After coarsening particle id 26 has first property 1.00000 and second property 1.00000 +DEAL:0:3d/3d::OK