void
TimerOutput::print_summary() const
{
- // we are going to change the
- // precision and width of output
- // below. store the old values so we
- // can restore it later on
+ // we are going to change the precision and width of output below. store the
+ // old values so we can restore it later on
const std::istream::fmtflags old_flags = out_stream.get_stream().flags();
const std::streamsize old_precision = out_stream.get_stream().precision();
const std::streamsize old_width = out_stream.get_stream().width();
+ // get the maximum width among all sections
+ unsigned int max_width = 0;
+ for (std::map<std::string, Section>::const_iterator i = sections.begin();
+ i != sections.end();
+ ++i)
+ max_width = std::max(max_width, (unsigned int)i->first.length());
+
+ // 32 is the default width until | character
+ max_width = std::max(max_width + 1, (unsigned int)32);
+ 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)
{
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.
+ // 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();
total_cpu_time = check_time;
// generate a nice table
- out_stream
- << "\n\n"
- << "+---------------------------------------------+------------"
- << "+------------+\n"
- << "| Total CPU time elapsed since start |";
+ 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
- << "| | "
- << "| |\n";
- out_stream << "| Section | no. calls |";
+ 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
- << "+---------------------------------+-----------+------------"
- << "+------------+";
+ 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
+ // 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(32, ' ');
+ name_out.resize(max_width, ' ');
out_stream << std::endl;
out_stream << "| " << name_out;
out_stream << "| ";
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
+ // 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 << 0.0 << "% |";
}
out_stream << std::endl
- << "+---------------------------------+-----------+"
+ << "+---------------------------------" << extra_dash
+ << "+-----------+"
<< "------------+------------+\n"
<< std::endl;
double total_wall_time = timer_all.wall_time();
// now generate a nice table
- out_stream
- << "\n\n"
- << "+---------------------------------------------+------------"
- << "+------------+\n"
- << "| Total wallclock time elapsed since start |";
+ 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
- << "| | "
- << "| |\n";
- out_stream << "| Section | no. calls |";
+ 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
- << "+---------------------------------+-----------+------------"
- << "+------------+";
+ 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
+ // 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(32, ' ');
+ name_out.resize(max_width, ' ');
out_stream << std::endl;
out_stream << "| " << name_out;
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
+ // 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 << 0.0 << "% |";
}
out_stream << std::endl
- << "+---------------------------------+-----------+"
+ << "+---------------------------------" << extra_dash
+ << "+-----------+"
<< "------------+------------+\n"
<< std::endl;
}
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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 large 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("hi");
+ burn(50);
+ t.leave_subsection("hi");
+
+ t.enter_subsection(
+ "this is a very long section name that previously did not work");
+ burn(50);
+ t.leave_subsection(
+ "this is a very long section name that previously did not work");
+ }
+
+ 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()
+{
+ std::ofstream logfile("output");
+ deallog.attach(logfile);
+
+ 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);
+}
--- /dev/null
+
+DEAL::cpu_times:
+DEAL::
+
++---------------------------------------------------------------------------+------------+------------+
+| Total CPU time elapsed since start | s | |
+| | | |
+| Section | no calls | CPU time | % of total |
++---------------------------------------------------------------+-----------+------------+------------+
+| hi | | s | % |
+| this is a very long section name that previously did not work | | s | % |
++---------------------------------------------------------------+-----------+------------+------------+
+
+
+DEAL::
+DEAL::wall_times:
+DEAL::
+
++---------------------------------------------------------------------------+------------+------------+
+| Total wallclock time elapsed since start | s | |
+| | | |
+| Section | no calls | wall time | % of total |
++---------------------------------------------------------------+-----------+------------+------------+
+| hi | | s | % |
+| this is a very long section name that previously did not work | | s | % |
++---------------------------------------------------------------+-----------+------------+------------+
+
+
+DEAL::
+DEAL::cpu_and_wall_times:
+DEAL::
+
++---------------------------------------------------------------------------+------------+------------+
+| Total CPU time elapsed since start | s | |
+| | | |
+| Section | no calls | CPU time | % of total |
++---------------------------------------------------------------+-----------+------------+------------+
+| hi | | s | % |
+| this is a very long section name that previously did not work | | s | % |
++---------------------------------------------------------------+-----------+------------+------------+
+
+
+
++---------------------------------------------------------------------------+------------+------------+
+| Total wallclock time elapsed since start | s | |
+| | | |
+| Section | no calls | wall time | % of total |
++---------------------------------------------------------------+-----------+------------+------------+
+| hi | | s | % |
+| this is a very long section name that previously did not work | | s | % |
++---------------------------------------------------------------+-----------+------------+------------+
+
+
+DEAL::