]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Rename ParticleHandler::locally_relevant_ids() to locally_owned_particles() 12119/head
authorPeter Munch <peterrmuench@gmail.com>
Sat, 1 May 2021 10:03:36 +0000 (12:03 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Mon, 3 May 2021 05:47:53 +0000 (07:47 +0200)
16 files changed:
doc/news/changes/incompatibilities/20210501Munch [new file with mode: 0644]
examples/step-70/step-70.cc
include/deal.II/particles/particle_handler.h
include/deal.II/particles/utilities.h
source/particles/data_out.cc
source/particles/particle_handler.cc
source/particles/utilities.cc
tests/particles/extract_index_set_01.cc
tests/particles/interpolate_on_particle_01.cc
tests/particles/particle_interpolation_01.cc
tests/particles/particle_interpolation_02.cc
tests/particles/particle_interpolation_03.cc
tests/particles/particle_interpolation_04.cc
tests/particles/particle_interpolation_05.cc
tests/particles/set_particle_positions_01.cc
tests/particles/set_particle_positions_02.cc

diff --git a/doc/news/changes/incompatibilities/20210501Munch b/doc/news/changes/incompatibilities/20210501Munch
new file mode 100644 (file)
index 0000000..87e2636
--- /dev/null
@@ -0,0 +1,4 @@
+Renamed: The function ParticleHandler::locally_relevant_ids() has been deprecated.
+Please use the new function ParticleHandler::locally_owned_particle_ids() instead.
+<br>
+(Peter Munch, 2021/05/01)
index 8ddfb1f216228c28c2701f020c9c98aece3b7d73..d5a2aced9bc3cd9dc83bcde76aa718bbb95ec769 100644 (file)
@@ -1062,7 +1062,7 @@ namespace Step70
     // successive vector elements (this is what the IndexSet::tensor_priduct()
     // function does).
     locally_owned_tracer_particle_coordinates =
-      tracer_particle_handler.locally_relevant_ids().tensor_product(
+      tracer_particle_handler.locally_owned_particle_ids().tensor_product(
         complete_index_set(spacedim));
 
     // At the beginning of the simulation, all particles are in their original
@@ -1865,7 +1865,7 @@ namespace Step70
           tracer_particle_velocities *= time_step;
 
           locally_relevant_tracer_particle_coordinates =
-            tracer_particle_handler.locally_relevant_ids().tensor_product(
+            tracer_particle_handler.locally_owned_particle_ids().tensor_product(
               complete_index_set(spacedim));
 
           relevant_tracer_particle_displacements.reinit(
index 267dfe30af791cbd80ecdecc3f93018b757d268b..87f8c1e444cdc4ff1ffb6321db77fba48d951c85 100644 (file)
@@ -444,10 +444,10 @@ namespace Particles
      *
      * The vector @p input_vector should have read access to the indices
      * created by extracting the locally relevant ids with
-     * locally_relevant_ids(), and taking its tensor
+     * locally_owned_particle_ids(), and taking its tensor
      * product with the index set representing the range `[0, spacedim)`, i.e.:
      * @code
-     * IndexSet ids = particle_handler.locally_relevant_ids().
+     * IndexSet ids = particle_handler.locally_owned_particle_ids().
      *  tensor_product(complete_index_set(spacedim));
      * @endcode
      *
@@ -497,7 +497,7 @@ namespace Particles
      * get_particle_positions() function, and then modify the resulting vector.
      *
      * @param [in] new_positions A vector of points of dimension
-     * particle_handler.n_locally_owned_particles()
+     * particle_handler.n_locally_owned_particle_ids()
      *
      * @param [in] displace_particles When true, this function adds the value
      * of the vector of points to the
@@ -641,7 +641,7 @@ namespace Particles
      * triangulation.
      */
     types::particle_index
-    n_locally_owned_particles() const;
+    n_locally_owned_particle_ids() const;
 
     /**
      * Return the next free particle index in the global set
@@ -651,6 +651,27 @@ namespace Particles
     types::particle_index
     get_next_free_particle_index() const;
 
+    /**
+     * Extract an IndexSet with global dimensions equal to
+     * get_next_free_particle_index(), containing the locally owned
+     * particle indices.
+     *
+     * This function can be used to construct distributed vectors and matrices
+     * to manipulate particles using linear algebra operations.
+     *
+     * Notice that it is the user's responsibility to guarantee that particle
+     * indices are unique, and no check is performed to verify that this is the
+     * case, nor that the union of all IndexSet objects on each mpi process is
+     * complete.
+     *
+     * @return An IndexSet of size get_next_free_particle_index(), containing
+     * n_locally_owned_particle() indices.
+     *
+     * @deprecated Use locally_owned_particle_ids() instead.
+     */
+    DEAL_II_DEPRECATED IndexSet
+                       locally_relevant_ids() const;
+
     /**
      * Extract an IndexSet with global dimensions equal to
      * get_next_free_particle_index(), containing the locally owned
@@ -668,7 +689,7 @@ namespace Particles
      * n_locally_owned_particle() indices.
      */
     IndexSet
-    locally_relevant_ids() const;
+    locally_owned_particle_ids() const;
 
     /**
      * Return the number of properties each particle has.
@@ -1146,6 +1167,15 @@ namespace Particles
       output_vector.compress(VectorOperation::insert);
   }
 
+
+
+  template <int dim, int spacedim>
+  inline IndexSet
+  ParticleHandler<dim, spacedim>::locally_relevant_ids() const
+  {
+    return this->locally_owned_particle_ids();
+  }
+
 } // namespace Particles
 
 DEAL_II_NAMESPACE_CLOSE
index 3b60163226f69f9926009d288cabc6b887aeb478..b6721504a8a8347be9f1dccb728dcde5b8252521 100644 (file)
@@ -191,7 +191,7 @@ namespace Particles
       OutputVectorType &                               interpolated_field,
       const ComponentMask &field_comps = ComponentMask())
     {
-      if (particle_handler.n_locally_owned_particles() == 0)
+      if (particle_handler.n_locally_owned_particle_ids() == 0)
         {
           interpolated_field.compress(VectorOperation::add);
           return; // nothing else to do here
index 9c016d5f0635d911e35c4d6deeef3312b1952554..b3de65897941c95fbd23a84caa366fa619b4319c 100644 (file)
@@ -40,7 +40,7 @@ namespace Particles
         "belong to a single vector or tensor."));
 
     if ((data_component_names.size() > 0) &&
-        (particles.n_locally_owned_particles() > 0))
+        (particles.n_locally_owned_particle_ids() > 0))
       {
         Assert(
           particles.begin()->has_properties(),
@@ -75,7 +75,7 @@ namespace Particles
     const unsigned int n_property_components = data_component_names.size();
     const unsigned int n_data_components     = dataset_names.size();
 
-    patches.resize(particles.n_locally_owned_particles());
+    patches.resize(particles.n_locally_owned_particle_ids());
 
     auto particle = particles.begin();
     for (unsigned int i = 0; particle != particles.end(); ++particle, ++i)
index 514ba68a181c2541eab11c2305f2317170973064..db3be2639566c285ed3e512062afa614b67afe31 100644 (file)
@@ -767,7 +767,7 @@ namespace Particles
 
   template <int dim, int spacedim>
   types::particle_index
-  ParticleHandler<dim, spacedim>::n_locally_owned_particles() const
+  ParticleHandler<dim, spacedim>::n_locally_owned_particle_ids() const
   {
     return particles.size();
   }
@@ -794,7 +794,7 @@ namespace Particles
 
   template <int dim, int spacedim>
   IndexSet
-  ParticleHandler<dim, spacedim>::locally_relevant_ids() const
+  ParticleHandler<dim, spacedim>::locally_owned_particle_ids() const
   {
     IndexSet set(get_next_free_particle_index());
     for (const auto &p : *this)
@@ -812,7 +812,7 @@ namespace Particles
     const bool                    add_to_output_vector)
   {
     // There should be one point per particle to gather
-    AssertDimension(positions.size(), n_locally_owned_particles());
+    AssertDimension(positions.size(), n_locally_owned_particle_ids());
 
     unsigned int i = 0;
     for (auto it = begin(); it != end(); ++it, ++i)
@@ -833,7 +833,7 @@ namespace Particles
     const bool                          displace_particles)
   {
     // There should be one point per particle to fix the new position
-    AssertDimension(new_positions.size(), n_locally_owned_particles());
+    AssertDimension(new_positions.size(), n_locally_owned_particle_ids());
 
     unsigned int i = 0;
     for (auto it = begin(); it != end(); ++it, ++i)
@@ -929,7 +929,7 @@ namespace Particles
     // processes around (with an invalid cell).
 
     std::vector<particle_iterator> particles_out_of_cell;
-    particles_out_of_cell.reserve(n_locally_owned_particles());
+    particles_out_of_cell.reserve(n_locally_owned_particle_ids());
 
     // Now update the reference locations of the moved particles
     std::vector<Point<spacedim>> real_locations;
index c90a04079f6d57a273947b13ec015638556bebaa..b8d32727b69f3cce20e9147a98f57389dfda50b0 100644 (file)
@@ -36,7 +36,7 @@ namespace Particles
       const AffineConstraints<number> &                constraints,
       const ComponentMask &                            space_comps)
     {
-      if (particle_handler.n_locally_owned_particles() == 0)
+      if (particle_handler.n_locally_owned_particle_ids() == 0)
         return; // nothing to do here
 
       const auto &tria = space_dh.get_triangulation();
@@ -119,7 +119,7 @@ namespace Particles
       const AffineConstraints<typename MatrixType::value_type> &constraints,
       const ComponentMask &                                     space_comps)
     {
-      if (particle_handler.n_locally_owned_particles() == 0)
+      if (particle_handler.n_locally_owned_particle_ids() == 0)
         {
           matrix.compress(VectorOperation::add);
           return; // nothing else to do here
index d341cee3b73b34c9b1fb71517aa23d4d8829817b..4ec1b119c9b11136d05977d27ca63c4bdb3dcba4 100644 (file)
@@ -13,7 +13,7 @@
 //
 // ---------------------------------------------------------------------
 
-// Test ParticleHandler::locally_relevant_ids().
+// Test ParticleHandler::locally_owned_particle_ids().
 
 #include <deal.II/base/mpi.h>
 #include <deal.II/base/utilities.h>
@@ -63,7 +63,7 @@ test()
   auto cpu_to_index =
     particle_handler.insert_global_particles(points, global_bounding_boxes);
 
-  const auto set = particle_handler.locally_relevant_ids();
+  const auto set = particle_handler.locally_owned_particle_ids();
   deallog << "Local set: ";
   set.print(deallog);
   deallog << std::endl;
index 1a49db8ce04a10f6a005ace5c0f0d83f46a9cc1e..837e77b370222bb89f920361bfb15e4fea349764 100644 (file)
@@ -78,7 +78,7 @@ test()
             << particle_handler.n_global_particles() << std::endl;
 
   const auto n_local_particles_dofs =
-    particle_handler.n_locally_owned_particles() * n_comps;
+    particle_handler.n_locally_owned_particle_ids() * n_comps;
 
   auto particle_sizes =
     Utilities::MPI::all_gather(MPI_COMM_WORLD, n_local_particles_dofs);
index abc6bd26f194112f8ba34b54c4f3e5c26ec4b78e..0d950ca732d47b04b9606dd3af08a69024fcd793 100644 (file)
@@ -99,7 +99,7 @@ test()
     space_dh, particle_handler, dsp, constraints, space_mask);
 
   const auto n_local_particles_dofs =
-    particle_handler.n_locally_owned_particles() * n_comps;
+    particle_handler.n_locally_owned_particle_ids() * n_comps;
 
   auto particle_sizes =
     Utilities::MPI::all_gather(MPI_COMM_WORLD, n_local_particles_dofs);
index 5cacdff2213271917f2a01ad720a4077cee8da23..f06433a15fdeffb6f5a61071ef4014b7f5abc70c 100644 (file)
@@ -103,7 +103,7 @@ test()
     space_dh, particle_handler, dsp, constraints, space_mask);
 
   const auto n_local_particles_dofs =
-    particle_handler.n_locally_owned_particles() * n_comps;
+    particle_handler.n_locally_owned_particle_ids() * n_comps;
 
   auto particle_sizes =
     Utilities::MPI::all_gather(MPI_COMM_WORLD, n_local_particles_dofs);
index b3edfb11ce6f3ddefc65cc71d7112701bb718459..b89d95bb2b557885b8579ef6ba65f30fa6242660 100644 (file)
@@ -102,7 +102,7 @@ test()
     space_dh, particle_handler, dsp, constraints, space_mask);
 
   const auto n_local_particles_dofs =
-    particle_handler.n_locally_owned_particles() * n_comps;
+    particle_handler.n_locally_owned_particle_ids() * n_comps;
 
   auto particle_sizes =
     Utilities::MPI::all_gather(MPI_COMM_WORLD, n_local_particles_dofs);
index 886e848890e5f398b7c570d7bc2033b971a594e4..c1f4a263b535a013e180f0cb09a4d39e663c9885 100644 (file)
@@ -97,7 +97,7 @@ test()
     space_dh, particle_handler, dsp, constraints, space_mask);
 
   const auto n_local_particles_dofs =
-    particle_handler.n_locally_owned_particles() * n_comps;
+    particle_handler.n_locally_owned_particle_ids() * n_comps;
 
   auto particle_sizes =
     Utilities::MPI::all_gather(MPI_COMM_WORLD, n_local_particles_dofs);
index 15c2c382f6333c20ec9b85efde15e85e59df8cd1..ca1ad396f2dcf309242472650e208169a89c815d 100644 (file)
@@ -99,7 +99,7 @@ test()
     space_dh, particle_handler, dsp, constraints, space_mask);
 
   const auto n_local_particles_dofs =
-    particle_handler.n_locally_owned_particles() * n_comps;
+    particle_handler.n_locally_owned_particle_ids() * n_comps;
 
   auto particle_sizes =
     Utilities::MPI::all_gather(MPI_COMM_WORLD, n_local_particles_dofs);
index 03393d72204cb12156a1ada1451e9f108e88c79e..7bec57ff5f44cb873c2890078332470e6c58d2f2 100644 (file)
@@ -65,7 +65,7 @@ test()
   auto cpu_to_index =
     particle_handler.insert_global_particles(points, global_bounding_boxes);
 
-  const auto set = particle_handler.locally_relevant_ids().tensor_product(
+  const auto set = particle_handler.locally_owned_particle_ids().tensor_product(
     complete_index_set(spacedim));
 
   TrilinosWrappers::MPI::Vector vector(set, MPI_COMM_WORLD);
@@ -81,8 +81,9 @@ test()
 
   // Make sure we have a new index set
 
-  const auto new_set = particle_handler.locally_relevant_ids().tensor_product(
-    complete_index_set(spacedim));
+  const auto new_set =
+    particle_handler.locally_owned_particle_ids().tensor_product(
+      complete_index_set(spacedim));
 
   TrilinosWrappers::MPI::Vector new_vector(new_set, MPI_COMM_WORLD);
 
index ae6853c9f72f978bbc599efce7e25bb1613e720f..c307a18d169b1955ddd242e12e6ca7f2b794836b 100644 (file)
@@ -72,10 +72,10 @@ test()
   Tensor<1, spacedim> displacement;
   displacement[0] = 0.1;
   std::vector<Point<spacedim>> new_particle_positions(
-    particle_handler.n_locally_owned_particles());
+    particle_handler.n_locally_owned_particle_ids());
 
   std::vector<Point<spacedim>> old_particle_positions(
-    particle_handler.n_locally_owned_particles());
+    particle_handler.n_locally_owned_particle_ids());
   particle_handler.get_particle_positions(old_particle_positions);
 
   for (unsigned int i = 0; i < old_particle_positions.size(); ++i)

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.