]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Make sure not to attach a PropertyPool if not asked for 6514/head
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Mon, 30 Apr 2018 00:07:17 +0000 (02:07 +0200)
committerMatthias Maier <tamiko@43-1.org>
Mon, 7 May 2018 20:24:48 +0000 (15:24 -0500)
include/deal.II/particles/particle.h
source/particles/particle.cc
source/particles/particle_handler.cc
tests/particles/particle_04.cc
tests/particles/particle_06.cc [new file with mode: 0644]
tests/particles/particle_06.output [new file with mode: 0644]

index 84bfd5793064bb9b10dc0b305dbc2f6f81538da9..db444913caee08b3554271dfe0c6ca5ec1b2da6f 100644 (file)
@@ -148,7 +148,7 @@ namespace Particles
      * by @p property_pool.
      */
     Particle (const void *&begin_data,
-              PropertyPool &property_pool);
+              PropertyPool *const = nullptr);
 
     /**
      * Move constructor for Particle, creates a particle from an existing
index 587a2573afe68461804b3ed184cedbffc665539f..600e90498224429746dac64ba87cee65fa67b1c2 100644 (file)
@@ -69,7 +69,7 @@ namespace Particles
 
   template <int dim, int spacedim>
   Particle<dim,spacedim>::Particle (const void *&data,
-                                    PropertyPool &new_property_pool)
+                                    PropertyPool *const new_property_pool)
   {
     const types::particle_index *id_data = static_cast<const types::particle_index *> (data);
     id = *id_data++;
@@ -81,14 +81,18 @@ namespace Particles
     for (unsigned int i = 0; i < dim; ++i)
       reference_location(i) = *pdata++;
 
-    property_pool = &new_property_pool;
-    properties = property_pool->allocate_properties_array();
+    property_pool = new_property_pool;
+    if (property_pool != nullptr)
+      properties = property_pool->allocate_properties_array();
+    else
+      properties = PropertyPool::invalid_handle;
 
     // See if there are properties to load
     if (has_properties())
       {
         const ArrayView<double> particle_properties = property_pool->get_properties(properties);
-        for (unsigned int i = 0; i < particle_properties.size(); ++i)
+        const unsigned int size = particle_properties.size();
+        for (unsigned int i = 0; i < size; ++i)
           particle_properties[i] = *pdata++;
       }
 
@@ -159,7 +163,7 @@ namespace Particles
   template <int dim, int spacedim>
   Particle<dim,spacedim>::~Particle ()
   {
-    if (properties != PropertyPool::invalid_handle)
+    if (property_pool != nullptr && properties != PropertyPool::invalid_handle)
       property_pool->deallocate_properties_array(properties);
   }
 
index e786eef43d9a3be5f582c407b0b4d34450ff3b14..6d365878132e27899fd72d6b05fdf2d36380faa7 100644 (file)
@@ -849,7 +849,7 @@ namespace Particles
 
         typename std::multimap<internal::LevelInd,Particle <dim,spacedim> >::iterator recv_particle =
           received_particles.insert(std::make_pair(internal::LevelInd(cell->level(),cell->index()),
-                                                   Particle<dim,spacedim>(recv_data_it,*property_pool)));
+                                                   Particle<dim,spacedim>(recv_data_it,property_pool.get())));
 
         if (load_callback)
           recv_data_it = load_callback(particle_iterator(received_particles,recv_particle),
@@ -1079,11 +1079,11 @@ namespace Particles
 #ifdef DEAL_II_WITH_CXX14
             position_hint = particles.emplace_hint(position_hint,
                                                    std::make_pair(cell->level(),cell->index()),
-                                                   Particle<dim,spacedim>(pdata,*property_pool));
+                                                   Particle<dim,spacedim>(pdata,property_pool.get()));
 #else
             position_hint = particles.insert(position_hint,
                                              std::make_pair(std::make_pair(cell->level(),cell->index()),
-                                                            Particle<dim,spacedim>(pdata,*property_pool)));
+                                                            Particle<dim,spacedim>(pdata,property_pool.get())));
 #endif
             ++position_hint;
           }
@@ -1101,11 +1101,11 @@ namespace Particles
 #ifdef DEAL_II_WITH_CXX14
             position_hint = particles.emplace_hint(position_hint,
                                                    std::make_pair(cell->level(),cell->index()),
-                                                   Particle<dim,spacedim>(pdata,*property_pool));
+                                                   Particle<dim,spacedim>(pdata,property_pool.get()));
 #else
             position_hint = particles.insert(position_hint,
                                              std::make_pair(std::make_pair(cell->level(),cell->index()),
-                                                            Particle<dim,spacedim>(pdata,*property_pool)));
+                                                            Particle<dim,spacedim>(pdata,property_pool.get())));
 #endif
             const Point<dim> p_unit = mapping->transform_real_to_unit_cell(cell, position_hint->second.get_location());
             position_hint->second.set_reference_location(p_unit);
@@ -1123,7 +1123,7 @@ namespace Particles
 
         for (unsigned int i = 0; i < *n_particles_in_cell_ptr; ++i)
           {
-            Particle<dim,spacedim> p (pdata,*property_pool);
+            Particle<dim,spacedim> p (pdata,property_pool.get());
 
             for (unsigned int child_index = 0; child_index < GeometryInfo<dim>::max_children_per_cell; ++child_index)
               {
index 49a6b6b0e070241a81db8ef50c5ba23f1251055d..182cecd4f02903e2ed115046081cd5fb8b76271c 100644 (file)
@@ -48,6 +48,7 @@ void test ()
     deallog << "Particle location: " << particle.get_location() << std::endl
             << "Particle reference location: " << particle.get_reference_location() << std::endl
             << "Particle index: " << particle.get_id() << std::endl;
+    Assert (!particle.has_properties(), ExcInternalError());
 
     std::vector<char> data(particle.serialized_size_in_bytes());
     void *write_pointer = static_cast<void *> (&data.front());
@@ -55,12 +56,12 @@ void test ()
     particle.write_data(write_pointer);
 
     const void *read_pointer = static_cast<const void *> (&data.front());
-    Particles::PropertyPool pool(0);
-    const Particles::Particle<dim,spacedim> new_particle(read_pointer,pool);
+    const Particles::Particle<dim,spacedim> new_particle(read_pointer);
 
     deallog << "Copy particle location: " << new_particle.get_location() << std::endl
             << "Copy particle reference location: " << new_particle.get_reference_location() << std::endl
             << "Copy particle index: " << new_particle.get_id() << std::endl;
+    Assert (!new_particle.has_properties(), ExcInternalError());
   }
 
   deallog << "OK" << std::endl;
diff --git a/tests/particles/particle_06.cc b/tests/particles/particle_06.cc
new file mode 100644 (file)
index 0000000..ad70ad2
--- /dev/null
@@ -0,0 +1,96 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Like particle_04, but also attach properties.
+
+#include "../tests.h"
+#include <deal.II/particles/particle.h>
+#include <deal.II/base/array_view.h>
+
+
+template <int dim, int spacedim>
+void test ()
+{
+  {
+    const unsigned int n_properties_per_particle = 3;
+    Particles::PropertyPool pool(n_properties_per_particle);
+
+    Point<spacedim> position;
+
+    position(0) = 0.3;
+    if (spacedim > 1)
+      position(1) = 0.5;
+    if (spacedim > 2)
+      position(2) = 0.7;
+
+    Point<dim> reference_position;
+    reference_position(0) = 0.2;
+    if (dim > 1)
+      reference_position(1) = 0.4;
+    if (dim > 2)
+      reference_position(2) = 0.6;
+
+    const types::particle_index index(7);
+
+    std::vector<double> properties = {0.15,0.45,0.75};
+
+    Particles::Particle<dim,spacedim> particle(position,reference_position,index);
+    particle.set_property_pool(pool);
+    particle.set_properties(ArrayView<double>(&properties[0],properties.size()));
+
+    deallog << "Particle location: " << particle.get_location() << std::endl
+            << "Particle reference location: " << particle.get_reference_location() << std::endl
+            << "Particle index: " << particle.get_id() << std::endl
+            << "Particle properties: " << std::vector<double>(particle.get_properties().begin(),
+                                                              particle.get_properties().end())
+            << std::endl;
+
+    std::vector<char> data(particle.serialized_size_in_bytes());
+    void *write_pointer = static_cast<void *> (&data.front());
+
+    particle.write_data(write_pointer);
+
+    const void *read_pointer = static_cast<const void *> (&data.front());
+    const Particles::Particle<dim,spacedim> new_particle(read_pointer, &pool);
+
+    deallog << "Copy particle location: " << new_particle.get_location() << std::endl
+            << "Copy particle reference location: " << new_particle.get_reference_location() << std::endl
+            << "Copy particle index: " << new_particle.get_id() << std::endl
+            << "Copy particle properties: " << std::vector<double>(new_particle.get_properties().begin(),
+                new_particle.get_properties().end())
+            << std::endl;
+  }
+
+  deallog << "OK" << std::endl;
+}
+
+
+
+int main ()
+{
+  initlog();
+
+  test<1,1>();
+  test<1,2>();
+  test<1,3>();
+
+  test<2,2>();
+  test<2,3>();
+
+  test<3,3>();
+
+}
diff --git a/tests/particles/particle_06.output b/tests/particles/particle_06.output
new file mode 100644 (file)
index 0000000..37135be
--- /dev/null
@@ -0,0 +1,55 @@
+
+DEAL::Particle location: 0.300000
+DEAL::Particle reference location: 0.200000
+DEAL::Particle index: 7
+DEAL::Particle properties: 0.150000 0.450000 0.750000
+DEAL::Copy particle location: 0.300000
+DEAL::Copy particle reference location: 0.200000
+DEAL::Copy particle index: 7
+DEAL::Copy particle properties: 0.150000 0.450000 0.750000
+DEAL::OK
+DEAL::Particle location: 0.300000 0.500000
+DEAL::Particle reference location: 0.200000
+DEAL::Particle index: 7
+DEAL::Particle properties: 0.150000 0.450000 0.750000
+DEAL::Copy particle location: 0.300000 0.500000
+DEAL::Copy particle reference location: 0.200000
+DEAL::Copy particle index: 7
+DEAL::Copy particle properties: 0.150000 0.450000 0.750000
+DEAL::OK
+DEAL::Particle location: 0.300000 0.500000 0.700000
+DEAL::Particle reference location: 0.200000
+DEAL::Particle index: 7
+DEAL::Particle properties: 0.150000 0.450000 0.750000
+DEAL::Copy particle location: 0.300000 0.500000 0.700000
+DEAL::Copy particle reference location: 0.200000
+DEAL::Copy particle index: 7
+DEAL::Copy particle properties: 0.150000 0.450000 0.750000
+DEAL::OK
+DEAL::Particle location: 0.300000 0.500000
+DEAL::Particle reference location: 0.200000 0.400000
+DEAL::Particle index: 7
+DEAL::Particle properties: 0.150000 0.450000 0.750000
+DEAL::Copy particle location: 0.300000 0.500000
+DEAL::Copy particle reference location: 0.200000 0.400000
+DEAL::Copy particle index: 7
+DEAL::Copy particle properties: 0.150000 0.450000 0.750000
+DEAL::OK
+DEAL::Particle location: 0.300000 0.500000 0.700000
+DEAL::Particle reference location: 0.200000 0.400000
+DEAL::Particle index: 7
+DEAL::Particle properties: 0.150000 0.450000 0.750000
+DEAL::Copy particle location: 0.300000 0.500000 0.700000
+DEAL::Copy particle reference location: 0.200000 0.400000
+DEAL::Copy particle index: 7
+DEAL::Copy particle properties: 0.150000 0.450000 0.750000
+DEAL::OK
+DEAL::Particle location: 0.300000 0.500000 0.700000
+DEAL::Particle reference location: 0.200000 0.400000 0.600000
+DEAL::Particle index: 7
+DEAL::Particle properties: 0.150000 0.450000 0.750000
+DEAL::Copy particle location: 0.300000 0.500000 0.700000
+DEAL::Copy particle reference location: 0.200000 0.400000 0.600000
+DEAL::Copy particle index: 7
+DEAL::Copy particle properties: 0.150000 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.