From: Luca Heltai Date: Thu, 11 Jun 2020 10:24:53 +0000 (+0200) Subject: Test also the version with Particles. X-Git-Tag: v9.3.0-rc1~1390^2~1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=bfa97b0d52f509fe7bacd64805673a00f214e3a4;p=dealii.git Test also the version with Particles. --- diff --git a/source/particles/particle_handler.cc b/source/particles/particle_handler.cc index ad1313d228..0f3f4274d3 100644 --- a/source/particles/particle_handler.cc +++ b/source/particles/particle_handler.cc @@ -739,7 +739,7 @@ namespace Particles positions.resize(particles.size()); ids.resize(particles.size()); if (n_properties_per_particle() > 0) - properties.resize(properties.size(), + properties.resize(particles.size(), std::vector(n_properties_per_particle())); unsigned int i = 0; diff --git a/tests/particles/particle_handler_17.cc b/tests/particles/particle_handler_17.cc new file mode 100644 index 0000000000..963d6c77b2 --- /dev/null +++ b/tests/particles/particle_handler_17.cc @@ -0,0 +1,141 @@ +// --------------------------------------------------------------------- +// +// 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. +// +// --------------------------------------------------------------------- + +// Test insert_global_particles twice. Make sure we don't lose particles +// along the way, and that the global ids don't overlap. Test that we can +// set ids arbitrarily. + +#include +#include +#include + +#include + +#include + +#include +#include +#include +#include + +#include + +#include + +#include + +#include "../tests.h" + +template +void +insert( + Particles::ParticleHandler & particle_handler, + unsigned int n_points, + unsigned starting_id, + const std::vector>> &global_bounding_boxes) +{ + const unsigned int my_cpu = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD); + const unsigned int n_cpus = Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD); + + std::vector> particles(n_points); + + std::vector> points(n_points); + std::vector> properties(n_points, + {my_cpu + 1000.0, + my_cpu + 1000.0}); + std::vector ids(n_points); + + for (unsigned int i = 0; i < n_points; ++i) + { + auto id = starting_id + my_cpu * n_points + i; + particles[i].set_property_pool(particle_handler.get_property_pool()); + particles[i].set_location(random_point()); + particles[i].set_id(id); + particles[i].get_properties()[0] = my_cpu + 1000.0; + particles[i].get_properties()[1] = id; + } + + auto cpu_to_index = + particle_handler.insert_global_particles(particles, global_bounding_boxes); + + for (const auto &c : cpu_to_index) + { + deallog << "From cpu: " << c.first << " I got : "; + c.second.print(deallog); + } + + if (cpu_to_index.find(my_cpu) != cpu_to_index.end()) + cpu_to_index.erase(cpu_to_index.find(my_cpu)); + auto received = Utilities::MPI::some_to_some(MPI_COMM_WORLD, cpu_to_index); + + for (const auto &c : received) + { + deallog << "To cpu : " << c.first << " I sent : "; + c.second.print(deallog); + } +} + + +template +void +test() +{ + parallel::distributed::Triangulation tr(MPI_COMM_WORLD); + + GridGenerator::hyper_cube(tr); + tr.refine_global(2); + + MappingQ mapping(1); + + Particles::ParticleHandler particle_handler(tr, mapping, 2); + + const unsigned int n_points = 3; + const unsigned int my_cpu = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD); + const unsigned int n_cpus = Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD); + Testing::srand(my_cpu + 1); + + + // Distribute the local points to the processor that owns them + // on the triangulation + auto my_bounding_box = GridTools::compute_mesh_predicate_bounding_box( + tr, IteratorFilters::LocallyOwnedCell()); + + auto global_bounding_boxes = + Utilities::MPI::all_gather(MPI_COMM_WORLD, my_bounding_box); + + insert(particle_handler, n_points, 100, global_bounding_boxes); + insert(particle_handler, n_points, 200, global_bounding_boxes); + + for (auto p : particle_handler) + { + deallog << "Particle : " << p.get_id() << ", properties: " + << static_cast(p.get_properties()[0]) << " - " + << static_cast(p.get_properties()[1]) << 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(); +} diff --git a/tests/particles/particle_handler_17.with_p4est=true.mpirun=2.output b/tests/particles/particle_handler_17.with_p4est=true.mpirun=2.output new file mode 100644 index 0000000000..79ac5649a5 --- /dev/null +++ b/tests/particles/particle_handler_17.with_p4est=true.mpirun=2.output @@ -0,0 +1,26 @@ + +DEAL:0:2d/2d::From cpu: 0 I got : {0, 2} +DEAL:0:2d/2d::From cpu: 1 I got : {[1,2]} +DEAL:0:2d/2d::To cpu : 1 I sent : {1} +DEAL:0:2d/2d::From cpu: 1 I got : {0, 2} +DEAL:0:2d/2d::To cpu : 1 I sent : {[0,2]} +DEAL:0:2d/2d::Particle : 104, properties: 1001 - 104 +DEAL:0:2d/2d::Particle : 105, properties: 1001 - 105 +DEAL:0:2d/2d::Particle : 203, properties: 1001 - 203 +DEAL:0:2d/2d::Particle : 102, properties: 1000 - 102 +DEAL:0:2d/2d::Particle : 100, properties: 1000 - 100 +DEAL:0:2d/2d::Particle : 205, properties: 1001 - 205 + +DEAL:1:2d/2d::From cpu: 0 I got : {1} +DEAL:1:2d/2d::From cpu: 1 I got : {0} +DEAL:1:2d/2d::To cpu : 0 I sent : {[1,2]} +DEAL:1:2d/2d::From cpu: 0 I got : {[0,2]} +DEAL:1:2d/2d::From cpu: 1 I got : {1} +DEAL:1:2d/2d::To cpu : 0 I sent : {0, 2} +DEAL:1:2d/2d::Particle : 201, properties: 1000 - 201 +DEAL:1:2d/2d::Particle : 202, properties: 1000 - 202 +DEAL:1:2d/2d::Particle : 200, properties: 1000 - 200 +DEAL:1:2d/2d::Particle : 204, properties: 1001 - 204 +DEAL:1:2d/2d::Particle : 103, properties: 1001 - 103 +DEAL:1:2d/2d::Particle : 101, properties: 1000 - 101 +