]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Particles as indexable. 10536/head
authorLuca Heltai <luca.heltai@sissa.it>
Mon, 15 Jun 2020 16:55:45 +0000 (18:55 +0200)
committerLuca Heltai <luca.heltai@sissa.it>
Mon, 15 Jun 2020 16:55:45 +0000 (18:55 +0200)
doc/news/changes/minor/20200615LucaHeltai [new file with mode: 0644]
include/deal.II/particles/particle.h
include/deal.II/particles/particle_accessor.h
tests/particles/particle_09.cc [new file with mode: 0644]
tests/particles/particle_09.output [new file with mode: 0644]
tests/particles/particle_handler_serial_08.cc [new file with mode: 0644]
tests/particles/particle_handler_serial_08.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20200615LucaHeltai b/doc/news/changes/minor/20200615LucaHeltai
new file mode 100644 (file)
index 0000000..2bb3031
--- /dev/null
@@ -0,0 +1,3 @@
+New: Particles::Particle and Particles::ParticleAccessor can now be used as
+indexable in boost::rtree objects. 
+<br> (Luca Heltai, 2020/06/15)
index dc239e1caeabec5345bc21f73a4334a08b28f8b5..9208c9ef941810a8f5e52ff752aedbfc7eb4671d 100644 (file)
@@ -432,4 +432,39 @@ namespace Particles
 
 DEAL_II_NAMESPACE_CLOSE
 
+
+namespace boost
+{
+  namespace geometry
+  {
+    namespace index
+    {
+      // Forward declaration of bgi::indexable
+      template <class T>
+      struct indexable;
+
+      /**
+       * Make sure we can construct an RTree of Particles::Particle objects.
+       */
+      template <int dim, int spacedim>
+      struct indexable<dealii::Particles::Particle<dim, spacedim>>
+      {
+        /**
+         * boost::rtree expects a const reference to an indexable object. For
+         * a Particles::Particle object, this is its reference location.
+         */
+        using result_type = const dealii::Point<spacedim> &;
+
+        result_type
+        operator()(
+          const dealii::Particles::Particle<dim, spacedim> &particle) const
+        {
+          return particle.get_location();
+        }
+      };
+
+    } // namespace index
+  }   // namespace geometry
+} // namespace boost
+
 #endif
index 415a69cddc631fdfc77c5a32a34f14c89c4463ef..50a38480fc07b8a03bcfe71df1dc21e7f939e7a7 100644 (file)
@@ -252,4 +252,38 @@ namespace Particles
 
 DEAL_II_NAMESPACE_CLOSE
 
+namespace boost
+{
+  namespace geometry
+  {
+    namespace index
+    {
+      // Forward declaration of bgi::indexable
+      template <class T>
+      struct indexable;
+
+      /**
+       * Make sure we can construct an RTree from Particles::ParticleAccessor
+       * objects.
+       */
+      template <int dim, int spacedim>
+      struct indexable<dealii::Particles::ParticleAccessor<dim, spacedim>>
+      {
+        /**
+         * boost::rtree expects a const reference to an indexable object. For
+         * a Particles::Particle object, this is its reference location.
+         */
+        using result_type = const dealii::Point<spacedim> &;
+
+        result_type
+        operator()(const dealii::Particles::ParticleAccessor<dim, spacedim>
+                     &accessor) const
+        {
+          return accessor.get_location();
+        }
+      };
+    } // namespace index
+  }   // namespace geometry
+} // namespace boost
+
 #endif
diff --git a/tests/particles/particle_09.cc b/tests/particles/particle_09.cc
new file mode 100644 (file)
index 0000000..e7ebf02
--- /dev/null
@@ -0,0 +1,54 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 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 constructing an rtree of particles.
+
+#include <deal.II/numerics/rtree.h>
+
+#include <deal.II/particles/particle.h>
+
+#include "../tests.h"
+
+namespace bgi = boost::geometry::index;
+
+template <int dim, int spacedim>
+void
+test()
+{
+  const int                                       n_particles = 10;
+  std::vector<Particles::Particle<dim, spacedim>> particles(n_particles);
+
+  unsigned int id = 0;
+  for (auto &p : particles)
+    {
+      p.set_location(random_point<spacedim>());
+      p.set_id(id++);
+    }
+
+  auto tree = pack_rtree(particles);
+
+  auto p = random_point<spacedim>();
+  for (const auto part : tree | bgi::adaptors::queried(bgi::nearest(p, 3)))
+    deallog << "Particle " << part.get_id() << " is close to " << p
+            << " (location = " << part.get_location() << ")" << std::endl;
+}
+
+int
+main()
+{
+  initlog();
+  test<2, 2>();
+}
diff --git a/tests/particles/particle_09.output b/tests/particles/particle_09.output
new file mode 100644 (file)
index 0000000..b7fe350
--- /dev/null
@@ -0,0 +1,4 @@
+
+DEAL::Particle 6 is close to 0.0163006 0.242887 (location = 0.364784 0.513401)
+DEAL::Particle 4 is close to 0.0163006 0.242887 (location = 0.277775 0.553970)
+DEAL::Particle 9 is close to 0.0163006 0.242887 (location = 0.141603 0.606969)
diff --git a/tests/particles/particle_handler_serial_08.cc b/tests/particles/particle_handler_serial_08.cc
new file mode 100644 (file)
index 0000000..7b12e1c
--- /dev/null
@@ -0,0 +1,62 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 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 constructing an rtree of particles from a ParticleHandler object.
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/numerics/rtree.h>
+
+#include <deal.II/particles/particle.h>
+#include <deal.II/particles/particle_handler.h>
+
+#include "../tests.h"
+
+namespace bgi = boost::geometry::index;
+
+template <int dim, int spacedim>
+void
+test()
+{
+  Triangulation<dim, spacedim> tria;
+  GridGenerator::hyper_cube(tria);
+  tria.refine_global(2);
+
+  Particles::ParticleHandler<dim, spacedim> particle_handler(
+    tria, StaticMappingQ1<dim, spacedim>::mapping);
+
+  const int                    n_particles = 10;
+  std::vector<Point<spacedim>> particles(n_particles);
+
+  for (auto &p : particles)
+    p = random_point<spacedim>();
+
+  particle_handler.insert_particles(particles);
+
+  auto tree = pack_rtree(particle_handler.begin(), particle_handler.end());
+
+  auto p = random_point<spacedim>();
+  for (const auto &part : tree | bgi::adaptors::queried(bgi::nearest(p, 3)))
+    deallog << "Particle " << part.get_id() << " is close to " << p
+            << " (location = " << part.get_location() << ")" << std::endl;
+}
+
+int
+main()
+{
+  initlog();
+  test<2, 2>();
+}
diff --git a/tests/particles/particle_handler_serial_08.output b/tests/particles/particle_handler_serial_08.output
new file mode 100644 (file)
index 0000000..3f1d499
--- /dev/null
@@ -0,0 +1,4 @@
+
+DEAL::Particle 6 is close to 0.0163006 0.242887 (location = 0.364784 0.513401)
+DEAL::Particle 9 is close to 0.0163006 0.242887 (location = 0.141603 0.606969)
+DEAL::Particle 4 is close to 0.0163006 0.242887 (location = 0.277775 0.553970)

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.