]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Modify usage and include in tests
authorRene Gassmoeller <rene.gassmoeller@mailbox.org>
Tue, 20 Feb 2024 08:00:15 +0000 (03:00 -0500)
committerRene Gassmoeller <rene.gassmoeller@mailbox.org>
Tue, 20 Feb 2024 18:32:26 +0000 (13:32 -0500)
examples/step-68/step-68.cc
include/deal.II/particles/particle.h
include/deal.II/particles/particle_accessor.h
include/deal.II/particles/particle_handler.h
include/deal.II/particles/property_pool.h
source/particles/particle_handler.cc
tests/particles/step-68.cc
tests/performance/timing_step_68.cc
tests/simplex/step-68.cc

index 149892bc5fe304c3330dab15dc0334cc95f2f930..12214e966149c308b7669a4be6e349b75b9118a3 100644 (file)
@@ -466,7 +466,7 @@ namespace Step68
       {
         // We calculate the velocity of the particles using their current
         // location.
-        Point<dim> particle_location = particle.get_location();
+        Point<dim> &particle_location = particle.get_location();
         velocity.vector_value(particle_location, particle_velocity);
 
         // This updates the position of the particles and sets the old position
@@ -474,8 +474,6 @@ namespace Step68
         for (int d = 0; d < dim; ++d)
           particle_location[d] += particle_velocity[d] * dt;
 
-        particle.set_location(particle_location);
-
         // We store the processor id (a scalar) and the particle velocity (a
         // vector) in the particle properties. In this example, this is done
         // purely for visualization purposes.
@@ -534,11 +532,10 @@ namespace Step68
         for (unsigned int particle_index = 0; particle != pic.end();
              ++particle, ++particle_index)
           {
-            Point<dim>            particle_location = particle->get_location();
+            Point<dim>            &particle_location = particle->get_location();
             const Tensor<1, dim> &particle_velocity =
               evaluator.get_value(particle_index);
             particle_location += particle_velocity * dt;
-            particle->set_location(particle_location);
 
             // Again, we store the particle velocity and the processor id in the
             // particle properties for visualization purposes.
index 600d6712935baba956c835af10791f57c099d730..1d60e894a4e6f56df39705b5e9f2a897710d2ffe 100644 (file)
@@ -287,12 +287,6 @@ namespace Particles
     const Point<dim> &
     get_reference_location() const;
 
-    /**
-     * Return the reference location of this particle in its current cell.
-     */
-    Point<dim> &
-    get_reference_location();
-
     /**
      * Return the ID number of this particle. The ID of a particle is intended
      * to be a property that is globally unique even in parallel computations
@@ -578,15 +572,6 @@ namespace Particles
 
 
 
-  template <int dim, int spacedim>
-  inline Point<dim> &
-  Particle<dim, spacedim>::get_reference_location()
-  {
-    return property_pool->get_reference_location(property_pool_handle);
-  }
-
-
-
   template <int dim, int spacedim>
   inline types::particle_index
   Particle<dim, spacedim>::get_id() const
index 1bce933fe6be46fdae56ba36823c851ac9b322c7..267d5c10cb68eee4f6efa7691170af0bf18fd92d 100644 (file)
@@ -149,7 +149,7 @@ namespace Particles
     set_location(const Point<spacedim> &new_location);
 
     /**
-     * Get the location of this particle.
+     * Get read-access to the location of this particle.
      *
      * @return The location of this particle.
      */
@@ -157,7 +157,23 @@ namespace Particles
     get_location() const;
 
     /**
-     * Get the location of this particle.
+     * Get read- and write-access to the location of this particle.
+     * Note that changing the location does not check
+     * whether this is a valid location in the simulation domain.
+     * 
+     * @note In parallel programs, the ParticleHandler class stores particles
+     *   on both the locally owned cells, as well as on ghost cells. The
+     *   particles on the latter are *copies* of particles owned on other
+     *   processors, and should therefore be treated in the same way as
+     *   ghost entries in
+     *   @ref GlossGhostedVector "vectors with ghost elements"
+     *   or
+     *   @ref GlossGhostCell "ghost cells":
+     *   In both cases, one should
+     *   treat the ghost elements or cells as `const` objects that shouldn't
+     *   be modified even if the objects allow for calls that modify
+     *   properties. Rather, properties should only be modified on processors
+     *   that actually *own* the particle.
      *
      * @return The location of this particle.
      */
@@ -193,12 +209,6 @@ namespace Particles
     const Point<dim> &
     get_reference_location() const;
 
-    /**
-     * Return the reference location of this particle in its current cell.
-     */
-    Point<dim> &
-    get_reference_location();
-
     /**
      * Set the ID number of this particle.
      */
@@ -698,17 +708,6 @@ namespace Particles
 
 
 
-  template <int dim, int spacedim>
-  inline Point<dim> &
-  ParticleAccessor<dim, spacedim>::get_reference_location()
-  {
-    Assert(state() == IteratorState::valid, ExcInternalError());
-
-    return property_pool->get_reference_location(get_handle());
-  }
-
-
-
   template <int dim, int spacedim>
   inline types::particle_index
   ParticleAccessor<dim, spacedim>::get_id() const
index 4a5c98d20d47c2d8c164333b6be6e89382f24c38..672a116c37c1b3d96c23e98ae877e5656f9d971a 100644 (file)
@@ -1341,12 +1341,15 @@ namespace Particles
                     get_next_free_particle_index() * spacedim);
     for (auto &p : *this)
       {
-        Point<spacedim> new_point(displace_particles ? p.get_location() :
-                                                       Point<spacedim>());
+        Point<spacedim> &position = p.get_location();
         const auto      id = p.get_id();
-        for (unsigned int i = 0; i < spacedim; ++i)
-          new_point[i] += input_vector[id * spacedim + i];
-        p.set_location(new_point);
+
+        if (displace_particles)
+          for (unsigned int i = 0; i < spacedim; ++i)
+            position[i] += input_vector[id * spacedim + i];
+        else
+          for (unsigned int i = 0; i < spacedim; ++i)
+            position[i] = input_vector[id * spacedim + i];
       }
     sort_particles_into_subdomains_and_cells();
   }
