]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Address comments, add test for ParticleIterator and ParticleAccessor 4729/head
authorRene Gassmoeller <rene.gassmoeller@mailbox.org>
Mon, 25 Sep 2017 17:35:32 +0000 (11:35 -0600)
committerRene Gassmoeller <rene.gassmoeller@mailbox.org>
Thu, 28 Sep 2017 08:14:36 +0000 (10:14 +0200)
include/deal.II/particles/particle_accessor.h
source/particles/particle_accessor.cc
tests/particles/particle_iterator_01.cc [new file with mode: 0644]
tests/particles/particle_iterator_01.output [new file with mode: 0644]

index 0336063e5e293c40c46714372760f1a16c973b61..bdb7aff9d10c934920feb813d40db2acb723076f 100644 (file)
@@ -18,7 +18,7 @@
 
 #include <deal.II/particles/particle.h>
 #include <deal.II/base/array_view.h>
-#include <deal.II/distributed/tria.h>
+#include <deal.II/grid/tria.h>
 
 DEAL_II_NAMESPACE_OPEN
 
@@ -141,8 +141,8 @@ namespace Particles
      * but the triangulation itself is not stored in the particle this
      * operation requires a reference to the triangulation.
      */
-    typename parallel::distributed::Triangulation<dim,spacedim>::cell_iterator
-    get_surrounding_cell (const parallel::distributed::Triangulation<dim,spacedim> &triangulation) const;
+    typename Triangulation<dim,spacedim>::cell_iterator
+    get_surrounding_cell (const Triangulation<dim,spacedim> &triangulation) const;
 
     /**
      * Serialize the contents of this class.
@@ -189,19 +189,18 @@ namespace Particles
      * A pointer to the container that stores the particles. Obviously,
      * this accessor is invalidated if the container changes.
      */
-    std::multimap<types::LevelInd, Particle<dim,spacedim> >            *map;
+    std::multimap<types::LevelInd, Particle<dim,spacedim> > *map;
 
     /**
      * An iterator into the container of particles. Obviously,
-       * this accessor is invalidated if the container changes.
+     * this accessor is invalidated if the container changes.
      */
-    typename std::multimap<types::LevelInd, Particle<dim,spacedim> >::iterator  particle;
+    typename std::multimap<types::LevelInd, Particle<dim,spacedim> >::iterator particle;
 
     /**
      * Make ParticleIterator a friend to allow it constructing ParticleAccessors.
      */
     template <int, int> friend class ParticleIterator;
-    template <int, int> friend class ParticleHandler;
   };
 }
 
index ea82fedfd16441f2fec827db3e2a6f08abefe451..42fabc146982ff20062785c992f655f6c0ee0410 100644 (file)
@@ -148,13 +148,13 @@ namespace Particles
 
 
   template <int dim, int spacedim>
-  typename parallel::distributed::Triangulation<dim,spacedim>::cell_iterator
-  ParticleAccessor<dim,spacedim>::get_surrounding_cell (const parallel::distributed::Triangulation<dim,spacedim> &triangulation) const
+  typename Triangulation<dim,spacedim>::cell_iterator
+  ParticleAccessor<dim,spacedim>::get_surrounding_cell (const Triangulation<dim,spacedim> &triangulation) const
   {
     Assert(particle != map->end(),
            ExcInternalError());
 
-    const typename parallel::distributed::Triangulation<dim,spacedim>::cell_iterator cell (&triangulation,
+    const typename Triangulation<dim,spacedim>::cell_iterator cell (&triangulation,
         particle->first.first,
         particle->first.second);
     return cell;
diff --git a/tests/particles/particle_iterator_01.cc b/tests/particles/particle_iterator_01.cc
new file mode 100644 (file)
index 0000000..67ccd45
--- /dev/null
@@ -0,0 +1,86 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Like particle_03, but tests the creation and use of a
+// particle iterator from the created particle.
+
+#include "../tests.h"
+#include <deal.II/particles/particle.h>
+#include <deal.II/particles/particle_iterator.h>
+
+#include <deal.II/base/array_view.h>
+
+
+template <int dim>
+void test ()
+{
+  {
+    const unsigned int n_properties_per_particle = 3;
+    Particles::PropertyPool pool(n_properties_per_particle);
+
+    Point<dim> position;
+    position(0) = 0.3;
+
+    if (dim>1)
+      position(1) = 0.5;
+
+    Point<dim> reference_position;
+    reference_position(0) = 0.2;
+
+    if (dim>1)
+      reference_position(1) = 0.4;
+
+    const Particles::types::particle_index index(7);
+
+    std::vector<double> properties = {0.15,0.45,0.75};
+
+    Particles::Particle<dim> particle(position,reference_position,index);
+    particle.set_property_pool(pool);
+    particle.set_properties(ArrayView<double>(&properties[0],properties.size()));
+
+    std::multimap<Particles::types::LevelInd,Particles::Particle<dim> > map;
+
+    Particles::types::LevelInd level_index = std::make_pair(0,0);
+    map.insert(std::make_pair(level_index,particle));
+
+    particle.get_properties()[0] = 0.05;
+    map.insert(std::make_pair(level_index,particle));
+
+    Particles::ParticleIterator<dim> particle_it (map,map.begin());
+    Particles::ParticleIterator<dim> particle_end (map,map.end());
+
+    for (; particle_it != particle_end; ++particle_it)
+      {
+        deallog << "Particle position: "
+                << (*particle_it).get_location()
+                << std::endl
+                << "Particle properties: "
+                << std::vector<double>(particle_it->get_properties().begin(),particle_it->get_properties().end())
+                << std::endl;
+      }
+  }
+
+  deallog << "OK" << std::endl;
+}
+
+
+
+int main ()
+{
+  initlog();
+  test<2>();
+}
diff --git a/tests/particles/particle_iterator_01.output b/tests/particles/particle_iterator_01.output
new file mode 100644 (file)
index 0000000..777a4e2
--- /dev/null
@@ -0,0 +1,6 @@
+
+DEAL::Particle position: 0.300000 0.500000
+DEAL::Particle properties: 0.150000 0.450000 0.750000
+DEAL::Particle position: 0.300000 0.500000
+DEAL::Particle properties: 0.0500000 0.450000 0.750000
+DEAL::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.