]> https://gitweb.dealii.org/ - dealii.git/commitdiff
update test 10679/head
authorTimo Heister <timo.heister@gmail.com>
Wed, 8 Jul 2020 20:02:20 +0000 (16:02 -0400)
committerTimo Heister <timo.heister@gmail.com>
Wed, 8 Jul 2020 20:22:19 +0000 (16:22 -0400)
tests/multithreading/bench_01.cc

index 756bdf2f171699eae02fd8e370aab6e00791c7b4..f315883a9d8696a61c22d3fba1436824481dfcf7 100644 (file)
@@ -48,7 +48,7 @@
 const unsigned int n_refinements = 3;
 // how many times to run each benchmark before averaging
 const unsigned int n_runs = 1;
-// maximum number of threads to test (to speed up test)
+// maximum number of threads to test (to speed up testing)
 const unsigned int n_max_threads = 4;
 
 
@@ -410,15 +410,29 @@ assemble()
 
   double reference_l2 = -1.;
   {
-    // reference
-    system_rhs = 0.;
-    WorkStream::internal::sequential::run(dof_handler.begin_active(),
-                                          dof_handler.end(),
-                                          worker,
-                                          copier,
-                                          AssemblyScratchData<dim>(fe),
-                                          AssemblyCopyData());
-    reference_l2 = system_rhs.l2_norm();
+    // reference sequential implementation
+
+    double avg = 0.;
+
+    for (unsigned int c = 0; c < n_runs; ++c)
+      {
+        system_rhs = 0.;
+        timer.reset();
+        timer.start();
+        WorkStream::internal::sequential::run(dof_handler.begin_active(),
+                                              dof_handler.end(),
+                                              worker,
+                                              copier,
+                                              AssemblyScratchData<dim>(fe),
+                                              AssemblyCopyData());
+        timer.stop();
+        const double time = timer.last_wall_time();
+        avg += time;
+        std::cout << time << " " << std::flush;
+        reference_l2 = system_rhs.l2_norm();
+      }
+    avg /= n_runs;
+    std::cout << " avg: " << avg << std::endl;
   }
 
   using Iterator = typename DoFHandler<dim>::active_cell_iterator;
@@ -523,8 +537,8 @@ assemble()
   }
 
   {
-    std::cout << "** TBB graph **" << std::endl;
-    for (unsigned int n_cores = 1 /*n_phys_cores*/; n_cores > 0; n_cores /= 2)
+    std::cout << "** TBB colored **" << std::endl;
+    for (unsigned int n_cores = n_phys_cores; n_cores > 0; n_cores /= 2)
       {
         MultithreadInfo::set_thread_limit(n_cores);
 
@@ -537,7 +551,44 @@ assemble()
             system_rhs = 0.;
             timer.reset();
             timer.start();
-            WorkStream::run(graph,
+            WorkStream::internal::tbb_colored::run(graph,
+                                                   worker,
+                                                   copier,
+                                                   AssemblyScratchData<dim>(fe),
+                                                   AssemblyCopyData(),
+                                                   8);
+
+            timer.stop();
+            const double time = timer.last_wall_time();
+            avg += time;
+            std::cout << time << " " << std::flush;
+            Assert(abs(reference_l2 - system_rhs.l2_norm()) < 1e-10,
+                   ExcInternalError());
+          }
+        avg /= n_runs;
+        std::cout << " avg: " << avg << std::endl;
+      }
+  }
+
+#endif
+
+  {
+    std::cout << "** WorkStream **" << std::endl;
+    for (unsigned int n_cores = n_phys_cores; n_cores > 0; n_cores /= 2)
+      {
+        MultithreadInfo::set_thread_limit(n_cores);
+
+        std::cout << "n_cores " << n_cores;
+        std::cout << ' ' << std::flush;
+        double avg = 0.;
+
+        for (unsigned int c = 0; c < n_runs; ++c)
+          {
+            system_rhs = 0.;
+            timer.reset();
+            timer.start();
+            WorkStream::run(dof_handler.begin_active(),
+                            dof_handler.end(),
                             worker,
                             copier,
                             AssemblyScratchData<dim>(fe),
@@ -555,7 +606,38 @@ assemble()
       }
   }
 
-#endif
+  {
+    std::cout << "** WorkStream colored **" << std::endl;
+    for (unsigned int n_cores = n_phys_cores; n_cores > 0; n_cores /= 2)
+      {
+        MultithreadInfo::set_thread_limit(n_cores);
+
+        std::cout << "n_cores " << n_cores;
+        std::cout << ' ' << std::flush;
+        double avg = 0.;
+
+        for (unsigned int c = 0; c < n_runs; ++c)
+          {
+            system_rhs = 0.;
+            timer.reset();
+            timer.start();
+            WorkStream::run(graph,
+                            worker,
+                            copier,
+                            AssemblyScratchData<dim>(fe),
+                            AssemblyCopyData());
+
+            timer.stop();
+            const double time = timer.last_wall_time();
+            avg += time;
+            std::cout << time << " " << std::flush;
+            Assert(abs(reference_l2 - system_rhs.l2_norm()) < 1e-10,
+                   ExcInternalError());
+          }
+        avg /= n_runs;
+        std::cout << " avg: " << avg << std::endl;
+      }
+  }
 }
 
 
@@ -564,6 +646,7 @@ int
 main()
 {
   initlog();
+  // remove the limit normally applied to all tests:
   MultithreadInfo::set_thread_limit();
 
   assemble<2>();

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.