]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Modify performance test to use FEPointEvaluator
authorBruno Blais <blais.bruno@gmail.com>
Tue, 4 Jul 2023 16:43:34 +0000 (12:43 -0400)
committerBruno Blais <blais.bruno@gmail.com>
Tue, 4 Jul 2023 16:43:34 +0000 (12:43 -0400)
examples/step-68/step-68.cc
tests/performance/timing_step_68.cc

index 9e23985f3ccfbfbab088d69dc8252f6d7a78d2a7..dc70cc20f120725d2df455fec2551e2d1fd16ea0 100644 (file)
@@ -543,8 +543,6 @@ namespace Step68
     Vector<double> local_dof_values(fluid_fe.dofs_per_cell);
 
     FEPointEvaluation<dim, dim> evaluator(mapping, fluid_fe, update_values);
-    const FEValuesExtractors::Vector velocity(0);
-
 
     // We loop over all the local particles. Although this could be achieved
     // directly by looping over all the cells, this would force us
index 03f7b0a519e0a06ef6959231221f07cf653f7970..5a0fff6dcc26d5ff3d51d8c48bada93f8e1c213f 100644 (file)
@@ -320,13 +320,13 @@ namespace Step68
     velocity_field.update_ghost_values();
   }
 
-
   template <int dim>
-  void
-  ParticleTracking<dim>::euler_step_interpolated(const double dt)
+  void ParticleTracking<dim>::euler_step_interpolated(const double dt)
   {
     Vector<double> local_dof_values(fluid_fe.dofs_per_cell);
 
+    FEPointEvaluation<dim, dim> evaluator(mapping, fluid_fe, update_values);
+
     auto particle = particle_handler.begin();
     while (particle != particle_handler.end())
       {
@@ -338,38 +338,37 @@ namespace Step68
 
         const auto pic = particle_handler.particles_in_cell(cell);
         Assert(pic.begin() == particle, ExcInternalError());
+        std::vector<Point<dim>> particle_positions;
         for (auto &p : pic)
-          {
-            const Point<dim> reference_location = p.get_reference_location();
-            Tensor<1, dim>   particle_velocity;
-            for (unsigned int j = 0; j < fluid_fe.dofs_per_cell; ++j)
-              {
-                const auto comp_j = fluid_fe.system_to_component_index(j);
-
-                particle_velocity[comp_j.first] +=
-                  fluid_fe.shape_value(j, reference_location) *
-                  local_dof_values[j];
-              }
-
-            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);
-
-            ArrayView<double> properties = p.get_properties();
-            for (int d = 0; d < dim; ++d)
-              properties[d] = particle_velocity[d];
-
-            properties[dim] =
-              Utilities::MPI::this_mpi_process(mpi_communicator);
-
-            ++particle;
-          }
+            particle_positions.push_back(p.get_reference_location());
+
+        evaluator.reinit(cell, particle_positions);
+        evaluator.evaluate(make_array_view(local_dof_values),
+                             EvaluationFlags::values);
+
+        for (unsigned int particle_index = 0;
+            particle != pic.end();
+            ++particle, ++particle_index)
+        {
+          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)
+            properties[d] = particle_velocity[d];
+
+          properties[dim] =
+            Utilities::MPI::this_mpi_process(mpi_communicator);
+        }
       }
   }
 
 
 
+
+
   template <int dim>
   void
   ParticleTracking<dim>::output_particles()

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.