]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Implement regular particle generation
authorRene Gassmoeller <rene.gassmoeller@mailbox.org>
Wed, 14 Aug 2019 20:44:34 +0000 (13:44 -0700)
committerRene Gassmoeller <rene.gassmoeller@mailbox.org>
Wed, 14 Aug 2019 21:49:34 +0000 (14:49 -0700)
include/deal.II/particles/particle_generator.h [new file with mode: 0644]
source/particles/CMakeLists.txt
source/particles/particle_generator.cc [new file with mode: 0644]
source/particles/particle_generator.inst.in [new file with mode: 0644]
tests/particles/particle_generator_01.cc [new file with mode: 0644]
tests/particles/particle_generator_01.with_p4est=true.output [new file with mode: 0644]

diff --git a/include/deal.II/particles/particle_generator.h b/include/deal.II/particles/particle_generator.h
new file mode 100644 (file)
index 0000000..e4cab73
--- /dev/null
@@ -0,0 +1,52 @@
+// ---------------------------------------------------------------------
+//
+// 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.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii_particles_particle_generator_h
+#define dealii_particles_particle_generator_h
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/fe/mapping.h>
+#include <deal.II/fe/mapping_q1.h>
+
+#include <deal.II/particles/particle_handler.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace Particles
+{
+  /**
+   * A namespace that contains all classes that are related to the particle
+   * generation.
+   */
+  namespace Generator
+  {
+    /**
+     * A function that generates particles in every cell at specified @p particle_reference_locations.
+     */
+    template <int dim, int spacedim = dim>
+    void
+    regular_reference_locations(
+      const parallel::distributed::Triangulation<dim, spacedim> &triangulation,
+      const std::vector<Point<dim>> & particle_reference_locations,
+      ParticleHandler<dim, spacedim> &particle_handler,
+      const Mapping<dim, spacedim> &  mapping = MappingQ1<dim, spacedim>());
+
+  } // namespace Generator
+} // namespace Particles
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
index d911d3ef1cd84ab6c6c654aa71e56ef1ab653121..0bbe6baaba855f70db6d2cbeea503b8a9b681be4 100644 (file)
@@ -21,6 +21,7 @@ SET(_src
   particle_accessor.cc
   particle_iterator.cc
   particle_handler.cc
+  particle_generator.cc
   property_pool.cc
   )
 
