]> https://gitweb.dealii.org/ - dealii.git/commitdiff
allow combining CPU and wallclock tables
authorDenis Davydov <davydden@gmail.com>
Tue, 29 May 2018 21:25:27 +0000 (23:25 +0200)
committerDenis Davydov <davydden@gmail.com>
Fri, 1 Jun 2018 10:18:31 +0000 (12:18 +0200)
doc/news/changes/minor/20180529DenisDavydov_b [new file with mode: 0644]
include/deal.II/base/timer.h
source/base/timer.cc
tests/base/timer_08_b.cc [new file with mode: 0644]
tests/base/timer_08_b.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20180529DenisDavydov_b b/doc/news/changes/minor/20180529DenisDavydov_b
new file mode 100644 (file)
index 0000000..8cddb37
--- /dev/null
@@ -0,0 +1,4 @@
+New: TimerOutput::cpu_and_wall_times_grouped will print CPU and
+wallclock time in the a single table.
+<br>
+(Denis Davydov 2018/05/29)
index 46422f6a97cc2bf81098c335b44465e95b087dd5..49bd330ae846b87a41572d953c9b3e407b0be837 100644 (file)
@@ -682,9 +682,13 @@ public:
      */
     wall_times,
     /**
-     * Output both CPU and wall clock times.
+     * Output both CPU and wall clock times in separate tables.
      */