index 7131149956045df41c095d8bab0b4f4feba5e401..5e4643553d1a1445b9372b6e7be25c58ed69d8ab 100644 (file)
@@ -177,13 +177,6 @@ namespace Particles
     const Point<dim> &
     get_reference_location(const Handle handle) const;
 
-    /**
-     * Return a writeable reference to the reference location of a particle
-     * identified by the given `handle`.
-     */
-    Point<dim> &
-    get_reference_location(const Handle handle);
-
     /**
      * Set the reference location of a particle identified by the given
      * `handle`.
@@ -391,29 +384,6 @@ namespace Particles
 
 
 
-  template <int dim, int spacedim>
-  inline Point<dim> &
-  PropertyPool<dim, spacedim>::get_reference_location(const Handle handle)
-  {
-    const std::vector<double>::size_type data_index =
-      (handle != invalid_handle) ? handle : 0;
-
-    // Ideally we would need to assert that 'handle' has not been deallocated
-    // by searching through 'currently_available_handles'. However, this
-    // is expensive and this function is performance critical, so instead
-    // just check against the array range, and rely on the fact
-    // that handles are invalidated when handed over to
-    // deallocate_properties_array().
-    Assert(data_index <= reference_locations.size() - 1,
-           ExcMessage("Invalid location handle. This can happen if the "
-                      "handle was duplicated and then one copy was deallocated "
-                      "before trying to access the properties."));
-
-    return reference_locations[data_index];
-  }
-
-
-
   template <int dim, int spacedim>
   inline void
   PropertyPool<dim, spacedim>::set_reference_location(
index df6d147959ebdfc322d2ecbcd04c6f9e4154877d..0d47946dc6c8928dbd604d7a07400bcb1212c7ea 100644 (file)
@@ -1127,10 +1127,11 @@ namespace Particles
     unsigned int i = 0;
     for (auto it = begin(); it != end(); ++it, ++i)
       {
+        Point<spacedim> &location = it->get_location();
         if (displace_particles)
-          it->set_location(it->get_location() + new_positions[i]);
+          location += new_positions[i];
         else
-          it->set_location(new_positions[i]);
+          location = new_positions[i];
       }
     sort_particles_into_subdomains_and_cells();
   }
@@ -1150,7 +1151,7 @@ namespace Particles
     Vector<double> new_position(spacedim);
     for (auto &particle : *this)
       {
-        Point<spacedim> particle_location = particle.get_location();
+        Point<spacedim> &particle_location = particle.get_location();
         function.vector_value(particle_location, new_position);
         if (displace_particles)
           for (unsigned int d = 0; d < spacedim; ++d)
@@ -1158,7 +1159,6 @@ namespace Particles
         else
           for (unsigned int d = 0; d < spacedim; ++d)
             particle_location[d] = new_position[d];
-        particle.set_location(particle_location);
       }
     sort_particles_into_subdomains_and_cells();
   }
index 572ad607670ffa57962a2a80feb6303d19cfa87a..24c1c41bd7d12bb4ac724f3ead338d4ae52585b0 100644 (file)
@@ -440,7 +440,7 @@ namespace Step68
       {
         // We calculate the velocity of the particles using their current
         // location.
-        Point<dim> particle_location = particle.get_location();
+        Point<dim> &particle_location = particle.get_location();
         velocity.vector_value(particle_location, particle_velocity);
 
         // This updates the position of the particles and sets the old position
@@ -448,8 +448,6 @@ namespace Step68
         for (int d = 0; d < dim; ++d)
           particle_location[d] += particle_velocity[d] * dt;
 
-        particle.set_location(particle_location);
-
         // We store the processor id (a scalar) and the particle velocity (a
         // vector) in the particle properties. In this example, this is done
         // purely for visualization purposes.
@@ -513,10 +511,9 @@ namespace Step68
                   local_dof_values[j];
               }
 
-            Point<dim> particle_location = particle->get_location();
+            Point<dim> &particle_location = particle->get_location();
             for (int d = 0; d < dim; ++d)
               particle_location[d] += particle_velocity[d] * dt;
-            p.set_location(particle_location);
 
             // Again, we store the particle velocity and the processor id in the
             // particle properties for visualization purposes.
index 4150ab19f750112aada3567e2486113ba49ae90d..0caa5178419186d94d2670c60d19e59863144e94 100644 (file)
@@ -317,11 +317,10 @@ namespace Step68
         for (unsigned int particle_index = 0; particle != pic.end();
              ++particle, ++particle_index)
           {
-            Point<dim>            particle_location = particle->get_location();
+            Point<dim>            &particle_location = particle->get_location();
             const Tensor<1, dim> &particle_velocity =
               evaluator.get_value(particle_index);
             particle_location += particle_velocity * dt;
-            particle->set_location(particle_location);
 
             ArrayView<double> properties = particle->get_properties();
             for (int d = 0; d < dim; ++d)
index ee645b6d2679473d4d43c890e6ddb8ce6dfc1983..524aac86961f2e3563607851101b35e3bc9595b6 100644 (file)
@@ -377,10 +377,9 @@ namespace Step68
                   local_dof_values[j];
               }
 
-            Point<dim> particle_location = particle->get_location();
+            Point<dim> &particle_location = particle->get_location();
             for (int d = 0; d < dim; ++d)
               particle_location[d] += particle_velocity[d] * dt;
-            p.set_location(particle_location);
 
             // Again, we store the particle velocity and the processor id in the
             // particle properties for visualization purposes.

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.