@@ -29,6 +30,7 @@ SET(_inst
   particle_accessor.inst.in
   particle_iterator.inst.in
   particle_handler.inst.in
+  particle_generator.inst.in
   )
 
 FILE(GLOB _header
diff --git a/source/particles/particle_generator.cc b/source/particles/particle_generator.cc
new file mode 100644 (file)
index 0000000..de923d0
--- /dev/null
@@ -0,0 +1,91 @@
+// ---------------------------------------------------------------------
+//
+// 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.
+//
+// ---------------------------------------------------------------------
+
+#include <deal.II/base/signaling_nan.h>
+
+#include <deal.II/particles/particle_generator.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace Particles
+{
+  namespace Generator
+  {
+    template <int dim, int spacedim = dim>
+    void
+    regular_reference_locations(
+      const parallel::distributed::Triangulation<dim, spacedim> &triangulation,
+      const std::vector<Point<dim>> & particle_reference_locations,
+      ParticleHandler<dim, spacedim> &particle_handler,
+      const Mapping<dim, spacedim> &  mapping)
+    {
+      types::particle_index n_particles_to_generate =
+        triangulation.n_locally_owned_active_cells() *
+        particle_reference_locations.size();
+
+      // The local start index is the number of all particles
+      // generated on lower MPI ranks.
+      types::particle_index local_start_index = 0;
+      MPI_Scan(&n_particles_to_generate,
+               &local_start_index,
+               1,
+               DEAL_II_PARTICLE_INDEX_MPI_TYPE,
+               MPI_SUM,
+               triangulation.get_communicator());
+      local_start_index -= n_particles_to_generate;
+
+      types::particle_index particle_index = local_start_index;
+
+      typename Triangulation<dim, spacedim>::active_cell_iterator
+        cell = triangulation.begin_active(),
+        endc = triangulation.end();
+
+      for (; cell != endc; cell++)
+        {
+          if (cell->is_locally_owned())
+            {
+              for (typename std::vector<Point<dim>>::const_iterator
+                     itr_particles_in_unit_cell =
+                       particle_reference_locations.begin();
+                   itr_particles_in_unit_cell !=
+                   particle_reference_locations.end();
+                   itr_particles_in_unit_cell++)
+                {
+                  const Point<spacedim> position_real =
+                    mapping.transform_unit_to_real_cell(
+                      cell, *itr_particles_in_unit_cell);
+
+                  const Particle<dim, spacedim> particle(
+                    position_real, *itr_particles_in_unit_cell, particle_index);
+                  particle_handler.insert_particle(particle, cell);
+                  ++particle_index;
+                }
+            }
+        }
+
+      particle_handler.update_cached_numbers();
+    }
+  } // namespace Generator
+} // namespace Particles
+
+
+
+DEAL_II_NAMESPACE_CLOSE
+
+DEAL_II_NAMESPACE_OPEN
+
+#include "particle_generator.inst"
+
+DEAL_II_NAMESPACE_CLOSE
diff --git a/source/particles/particle_generator.inst.in b/source/particles/particle_generator.inst.in
new file mode 100644 (file)
index 0000000..c642d28
--- /dev/null
@@ -0,0 +1,37 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 - 2018 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.
+//
+// ---------------------------------------------------------------------
+
+
+for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
+  {
+#if deal_II_dimension <= deal_II_space_dimension
+    namespace Particles
+    \{
+      namespace Generator
+      \{
+        template void
+        regular_reference_locations<deal_II_dimension, deal_II_space_dimension>(
+          const parallel::distributed::Triangulation<deal_II_dimension,
+                                                     deal_II_space_dimension>
+            &triangulation,
+          const std::vector<Point<deal_II_dimension>>
+            &particle_reference_locations,
+          ParticleHandler<deal_II_dimension, deal_II_space_dimension>
+            &particle_handler,
+          const Mapping<deal_II_dimension, deal_II_space_dimension> &mapping);
+      \}
+    \}
+#endif
+  }
diff --git a/tests/particles/particle_generator_01.cc b/tests/particles/particle_generator_01.cc
new file mode 100644 (file)
index 0000000..4e1265c
--- /dev/null
@@ -0,0 +1,82 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 - 2018 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// check the creation and destruction of particle within the particle handler
+// class using a particle generator.
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/fe/mapping_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_tools.h>
+
+#include <deal.II/particles/particle_generator.h>
+#include <deal.II/particles/particle_handler.h>
+
+#include "../tests.h"
+
+template <int dim, int spacedim>
+void
+test()
+{
+  {
+    parallel::distributed::Triangulation<dim, spacedim> tr(MPI_COMM_WORLD);
+
+    GridGenerator::hyper_cube(tr);
+    MappingQ<dim, spacedim> mapping(1);
+
+    Particles::ParticleHandler<dim, spacedim> particle_handler(tr, mapping);
+
+    std::vector<Point<dim>> particle_reference_locations(1, Point<dim>());
+    Particles::Generator::regular_reference_locations(
+      tr, particle_reference_locations, particle_handler);
+
+    deallog << "Particle number: " << particle_handler.n_global_particles()
+            << std::endl;
+
+    for (const auto &particle : particle_handler)
+      {
+        deallog << "Particle location: " << particle.get_location()
+                << std::endl;
+        deallog << "Particle reference location: "
+                << particle.get_reference_location() << std::endl;
+      }
+  }
+
+  deallog << "OK" << std::endl;
+}
+
+
+
+int
+main(int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+
+  initlog();
+
+  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_generator_01.with_p4est=true.output b/tests/particles/particle_generator_01.with_p4est=true.output
new file mode 100644 (file)
index 0000000..89addf5
--- /dev/null
@@ -0,0 +1,13 @@
+
+DEAL:2d/2d::Particle number: 1
+DEAL:2d/2d::Particle location: 0.00000 0.00000
+DEAL:2d/2d::Particle reference location: 0.00000 0.00000
+DEAL:2d/2d::OK
+DEAL:2d/3d::Particle number: 1
+DEAL:2d/3d::Particle location: 0.00000 0.00000 0.00000
+DEAL:2d/3d::Particle reference location: 0.00000 0.00000
+DEAL:2d/3d::OK
+DEAL:3d/3d::Particle number: 1
+DEAL:3d/3d::Particle location: 0.00000 0.00000 0.00000
+DEAL:3d/3d::Particle reference location: 0.00000 0.00000 0.00000
+DEAL:3d/3d::OK

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.