-    cpu_and_wall_times
+    cpu_and_wall_times,
+    /**
+     * Output both CPU and wall clock times in a single table.
+     */
+    cpu_and_wall_times_grouped
   };
 
   /**
index b63a08afb0585e15683e23bd96c6b1df5f9f5fee..df1b6505eca0f561d306272e1c141696e0743dc0 100644 (file)
@@ -558,15 +558,187 @@ TimerOutput::print_summary() const
   const std::string extra_dash  = std::string(max_width - 32, '-');
   const std::string extra_space = std::string(max_width - 32, ' ');
 
-  // in case we want to write CPU times
-  if (output_type != wall_times)
+  if (output_type != cpu_and_wall_times_grouped)
     {
-      double total_cpu_time =
+      // in case we want to write CPU times
+      if (output_type != wall_times)
+        {
+          double total_cpu_time =
+            Utilities::MPI::sum(timer_all(), mpi_communicator);
+
+          // check that the sum of all times is less or equal than the total
+          // time. otherwise, we might have generated a lot of overhead in this
+          // function.
+          double check_time = 0.;
+          for (std::map<std::string, Section>::const_iterator i =
+                 sections.begin();
+               i != sections.end();
+               ++i)
+            check_time += i->second.total_cpu_time;
+
+          const double time_gap = check_time - total_cpu_time;
+          if (time_gap > 0.0)
+            total_cpu_time = check_time;
+
+          // generate a nice table
+          out_stream << "\n\n"
+                     << "+---------------------------------------------"
+                     << extra_dash << "+------------"
+                     << "+------------+\n"
+                     << "| Total CPU time elapsed since start          "
+                     << extra_space << "|";
+          out_stream << std::setw(10) << std::setprecision(3) << std::right;
+          out_stream << total_cpu_time << "s |            |\n";
+          out_stream << "|                                             "
+                     << extra_space << "|            "
+                     << "|            |\n";
+          out_stream << "| Section                         " << extra_space
+                     << "| no. calls |";
+          out_stream << std::setw(10);
+          out_stream << std::setprecision(3);
+          out_stream << "  CPU time "
+                     << " | % of total |\n";
+          out_stream << "+---------------------------------" << extra_dash
+                     << "+-----------+------------"
+                     << "+------------+";
+          for (std::map<std::string, Section>::const_iterator i =
+                 sections.begin();
+               i != sections.end();
+               ++i)
+            {
+              std::string name_out = i->first;
+
+              // resize the array so that it is always of the same size
+              unsigned int pos_non_space = name_out.find_first_not_of(' ');
+              name_out.erase(0, pos_non_space);
+              name_out.resize(max_width, ' ');
+              out_stream << std::endl;
+              out_stream << "| " << name_out;
+              out_stream << "| ";
+              out_stream << std::setw(9);
+              out_stream << i->second.n_calls << " |";
+              out_stream << std::setw(10);
+              out_stream << std::setprecision(3);
+              out_stream << i->second.total_cpu_time << "s |";
+              out_stream << std::setw(10);
+              if (total_cpu_time != 0)
+                {
+                  // if run time was less than 0.1%, just print a zero to avoid
+                  // printing silly things such as "2.45e-6%". otherwise print
+                  // the actual percentage
+                  const double fraction =
+                    i->second.total_cpu_time / total_cpu_time;
+                  if (fraction > 0.001)
+                    {
+                      out_stream << std::setprecision(2);
+                      out_stream << fraction * 100;
+                    }
+                  else
+                    out_stream << 0.0;
+
+                  out_stream << "% |";
+                }
+              else
+                out_stream << 0.0 << "% |";
+            }
+          out_stream << std::endl
+                     << "+---------------------------------" << extra_dash
+                     << "+-----------+"
+                     << "------------+------------+\n"
+                     << std::endl;
+
+          if (time_gap > 0.0)
+            out_stream
+              << std::endl
+              << "Note: The sum of counted times is " << time_gap
+              << " seconds larger than the total time.\n"
+              << "(Timer function may have introduced too much overhead, or different\n"
+              << "section timers may have run at the same time.)" << std::endl;
+        }
+
+      // in case we want to write out wallclock times
+      if (output_type != cpu_times)
+        {
+          double total_wall_time = timer_all.wall_time();
+
+          // now generate a nice table
+          out_stream << "\n\n"
+                     << "+---------------------------------------------"
+                     << extra_dash << "+------------"
+                     << "+------------+\n"
+                     << "| Total wallclock time elapsed since start    "
+                     << extra_space << "|";
+          out_stream << std::setw(10) << std::setprecision(3) << std::right;
+          out_stream << total_wall_time << "s |            |\n";
+          out_stream << "|                                             "
+                     << extra_space << "|            "
+                     << "|            |\n";
+          out_stream << "| Section                         " << extra_space
+                     << "| no. calls |";
+          out_stream << std::setw(10);
+          out_stream << std::setprecision(3);
+          out_stream << "  wall time | % of total |\n";
+          out_stream << "+---------------------------------" << extra_dash
+                     << "+-----------+------------"
+                     << "+------------+";
+          for (std::map<std::string, Section>::const_iterator i =
+                 sections.begin();
+               i != sections.end();
+               ++i)
+            {
+              std::string name_out = i->first;
+
+              // resize the array so that it is always of the same size
+              unsigned int pos_non_space = name_out.find_first_not_of(' ');
+              name_out.erase(0, pos_non_space);
+              name_out.resize(max_width, ' ');
+              out_stream << std::endl;
+              out_stream << "| " << name_out;
+              out_stream << "| ";
+              out_stream << std::setw(9);
+              out_stream << i->second.n_calls << " |";
+              out_stream << std::setw(10);
+              out_stream << std::setprecision(3);
+              out_stream << i->second.total_wall_time << "s |";
+              out_stream << std::setw(10);
+
+              if (total_wall_time != 0)
+                {
+                  // if run time was less than 0.1%, just print a zero to avoid
+                  // printing silly things such as "2.45e-6%". otherwise print
+                  // the actual percentage
+                  const double fraction =
+                    i->second.total_wall_time / total_wall_time;
+                  if (fraction > 0.001)
+                    {
+                      out_stream << std::setprecision(2);
+                      out_stream << fraction * 100;
+                    }
+                  else
+                    out_stream << 0.0;
+
+                  out_stream << "% |";
+                }
+              else
+                out_stream << 0.0 << "% |";
+            }
+          out_stream << std::endl
+                     << "+---------------------------------" << extra_dash
+                     << "+-----------+"
+                     << "------------+------------+\n"
+                     << std::endl;
+        }
+    }
+  else
+    {
+      const double total_wall_time = timer_all.wall_time();
+      double       total_cpu_time =
         Utilities::MPI::sum(timer_all(), mpi_communicator);
 
       // check that the sum of all times is less or equal than the total time.
       // otherwise, we might have generated a lot of overhead in this function.
       double check_time = 0.;
+
       for (std::map<std::string, Section>::const_iterator i = sections.begin();
            i != sections.end();
            ++i)
@@ -577,26 +749,66 @@ TimerOutput::print_summary() const
         total_cpu_time = check_time;
 
       // generate a nice table
-      out_stream << "\n\n"
-                 << "+---------------------------------------------"
-                 << extra_dash << "+------------"
-                 << "+------------+\n"
-                 << "| Total CPU time elapsed since start          "
+      out_stream << "\n\n+---------------------------------------------"
+                 << extra_dash << "+";
+
+      if (output_type != wall_times)
+        out_stream << "------------+------------+";
+
+      if (output_type != cpu_times)
+        out_stream << "------------+------------+";
+
+      out_stream << "\n";
+
+      if (output_type == cpu_times)
+        out_stream << "| Total CPU time elapsed since start          ";
+      else if (output_type == wall_times)
+        out_stream << "| Total wallclock time elapsed since start    ";
+      else
+        out_stream << "| Total CPU/wall time elapsed since start     ";
+
+      out_stream << extra_space << "|";
+
+      if (output_type != wall_times)
+        out_stream << std::setw(10) << std::setprecision(3) << std::right
+                   << total_cpu_time << "s |            |";
+
+      if (output_type != cpu_times)
+        out_stream << std::setw(10) << std::setprecision(3) << std::right
+                   << total_wall_time << "s |            |";
+
+      out_stream << "\n|                                             "
                  << extra_space << "|";
-      out_stream << std::setw(10) << std::setprecision(3) << std::right;
-      out_stream << total_cpu_time << "s |            |\n";
-      out_stream << "|                                             "
-                 << extra_space << "|            "
-                 << "|            |\n";
-      out_stream << "| Section                         " << extra_space
+
+      if (output_type != wall_times)
+        out_stream << "            |            |";
+
+      if (output_type != cpu_times)
+        out_stream << "            |            |";
+
+      out_stream << "\n| Section                         " << extra_space
                  << "| no. calls |";
-      out_stream << std::setw(10);
-      out_stream << std::setprecision(3);
-      out_stream << "  CPU time "
-                 << " | % of total |\n";
-      out_stream << "+---------------------------------" << extra_dash
-                 << "+-----------+------------"
-                 << "+------------+";
+
+      // FIXME: do we need this?
+      out_stream << std::setw(10) << std::setprecision(3);
+
+      if (output_type != wall_times)
+        out_stream << "  CPU time  | % of total |";
+
+      if (output_type != cpu_times)
+        out_stream << "  wall time | % of total |";
+
+      out_stream << "\n+---------------------------------" << extra_dash
+                 << "+-----------+";
+
+      if (output_type != wall_times)
+        out_stream << "------------+------------+";
+
+      if (output_type != cpu_times)
+        out_stream << "------------+------------+";
+
+      out_stream << std::endl;
+
       for (std::map<std::string, Section>::const_iterator i = sections.begin();
            i != sections.end();
            ++i)
@@ -607,41 +819,80 @@ TimerOutput::print_summary() const
           unsigned int pos_non_space = name_out.find_first_not_of(' ');
           name_out.erase(0, pos_non_space);
           name_out.resize(max_width, ' ');
-          out_stream << std::endl;
-          out_stream << "| " << name_out;
-          out_stream << "| ";
+          out_stream << "| " << name_out << "| ";
+
           out_stream << std::setw(9);
           out_stream << i->second.n_calls << " |";
-          out_stream << std::setw(10);
-          out_stream << std::setprecision(3);
-          out_stream << i->second.total_cpu_time << "s |";
-          out_stream << std::setw(10);
-          if (total_cpu_time != 0)
+
+          if (output_type != wall_times)
             {
-              // if run time was less than 0.1%, just print a zero to avoid
-              // printing silly things such as "2.45e-6%". otherwise print the
-              // actual percentage
-              const double fraction = i->second.total_cpu_time / total_cpu_time;
-              if (fraction > 0.001)
+              out_stream << std::setw(10);
+              out_stream << std::setprecision(3);
+              out_stream << i->second.total_cpu_time << "s |";
+              out_stream << std::setw(10);
+              if (total_cpu_time != 0)
                 {
-                  out_stream << std::setprecision(2);
-                  out_stream << fraction * 100;
+                  // if run time was less than 0.1%, just print a zero to avoid
+                  // printing silly things such as "2.45e-6%". otherwise print
+                  // the actual percentage
+                  const double fraction =
+                    i->second.total_cpu_time / total_cpu_time;
+                  if (fraction > 0.001)
+                    {
+                      out_stream << std::setprecision(2);
+                      out_stream << fraction * 100;
+                    }
+                  else
+                    out_stream << 0.0;
+
+                  out_stream << "% |";
                 }
               else
-                out_stream << 0.0;
+                out_stream << 0.0 << "% |";
+            }
+
+          if (output_type != cpu_times)
+            {
+              out_stream << std::setw(10);
+              out_stream << std::setprecision(3);
+              out_stream << i->second.total_wall_time << "s |";
+              out_stream << std::setw(10);
 
-              out_stream << "% |";
+              if (total_wall_time != 0)
+                {
+                  // if run time was less than 0.1%, just print a zero to avoid
+                  // printing silly things such as "2.45e-6%". otherwise print
+                  // the actual percentage
+                  const double fraction =
+                    i->second.total_wall_time / total_wall_time;
+                  if (fraction > 0.001)
+                    {
+                      out_stream << std::setprecision(2);
+                      out_stream << fraction * 100;
+                    }
+                  else
+                    out_stream << 0.0;
+
+                  out_stream << "% |";
+                }
+              else
+                out_stream << 0.0 << "% |";
             }
-          else
-            out_stream << 0.0 << "% |";
+          out_stream << std::endl;
         }
-      out_stream << std::endl
-                 << "+---------------------------------" << extra_dash
-                 << "+-----------+"
-                 << "------------+------------+\n"
-                 << std::endl;
 
-      if (time_gap > 0.0)
+      out_stream << "+---------------------------------" << extra_dash
+                 << "+-----------+";
+
+      if (output_type != wall_times)
+        out_stream << "------------+------------+";
+
+      if (output_type != cpu_times)
+        out_stream << "------------+------------+";
+
+      out_stream << std::endl << std::endl;
+
+      if (output_type != wall_times && time_gap > 0.0)
         out_stream
           << std::endl
           << "Note: The sum of counted times is " << time_gap
@@ -650,78 +901,6 @@ TimerOutput::print_summary() const
           << "section timers may have run at the same time.)" << std::endl;
     }
 
-  // in case we want to write out wallclock times
-  if (output_type != cpu_times)
-    {
-      double total_wall_time = timer_all.wall_time();
-
-      // now generate a nice table
-      out_stream << "\n\n"
-                 << "+---------------------------------------------"
-                 << extra_dash << "+------------"
-                 << "+------------+\n"
-                 << "| Total wallclock time elapsed since start    "
-                 << extra_space << "|";
-      out_stream << std::setw(10) << std::setprecision(3) << std::right;
-      out_stream << total_wall_time << "s |            |\n";
-      out_stream << "|                                             "
-                 << extra_space << "|            "
-                 << "|            |\n";
-      out_stream << "| Section                         " << extra_space
-                 << "| no. calls |";
-      out_stream << std::setw(10);
-      out_stream << std::setprecision(3);
-      out_stream << "  wall time | % of total |\n";
-      out_stream << "+---------------------------------" << extra_dash
-                 << "+-----------+------------"
-                 << "+------------+";
-      for (std::map<std::string, Section>::const_iterator i = sections.begin();
-           i != sections.end();
-           ++i)
-        {
-          std::string name_out = i->first;
-
-          // resize the array so that it is always of the same size
-          unsigned int pos_non_space = name_out.find_first_not_of(' ');
-          name_out.erase(0, pos_non_space);
-          name_out.resize(max_width, ' ');
-          out_stream << std::endl;
-          out_stream << "| " << name_out;
-          out_stream << "| ";
-          out_stream << std::setw(9);
-          out_stream << i->second.n_calls << " |";
-          out_stream << std::setw(10);
-          out_stream << std::setprecision(3);
-          out_stream << i->second.total_wall_time << "s |";
-          out_stream << std::setw(10);
-
-          if (total_wall_time != 0)
-            {
-              // if run time was less than 0.1%, just print a zero to avoid
-              // printing silly things such as "2.45e-6%". otherwise print the
-              // actual percentage
-              const double fraction =
-                i->second.total_wall_time / total_wall_time;
-              if (fraction > 0.001)
-                {
-                  out_stream << std::setprecision(2);
-                  out_stream << fraction * 100;
-                }
-              else
-                out_stream << 0.0;
-
-              out_stream << "% |";
-            }
-          else
-            out_stream << 0.0 << "% |";
-        }
-      out_stream << std::endl
-                 << "+---------------------------------" << extra_dash
-                 << "+-----------+"
-                 << "------------+------------+\n"
-                 << std::endl;
-    }
-
   // restore previous precision and width
   out_stream.get_stream().precision(old_precision);
   out_stream.get_stream().width(old_width);
diff --git a/tests/base/timer_08_b.cc b/tests/base/timer_08_b.cc
new file mode 100644 (file)
index 0000000..de3a034
--- /dev/null
@@ -0,0 +1,75 @@
+// ---------------------------------------------------------------------
+//
+// 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.
+//
+// ---------------------------------------------------------------------
+
+
+// test TimerOutput with standard section names
+
+#include <deal.II/base/timer.h>
+
+#include <algorithm>
+#include <sstream>
+
+#include "../tests.h"
+
+// burn computer time
+double s = 0.;
+void
+burn(unsigned int n)
+{
+  for (unsigned int i = 0; i < n; ++i)
+    {
+      for (unsigned int j = 1; j < 100000; ++j)
+        {
+          s += 1. / j * i;
+        }
+    }
+}
+
+void
+test(TimerOutput::OutputType output_type)
+{
+  std::stringstream ss;
+  {
+    TimerOutput t(ss, TimerOutput::summary, output_type);
+
+    t.enter_subsection("Hello? Hello? Hello?");
+    burn(50);
+    t.leave_subsection("Hello? Hello? Hello?");
+
+    t.enter_subsection("Is there anybody in there?");
+    burn(50);
+    t.leave_subsection("Is there anybody in there?");
+  }
+
+  std::string s = ss.str();
+  std::replace_if(s.begin(), s.end(), ::isdigit, ' ');
+  std::replace_if(s.begin(), s.end(), [](char x) { return x == '.'; }, ' ');
+  deallog << s << std::endl << std::endl;
+}
+
+int
+main()
+{
+  initlog();
+
+  deallog << "cpu_times:" << std::endl;
+  test(TimerOutput::cpu_times);
+  deallog << "wall_times:" << std::endl;
+  test(TimerOutput::wall_times);
+  deallog << "cpu_and_wall_times:" << std::endl;
+  test(TimerOutput::cpu_and_wall_times);
+  deallog << "cpu_and_wall_times_grouped:" << std::endl;
+  test(TimerOutput::cpu_and_wall_times_grouped);
+}
diff --git a/tests/base/timer_08_b.output b/tests/base/timer_08_b.output
new file mode 100644 (file)
index 0000000..a0df081
--- /dev/null
@@ -0,0 +1,68 @@
+
+DEAL::cpu_times:
+DEAL::
+
++---------------------------------------------+------------+------------+
+| Total CPU time elapsed since start          |          s |            |
+|                                             |            |            |
+| Section                         | no  calls |  CPU time  | % of total |
++---------------------------------+-----------+------------+------------+
+| Hello? Hello? Hello?            |           |          s |          % |
+| Is there anybody in there?      |           |          s |          % |
++---------------------------------+-----------+------------+------------+
+
+
+DEAL::
+DEAL::wall_times:
+DEAL::
+
++---------------------------------------------+------------+------------+
+| Total wallclock time elapsed since start    |          s |            |
+|                                             |            |            |
+| Section                         | no  calls |  wall time | % of total |
++---------------------------------+-----------+------------+------------+
+| Hello? Hello? Hello?            |           |          s |          % |
+| Is there anybody in there?      |           |          s |          % |
++---------------------------------+-----------+------------+------------+
+
+
+DEAL::
+DEAL::cpu_and_wall_times:
+DEAL::
+
++---------------------------------------------+------------+------------+
+| Total CPU time elapsed since start          |          s |            |
+|                                             |            |            |
+| Section                         | no  calls |  CPU time  | % of total |
++---------------------------------+-----------+------------+------------+
+| Hello? Hello? Hello?            |           |          s |          % |
+| Is there anybody in there?      |           |          s |          % |
++---------------------------------+-----------+------------+------------+
+
+
+
++---------------------------------------------+------------+------------+
+| Total wallclock time elapsed since start    |          s |            |
+|                                             |            |            |
+| Section                         | no  calls |  wall time | % of total |
++---------------------------------+-----------+------------+------------+
+| Hello? Hello? Hello?            |           |          s |          % |
+| Is there anybody in there?      |           |          s |          % |
++---------------------------------+-----------+------------+------------+
+
+
+DEAL::
+DEAL::cpu_and_wall_times_grouped:
+DEAL::
+
++---------------------------------------------+------------+------------+------------+------------+
+| Total CPU/wall time elapsed since start     |          s |            |          s |            |
+|                                             |            |            |            |            |
+| Section                         | no  calls |  CPU time  | % of total |  wall time | % of total |
++---------------------------------+-----------+------------+------------+------------+------------+
+| Hello? Hello? Hello?            |           |          s |          % |          s |          % |
+| Is there anybody in there?      |           |          s |          % |          s |          % |
++---------------------------------+-----------+------------+------------+------------+------------+
+
+
+DEAL::

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.