]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Test also the version with Particles.
authorLuca Heltai <luca.heltai@sissa.it>
Thu, 11 Jun 2020 10:24:53 +0000 (12:24 +0200)
committerLuca Heltai <luca.heltai@sissa.it>
Thu, 11 Jun 2020 10:24:53 +0000 (12:24 +0200)
source/particles/particle_handler.cc
tests/particles/particle_handler_17.cc [new file with mode: 0644]
tests/particles/particle_handler_17.with_p4est=true.mpirun=2.output [new file with mode: 0644]

index ad1313d2287f4ce14f616ffc434930841bbc2b88..0f3f4274d371cac5c3a4ab731cb53240507414fa 100644 (file)
@@ -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<double>(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 (file)
index 0000000..963d6c7
--- /dev/null
@@ -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 <deal.II/base/mpi.h>
+#include <deal.II/base/std_cxx20/iota_view.h>
+#include <deal.II/base/utilities.h>
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/fe/mapping_q.h>
+
+#include <deal.II/grid/filtered_iterator.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/grid_tools_cache.h>
+
+#include <deal.II/particles/particle_handler.h>
+
+#include <unistd.h>
+
+#include <iostream>
+
+#include "../tests.h"
+
+template <int dim, int spacedim>
+void
+insert(
+  Particles::ParticleHandler<dim, spacedim> &            particle_handler,
+  unsigned int                                           n_points,
+  unsigned                                               starting_id,
+  const std::vector<std::vector<BoundingBox<spacedim>>> &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::Particle<dim, spacedim>> particles(n_points);
+
+  std::vector<Point<spacedim>>       points(n_points);
+  std::vector<std::vector<double>>   properties(n_points,
+                                              {my_cpu + 1000.0,
+                                               my_cpu + 1000.0});
+  std::vector<types::particle_index> 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<spacedim>());
+      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 <int dim, int spacedim>
+void
+test()
+{
+  parallel::distributed::Triangulation<dim, spacedim> tr(MPI_COMM_WORLD);
+
+  GridGenerator::hyper_cube(tr);
+  tr.refine_global(2);
+
+  MappingQ<dim, spacedim> mapping(1);
+
+  Particles::ParticleHandler<dim, spacedim> 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<unsigned int>(p.get_properties()[0]) << " - "
+              << static_cast<unsigned int>(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 (file)
index 0000000..79ac564
--- /dev/null
@@ -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
+

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.