From: Rene Gassmoeller Date: Tue, 22 Jun 2021 14:31:27 +0000 (-0400) Subject: Clarify locally owned particles X-Git-Tag: v9.4.0-rc1~1200^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=445722bb53fe892ef2a60b1aaf4a4393447b8c66;p=dealii.git Clarify locally owned particles --- diff --git a/include/deal.II/particles/particle_handler.h b/include/deal.II/particles/particle_handler.h index a674d73d54..865d94b205 100644 --- a/include/deal.II/particles/particle_handler.h +++ b/include/deal.II/particles/particle_handler.h @@ -830,10 +830,10 @@ namespace Particles types::particle_index global_number_of_particles; /** - * This variable stores how many particles are stored locally. It is + * This variable stores how many particles are locally owned. It is * calculated by update_cached_numbers(). */ - types::particle_index local_number_of_particles; + types::particle_index number_of_locally_owned_particles; /** * The maximum number of particles per cell in the global domain. This diff --git a/include/deal.II/particles/property_pool.h b/include/deal.II/particles/property_pool.h index 3f74a8167c..944d2ca643 100644 --- a/include/deal.II/particles/property_pool.h +++ b/include/deal.II/particles/property_pool.h @@ -211,6 +211,12 @@ namespace Particles unsigned int n_properties_per_slot() const; + /** + * Return how many slots are currently registered in the pool. + */ + unsigned int + n_registered_slots() const; + /** * This function makes sure that all internally stored memory blocks * are sorted in the same order as one would loop over the @p handles_to_sort diff --git a/source/particles/particle_handler.cc b/source/particles/particle_handler.cc index 2269b73c74..a8cb1e29b3 100644 --- a/source/particles/particle_handler.cc +++ b/source/particles/particle_handler.cc @@ -58,7 +58,7 @@ namespace Particles , property_pool(std::make_unique>(0)) , particles() , global_number_of_particles(0) - , local_number_of_particles(0) + , number_of_locally_owned_particles(0) , global_max_particles_per_cell(0) , next_free_particle_index(0) , size_callback() @@ -79,7 +79,7 @@ namespace Particles , property_pool(std::make_unique>(n_properties)) , particles(triangulation.n_active_cells()) , global_number_of_particles(0) - , local_number_of_particles(0) + , number_of_locally_owned_particles(0) , global_max_particles_per_cell(0) , next_free_particle_index(0) , size_callback() @@ -167,7 +167,9 @@ namespace Particles // copy static members global_number_of_particles = particle_handler.global_number_of_particles; - local_number_of_particles = particle_handler.local_number_of_particles; + number_of_locally_owned_particles = + particle_handler.number_of_locally_owned_particles; + global_max_particles_per_cell = particle_handler.global_max_particles_per_cell; next_free_particle_index = particle_handler.next_free_particle_index; @@ -185,10 +187,10 @@ namespace Particles ParticleHandler::clear() { clear_particles(); - global_number_of_particles = 0; - local_number_of_particles = 0; - next_free_particle_index = 0; - global_max_particles_per_cell = 0; + global_number_of_particles = 0; + number_of_locally_owned_particles = 0; + next_free_particle_index = 0; + global_max_particles_per_cell = 0; } @@ -217,18 +219,22 @@ namespace Particles void ParticleHandler::update_cached_numbers() { - local_number_of_particles = 0; + number_of_locally_owned_particles = 0; types::particle_index local_max_particle_index = 0; types::particle_index local_max_particles_per_cell = 0; - for (const auto &particles_in_cell : particles) + for (const auto &cell : triangulation->active_cell_iterators()) { + const auto &particles_in_cell = particles[cell->active_cell_index()]; + const types::particle_index n_particles_in_cell = particles_in_cell.size(); local_max_particles_per_cell = std::max(local_max_particles_per_cell, n_particles_in_cell); - local_number_of_particles += n_particles_in_cell; + + if (cell->is_locally_owned()) + number_of_locally_owned_particles += n_particles_in_cell; for (const auto &particle : particles_in_cell) { @@ -243,7 +249,7 @@ namespace Particles &*triangulation)) { global_number_of_particles = dealii::Utilities::MPI::sum( - local_number_of_particles, + number_of_locally_owned_particles, parallel_triangulation->get_communicator()); if (global_number_of_particles == 0) @@ -267,7 +273,7 @@ namespace Particles } else { - global_number_of_particles = local_number_of_particles; + global_number_of_particles = number_of_locally_owned_particles; next_free_particle_index = global_number_of_particles == 0 ? 0 : local_max_particle_index + 1; global_max_particles_per_cell = local_max_particles_per_cell; @@ -368,6 +374,12 @@ namespace Particles property_pool->deregister_particle(handle); } + // need to reduce the cached number before deleting, because the iterator + // may be invalid after removing the particle even if only + // accessing the cell + if (particle->get_surrounding_cell()->is_locally_owned()) + --number_of_locally_owned_particles; + if (particles_on_cell.size() > 1) { particles_on_cell[particle->particle_index_within_cell] = @@ -378,8 +390,6 @@ namespace Particles { particles_on_cell.clear(); } - - --local_number_of_particles; } @@ -480,7 +490,7 @@ namespace Particles data = particle_it->read_particle_data_from_memory(data); - ++local_number_of_particles; + ++number_of_locally_owned_particles; return particle_it; } @@ -519,7 +529,7 @@ namespace Particles if (properties.size() != 0) particle_it->set_properties(properties); - ++local_number_of_particles; + ++number_of_locally_owned_particles; return particle_it; } @@ -894,7 +904,7 @@ namespace Particles types::particle_index ParticleHandler::n_locally_owned_particles() const { - return local_number_of_particles; + return number_of_locally_owned_particles; } @@ -1296,7 +1306,7 @@ namespace Particles // now make sure particle data is sorted in order of iteration std::vector::Handle> unsorted_handles; - unsorted_handles.reserve(local_number_of_particles); + unsorted_handles.reserve(property_pool->n_registered_slots()); typename PropertyPool::Handle sorted_handle = 0; for (auto &particles_in_cell : particles) @@ -1354,8 +1364,7 @@ namespace Particles parallel_triangulation->ghost_owners(); for (const auto ghost_owner : ghost_owners) ghost_particles_cache.ghost_particles_by_domain[ghost_owner].reserve( - static_cast::size_type>( - local_number_of_particles * 0.25)); + n_locally_owned_particles() / 4); const std::vector> vertex_to_neighbor_subdomain = triangulation_cache->get_vertex_to_neighbor_subdomain(); diff --git a/source/particles/property_pool.cc b/source/particles/property_pool.cc index 1189a664b1..9a03f3195f 100644 --- a/source/particles/property_pool.cc +++ b/source/particles/property_pool.cc @@ -167,6 +167,29 @@ namespace Particles return n_properties; } + + + template + unsigned int + PropertyPool::n_registered_slots() const + { + Assert(locations.size() == reference_locations.size(), + ExcMessage("Number of registered locations is not equal to number " + "of registered reference locations.")); + + Assert(locations.size() == ids.size(), + ExcMessage("Number of registered locations is not equal to number " + "of registered ids.")); + + Assert(locations.size() * n_properties == properties.size(), + ExcMessage("Number of registered locations is not equal to number " + "of registered property slots.")); + + return locations.size() - currently_available_handles.size(); + } + + + template void PropertyPool::sort_memory_slots( diff --git a/tests/particles/exchange_ghosts_01.cc b/tests/particles/exchange_ghosts_01.cc new file mode 100644 index 0000000000..a5dbe169b7 --- /dev/null +++ b/tests/particles/exchange_ghosts_01.cc @@ -0,0 +1,124 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2021 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_06, but tests that the particle handlers function +// n_locally_owned_particles only counts locally owned particles. + +#include + +#include + +#include +#include + +#include + +#include "../tests.h" + +template +void +test() +{ + { + parallel::distributed::Triangulation tr(MPI_COMM_WORLD); + + GridGenerator::hyper_cube(tr); + tr.refine_global(2); + MappingQ mapping(1); + + // both processes create a particle handler with a particle in a + // position that is in a ghost cell to the other process + Particles::ParticleHandler particle_handler(tr, mapping); + + Point position; + Point reference_position; + + if (Utilities::MPI::this_mpi_process(tr.get_communicator()) == 0) + for (unsigned int i = 0; i < dim; ++i) + position(i) = 0.475; + else + for (unsigned int i = 0; i < dim; ++i) + position(i) = 0.525; + + Particles::Particle particle( + position, + reference_position, + Utilities::MPI::this_mpi_process(tr.get_communicator())); + + // We give a local random cell hint to check that sorting and + // transferring ghost particles works. + typename Triangulation::active_cell_iterator cell = + tr.begin_active(); + while (!cell->is_locally_owned()) + ++cell; + + particle_handler.insert_particle(particle, cell); + + particle_handler.sort_particles_into_subdomains_and_cells(); + + deallog << "Before ghost exchange: " + << particle_handler.n_locally_owned_particles() + << " locally owned particles on process " + << Utilities::MPI::this_mpi_process(tr.get_communicator()) + << std::endl; + + deallog << "Before ghost exchange: " + << particle_handler.get_property_pool().n_registered_slots() + << " stored particles on process " + << Utilities::MPI::this_mpi_process(tr.get_communicator()) + << std::endl; + + particle_handler.exchange_ghost_particles(); + + // Usually this update is unnecessary, because exchanging ghost particles + // does not change anything about the cache. However, there was a bug + // that ghost particles were counted as locally owned, so make sure + // this does not happen again. + particle_handler.update_cached_numbers(); + + deallog << "After ghost exchange: " + << particle_handler.n_locally_owned_particles() + << " locally owned particles on process " + << Utilities::MPI::this_mpi_process(tr.get_communicator()) + << std::endl; + + deallog << "After ghost exchange: " + << particle_handler.get_property_pool().n_registered_slots() + << " stored particles on process " + << Utilities::MPI::this_mpi_process(tr.get_communicator()) + << 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("3d/3d"); + test<3, 3>(); + deallog.pop(); +} diff --git a/tests/particles/exchange_ghosts_01.with_p4est=true.mpirun=2.output b/tests/particles/exchange_ghosts_01.with_p4est=true.mpirun=2.output new file mode 100644 index 0000000000..8e729f83b9 --- /dev/null +++ b/tests/particles/exchange_ghosts_01.with_p4est=true.mpirun=2.output @@ -0,0 +1,23 @@ + +DEAL:0:2d/2d::Before ghost exchange: 1 locally owned particles on process 0 +DEAL:0:2d/2d::Before ghost exchange: 1 stored particles on process 0 +DEAL:0:2d/2d::After ghost exchange: 1 locally owned particles on process 0 +DEAL:0:2d/2d::After ghost exchange: 2 stored particles on process 0 +DEAL:0:2d/2d::OK +DEAL:0:3d/3d::Before ghost exchange: 1 locally owned particles on process 0 +DEAL:0:3d/3d::Before ghost exchange: 1 stored particles on process 0 +DEAL:0:3d/3d::After ghost exchange: 1 locally owned particles on process 0 +DEAL:0:3d/3d::After ghost exchange: 2 stored particles on process 0 +DEAL:0:3d/3d::OK + +DEAL:1:2d/2d::Before ghost exchange: 1 locally owned particles on process 1 +DEAL:1:2d/2d::Before ghost exchange: 1 stored particles on process 1 +DEAL:1:2d/2d::After ghost exchange: 1 locally owned particles on process 1 +DEAL:1:2d/2d::After ghost exchange: 2 stored particles on process 1 +DEAL:1:2d/2d::OK +DEAL:1:3d/3d::Before ghost exchange: 1 locally owned particles on process 1 +DEAL:1:3d/3d::Before ghost exchange: 1 stored particles on process 1 +DEAL:1:3d/3d::After ghost exchange: 1 locally owned particles on process 1 +DEAL:1:3d/3d::After ghost exchange: 2 stored particles on process 1 +DEAL:1:3d/3d::OK +