Reference documentation for deal.II version 9.0.0
timer.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/exceptions.h>
17 #include <deal.II/base/mpi.h>
18 #include <deal.II/base/signaling_nan.h>
19 #include <deal.II/base/thread_management.h>
20 #include <deal.II/base/timer.h>
21 #include <deal.II/base/utilities.h>
22 
23 #include <algorithm>
24 #include <chrono>
25 #include <iostream>
26 #include <iomanip>
27 #include <map>
28 #include <sstream>
29 #include <string>
30 #include <type_traits>
31 
32 #ifdef DEAL_II_HAVE_SYS_RESOURCE_H
33 # include <sys/resource.h>
34 #endif
35 
36 #ifdef DEAL_II_MSVC
37 # include <windows.h>
38 #endif
39 
40 
41 
42 DEAL_II_NAMESPACE_OPEN
43 
44 namespace internal
45 {
46  namespace TimerImplementation
47  {
48  namespace
49  {
53  template <typename T>
54  struct is_duration : std::false_type {};
55 
59  template <typename Rep, typename Period>
60  struct is_duration<std::chrono::duration<Rep, Period>> : std::true_type {};
61 
67  template <typename T>
68  T
69  from_seconds(const double time)
70  {
71  static_assert(is_duration<T>::value,
72  "The template type should be a duration type.");
73  return T(std::lround(T::period::den*(time/T::period::num)));
74  }
75 
80  template <typename Rep, typename Period>
81  double
82  to_seconds(const std::chrono::duration<Rep, Period> duration)
83  {
84  return Period::num*double(duration.count())/Period::den;
85  }
86 
90  void
91  clear_timing_data(Utilities::MPI::MinMaxAvg &data)
92  {
93  data.sum = numbers::signaling_nan<double>();
94  data.min = numbers::signaling_nan<double>();
95  data.max = numbers::signaling_nan<double>();
96  data.avg = numbers::signaling_nan<double>();
99  }
100  }
101  }
102 }
103 
104 
105 
107 {
108  double system_cpu_duration = 0.0;
109 #ifdef DEAL_II_MSVC
110  FILETIME cpuTime, sysTime, createTime, exitTime;
111  const auto succeeded = GetProcessTimes
112  (GetCurrentProcess(), &createTime, &exitTime, &sysTime, &cpuTime);
113  if (succeeded)
114  {
115  system_cpu_duration = (double)
116  (((unsigned long long)cpuTime.dwHighDateTime << 32)
117  | cpuTime.dwLowDateTime) / 1e7;
118  }
119  // keep the zero value if GetProcessTimes didn't work
120 #elif defined(DEAL_II_HAVE_SYS_RESOURCE_H)
121  rusage usage;
122  getrusage (RUSAGE_SELF, &usage);
123  system_cpu_duration = usage.ru_utime.tv_sec + 1.e-6 * usage.ru_utime.tv_usec;
124 #else
125 # warning "Unsupported platform. Porting not finished."
126 #endif
127  return time_point(internal::TimerImplementation::from_seconds<duration>(system_cpu_duration));
128 }
129 
130 
131 
132 template <typename clock_type_>
134  :
135  current_lap_start_time(clock_type::now()),
136  accumulated_time(duration_type::zero()),
137  last_lap_time(duration_type::zero())
138 {}
139 
140 
141 
142 template <typename clock_type_>
143 void
145 {
146  current_lap_start_time = clock_type::now();
147  accumulated_time = duration_type::zero();
148  last_lap_time = duration_type::zero();
149 }
150 
151 
152 
154  :
155  Timer(MPI_COMM_SELF, /*sync_lap_times=*/false)
156 {}
157 
158 
159 
160 Timer::Timer(MPI_Comm mpi_communicator,
161  const bool sync_lap_times_)
162  :
163  running (false),
164  mpi_communicator (mpi_communicator),
165  sync_lap_times(sync_lap_times_)
166 {
167  reset();
168  start();
169 }
170 
171 
172 
174 {
175  running = true;
176 #ifdef DEAL_II_WITH_MPI
177  if (sync_lap_times)
178  {
179  const int ierr = MPI_Barrier(mpi_communicator);
180  AssertThrowMPI(ierr);
181  }
182 #endif
183  wall_times.current_lap_start_time = wall_clock_type::now();
184  cpu_times.current_lap_start_time = cpu_clock_type::now();
185 }
186 
187 
188 
189 double Timer::stop ()
190 {
191  if (running)
192  {
193  running = false;
194 
195  wall_times.last_lap_time = wall_clock_type::now() - wall_times.current_lap_start_time;
196  cpu_times.last_lap_time = cpu_clock_type::now() - cpu_times.current_lap_start_time;
197 
199  (internal::TimerImplementation::to_seconds(wall_times.last_lap_time),
201  if (sync_lap_times)
202  {
203  wall_times.last_lap_time = internal::TimerImplementation::from_seconds<decltype(wall_times)::duration_type>
205  cpu_times.last_lap_time = internal::TimerImplementation::from_seconds<decltype(cpu_times)::duration_type>
207  (internal::TimerImplementation::to_seconds(cpu_times.last_lap_time),
209  }
211  cpu_times.accumulated_time += cpu_times.last_lap_time;
213  (internal::TimerImplementation::to_seconds(wall_times.accumulated_time),
215  }
216  return internal::TimerImplementation::to_seconds(cpu_times.accumulated_time);
217 }
218 
219 
220 
221 double Timer::cpu_time() const
222 {
223  if (running)
224  {
225  const double running_time = internal::TimerImplementation::to_seconds
227  - cpu_times.current_lap_start_time
228  + cpu_times.accumulated_time);
229  return Utilities::MPI::sum (running_time, mpi_communicator);
230  }
231  else
232  {
233  return Utilities::MPI::sum (internal::TimerImplementation::to_seconds(cpu_times.accumulated_time),
235  }
236 }
237 
238 
239 
240 double Timer::last_cpu_time() const
241 {
242  return internal::TimerImplementation::to_seconds(cpu_times.last_lap_time);
243 }
244 
245 
246 
247 double Timer::get_lap_time() const
248 {
249  return internal::TimerImplementation::to_seconds(wall_times.last_lap_time);
250 }
251 
252 
253 
254 double Timer::operator() () const
255 {
256  return cpu_time();
257 }
258 
259 
260 
261 double Timer::wall_time () const
262 {
263  wall_clock_type::duration current_elapsed_wall_time;
264  if (running)
265  current_elapsed_wall_time = wall_clock_type::now()
268  else
269  current_elapsed_wall_time = wall_times.accumulated_time;
270 
271  return internal::TimerImplementation::to_seconds(current_elapsed_wall_time);
272 }
273 
274 
275 
276 double Timer::last_wall_time () const
277 {
278  return internal::TimerImplementation::to_seconds(wall_times.last_lap_time);
279 }
280 
281 
282 
284 {
285  wall_times.reset();
286  cpu_times.reset();
287  running = false;
288  internal::TimerImplementation::clear_timing_data(last_lap_wall_time_data);
289  internal::TimerImplementation::clear_timing_data(accumulated_wall_time_data);
290 }
291 
292 
293 
294 /* ---------------------------- TimerOutput -------------------------- */
295 
296 TimerOutput::TimerOutput (std::ostream &stream,
297  const OutputFrequency output_frequency,
298  const OutputType output_type)
299  :
300  output_frequency (output_frequency),
301  output_type (output_type),
302  out_stream (stream, true),
303  output_is_enabled (true),
304  mpi_communicator (MPI_COMM_SELF)
305 {}
306 
307 
308 
310  const OutputFrequency output_frequency,
311  const OutputType output_type)
312  :
313  output_frequency (output_frequency),
314  output_type (output_type),
315  out_stream (stream),
316  output_is_enabled (true),
317  mpi_communicator (MPI_COMM_SELF)
318 {}
319 
320 
321 
322 TimerOutput::TimerOutput (MPI_Comm mpi_communicator,
323  std::ostream &stream,
324  const OutputFrequency output_frequency,
325  const OutputType output_type)
326  :
327  output_frequency (output_frequency),
328  output_type (output_type),
329  out_stream (stream, true),
330  output_is_enabled (true),
331  mpi_communicator (mpi_communicator)
332 {}
333 
334 
335 
336 TimerOutput::TimerOutput (MPI_Comm mpi_communicator,
337  ConditionalOStream &stream,
338  const OutputFrequency output_frequency,
339  const OutputType output_type)
340  :
341  output_frequency (output_frequency),
342  output_type (output_type),
343  out_stream (stream),
344  output_is_enabled (true),
345  mpi_communicator (mpi_communicator)
346 {}
347 
348 
349 
351 {
352  auto do_exit = [this]()
353  {
354  try
355  {
356  while (active_sections.size() > 0)
358  // don't print unless we leave all subsections
359  if ((output_frequency == summary ||
361  && output_is_enabled == true)
362  print_summary();
363  }
364  catch (...)
365  {}
366  };
367 
368  // avoid communicating with other processes if there is an uncaught
369  // exception
370 #ifdef DEAL_II_WITH_MPI
371  if (std::uncaught_exception() && mpi_communicator != MPI_COMM_SELF)
372  {
373  std::cerr << "---------------------------------------------------------"
374  << std::endl
375  << "TimerOutput objects finalize timed values printed to the"
376  << std::endl
377  << "screen by communicating over MPI in their destructors."
378  << std::endl
379  << "Since an exception is currently uncaught, this"
380  << std::endl
381  << "synchronization (and subsequent output) will be skipped to"
382  << std::endl
383  << "avoid a possible deadlock."
384  << std::endl
385  << "---------------------------------------------------------"
386  << std::endl;
387  }
388  else
389  {
390  do_exit();
391  }
392 #else
393  do_exit();
394 #endif
395 }
396 
397 
398 
399 void
400 TimerOutput::enter_subsection (const std::string &section_name)
401 {
403 
404  Assert (section_name.empty() == false,
405  ExcMessage ("Section string is empty."));
406 
407  Assert (std::find (active_sections.begin(), active_sections.end(),
408  section_name) == active_sections.end(),
409  ExcMessage (std::string("Cannot enter the already active section <")
410  + section_name + ">."));
411 
412  if (sections.find (section_name) == sections.end())
413  {
414  if (mpi_communicator != MPI_COMM_SELF)
415  {
416  // create a new timer for this section. the second argument
417  // will ensure that we have an MPI barrier before starting
418  // and stopping a timer, and this ensures that we get the
419  // maximum run time for this section over all processors.
420  // The mpi_communicator from TimerOutput is passed to the
421  // Timer here, so this Timer will collect timing information
422  // among all processes inside mpi_communicator.
423  sections[section_name].timer = Timer(mpi_communicator, true);
424  }
425 
426 
427  sections[section_name].total_cpu_time = 0;
428  sections[section_name].total_wall_time = 0;
429  sections[section_name].n_calls = 0;
430  }
431 
432  sections[section_name].timer.reset();
433  sections[section_name].timer.start();
434  sections[section_name].n_calls++;
435 
436  active_sections.push_back (section_name);
437 }
438 
439 
440 
441 void
442 TimerOutput::leave_subsection (const std::string &section_name)
443 {
444  Assert (!active_sections.empty(),
445  ExcMessage("Cannot exit any section because none has been entered!"));
446 
448 
449  if (section_name != "")
450  {
451  Assert (sections.find (section_name) != sections.end(),
452  ExcMessage ("Cannot delete a section that was never created."));
453  Assert (std::find (active_sections.begin(), active_sections.end(),
454  section_name) != active_sections.end(),
455  ExcMessage ("Cannot delete a section that has not been entered."));
456  }
457 
458  // if no string is given, exit the last
459  // active section.
460  const std::string actual_section_name = (section_name == "" ?
461  active_sections.back () :
462  section_name);
463 
464  sections[actual_section_name].timer.stop();
465  sections[actual_section_name].total_wall_time
466  += sections[actual_section_name].timer.last_wall_time();
467 
468  // Get cpu time. On MPI systems, if constructed with an mpi_communicator
469  // like MPI_COMM_WORLD, then the Timer will sum up the CPU time between
470  // processors among the provided mpi_communicator. Therefore, no
471  // communication is needed here.
472  const double cpu_time = sections[actual_section_name].timer.last_cpu_time();
473  sections[actual_section_name].total_cpu_time += cpu_time;
474 
475  // in case we have to print out something, do that here...
477  && output_is_enabled == true)
478  {
479  std::string output_time;
480  std::ostringstream cpu;
481  cpu << cpu_time << "s";
482  std::ostringstream wall;
483  wall << sections[actual_section_name].timer.last_wall_time() << "s";
484  if (output_type == cpu_times)
485  output_time = ", CPU time: " + cpu.str();
486  else if (output_type == wall_times)
487  output_time = ", wall time: " + wall.str() + ".";
488  else
489  output_time = ", CPU/wall time: " + cpu.str() + " / " + wall.str() + ".";
490 
491  out_stream << actual_section_name << output_time
492  << std::endl;
493  }
494 
495  // delete the index from the list of
496  // active ones
497  active_sections.erase (std::find (active_sections.begin(), active_sections.end(),
498  actual_section_name));
499 }
500 
501 
502 
503 std::map<std::string, double>
505 {
506  std::map<std::string, double> output;
507  for (const auto &section : sections)
508  {
509  switch (kind)
510  {
511  case TimerOutput::OutputData::total_cpu_time:
512  output[section.first] = section.second.total_cpu_time;
513  break;
514  case TimerOutput::OutputData::total_wall_time:
515  output[section.first] = section.second.total_wall_time;
516  break;
517  case TimerOutput::OutputData::n_calls:
518  output[section.first] = section.second.n_calls;
519  break;
520  default:
521  Assert(false, ExcNotImplemented());
522  }
523  }
524  return output;
525 }
526 
527 
528 
529 void
531 {
532  // we are going to change the
533  // precision and width of output
534  // below. store the old values so we
535  // can restore it later on
536  const std::istream::fmtflags old_flags = out_stream.get_stream().flags();
537  const std::streamsize old_precision = out_stream.get_stream().precision ();
538  const std::streamsize old_width = out_stream.get_stream().width ();
539 
540  // in case we want to write CPU times
541  if (output_type != wall_times)
542  {
544 
545  // check that the sum of all times is
546  // less or equal than the total
547  // time. otherwise, we might have
548  // generated a lot of overhead in this
549  // function.
550  double check_time = 0.;
551  for (std::map<std::string, Section>::const_iterator
552  i = sections.begin(); i!=sections.end(); ++i)
553  check_time += i->second.total_cpu_time;
554 
555  const double time_gap = check_time-total_cpu_time;
556  if (time_gap > 0.0)
557  total_cpu_time = check_time;
558 
559  // generate a nice table
560  out_stream << "\n\n"
561  << "+---------------------------------------------+------------"
562  << "+------------+\n"
563  << "| Total CPU time elapsed since start |";
564  out_stream << std::setw(10) << std::setprecision(3) << std::right;
565  out_stream << total_cpu_time << "s | |\n";
566  out_stream << "| | "
567  << "| |\n";
568  out_stream << "| Section | no. calls |";
569  out_stream << std::setw(10);
570  out_stream << std::setprecision(3);
571  out_stream << " CPU time " << " | % of total |\n";
572  out_stream << "+---------------------------------+-----------+------------"
573  << "+------------+";
574  for (std::map<std::string, Section>::const_iterator
575  i = sections.begin(); i!=sections.end(); ++i)
576  {
577  std::string name_out = i->first;
578 
579  // resize the array so that it is always
580  // of the same size
581  unsigned int pos_non_space = name_out.find_first_not_of (' ');
582  name_out.erase(0, pos_non_space);
583  name_out.resize (32, ' ');
584  out_stream << std::endl;
585  out_stream << "| " << name_out;
586  out_stream << "| ";
587  out_stream << std::setw(9);
588  out_stream << i->second.n_calls << " |";
589  out_stream << std::setw(10);
590  out_stream << std::setprecision(3);
591  out_stream << i->second.total_cpu_time << "s |";
592  out_stream << std::setw(10);
593  if (total_cpu_time != 0)
594  {
595  // if run time was less than 0.1%, just print a zero to avoid
596  // printing silly things such as "2.45e-6%". otherwise print
597  // the actual percentage
598  const double fraction = i->second.total_cpu_time/total_cpu_time;
599  if (fraction > 0.001)
600  {
601  out_stream << std::setprecision(2);
602  out_stream << fraction * 100;
603  }
604  else
605  out_stream << 0.0;
606 
607  out_stream << "% |";
608  }
609  else
610  out_stream << 0.0 << "% |";
611  }
612  out_stream << std::endl
613  << "+---------------------------------+-----------+"
614  << "------------+------------+\n"
615  << std::endl;
616 
617  if (time_gap > 0.0)
618  out_stream << std::endl
619  << "Note: The sum of counted times is " << time_gap
620  << " seconds larger than the total time.\n"
621  << "(Timer function may have introduced too much overhead, or different\n"
622  << "section timers may have run at the same time.)" << std::endl;
623  }
624 
625  // in case we want to write out wallclock times
626  if (output_type != cpu_times)
627  {
629 
630  // now generate a nice table
631  out_stream << "\n\n"
632  << "+---------------------------------------------+------------"
633  << "+------------+\n"
634  << "| Total wallclock time elapsed since start |";
635  out_stream << std::setw(10) << std::setprecision(3) << std::right;
636  out_stream << total_wall_time << "s | |\n";
637  out_stream << "| | "
638  << "| |\n";
639  out_stream << "| Section | no. calls |";
640  out_stream << std::setw(10);
641  out_stream << std::setprecision(3);
642  out_stream << " wall time | % of total |\n";
643  out_stream << "+---------------------------------+-----------+------------"
644  << "+------------+";
645  for (std::map<std::string, Section>::const_iterator
646  i = sections.begin(); i!=sections.end(); ++i)
647  {
648  std::string name_out = i->first;
649 
650  // resize the array so that it is always
651  // of the same size
652  unsigned int pos_non_space = name_out.find_first_not_of (' ');
653  name_out.erase(0, pos_non_space);
654  name_out.resize (32, ' ');
655  out_stream << std::endl;
656  out_stream << "| " << name_out;
657  out_stream << "| ";
658  out_stream << std::setw(9);
659  out_stream << i->second.n_calls << " |";
660  out_stream << std::setw(10);
661  out_stream << std::setprecision(3);
662  out_stream << i->second.total_wall_time << "s |";
663  out_stream << std::setw(10);
664 
665  if (total_wall_time != 0)
666  {
667  // if run time was less than 0.1%, just print a zero to avoid
668  // printing silly things such as "2.45e-6%". otherwise print
669  // the actual percentage
670  const double fraction = i->second.total_wall_time/total_wall_time;
671  if (fraction > 0.001)
672  {
673  out_stream << std::setprecision(2);
674  out_stream << fraction * 100;
675  }
676  else
677  out_stream << 0.0;
678 
679  out_stream << "% |";
680  }
681  else
682  out_stream << 0.0 << "% |";
683  }
684  out_stream << std::endl
685  << "+---------------------------------+-----------+"
686  << "------------+------------+\n"
687  << std::endl;
688  }
689 
690  // restore previous precision and width
691  out_stream.get_stream().precision (old_precision);
692  out_stream.get_stream().width (old_width);
693  out_stream.get_stream().flags (old_flags);
694 }
695 
696 
697 
698 void
700 {
701  output_is_enabled = false;
702 }
703 
704 
705 
706 void
708 {
709  output_is_enabled = true;
710 }
711 
712 void
714 {
716  sections.clear();
717  active_sections.clear();
718  timer_all.restart();
719 }
720 
722 {
723  try
724  {
725  stop();
726  }
727  catch (...)
728  {}
729 }
730 
731 
732 DEAL_II_NAMESPACE_CLOSE
void print_summary() const
Definition: timer.cc:530
void reset()
Definition: timer.cc:713
static const unsigned int invalid_unsigned_int
Definition: types.h:173
double last_wall_time() const
Definition: timer.cc:276
double get_lap_time() const
Definition: timer.cc:247
Definition: timer.h:37
MPI_Comm mpi_communicator
Definition: timer.h:874
void reset()
Definition: timer.cc:283
MPI_Comm mpi_communicator
Definition: timer.h:379
STL namespace.
static time_point now() noexcept
Definition: timer.cc:106
double stop()
Definition: timer.cc:189
clock_type::duration duration_type
Definition: timer.h:318
Timer timer_all
Definition: timer.h:833
Threads::Mutex mutex
Definition: timer.h:880
OutputType output_type
Definition: timer.h:826
bool output_is_enabled
Definition: timer.h:861
Utilities::MPI::MinMaxAvg last_lap_wall_time_data
Definition: timer.h:392
bool sync_lap_times
Definition: timer.h:385
static ::ExceptionBase & ExcMessage(std::string arg1)
T sum(const T &t, const MPI_Comm &mpi_communicator)
double wall_time() const
Definition: timer.cc:261
#define Assert(cond, exc)
Definition: exceptions.h:1142
double operator()() const
Definition: timer.cc:254
duration_type accumulated_time
Definition: timer.h:329
Timer()
Definition: timer.cc:153
OutputFrequency
Definition: timer.h:611
std::chrono::time_point< CPUClock, duration > time_point
Definition: timer.h:59
ClockMeasurements< wall_clock_type > wall_times
Definition: timer.h:362
void disable_output()
Definition: timer.cc:699
bool running
Definition: timer.h:372
duration_type last_lap_time
Definition: timer.h:334
std::list< std::string > active_sections
Definition: timer.h:869
void start()
Definition: timer.cc:173
double cpu_time() const
Definition: timer.cc:221
std::map< std::string, double > get_summary_data(const OutputData kind) const
Definition: timer.cc:504
void restart()
Definition: timer.h:889
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1314
double last_cpu_time() const
Definition: timer.cc:240
TimerOutput(std::ostream &stream, const OutputFrequency output_frequency, const OutputType output_type)
Definition: timer.cc:296
std::map< std::string, Section > sections
Definition: timer.h:850
ConditionalOStream out_stream
Definition: timer.h:855
std::ostream & get_stream() const
void enable_output()
Definition: timer.cc:707
static ::ExceptionBase & ExcNotImplemented()
Definition: timer.h:117
void enter_subsection(const std::string &section_name)
Definition: timer.cc:400
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm &mpi_communicator)
Definition: mpi.cc:275
OutputFrequency output_frequency
Definition: timer.h:821
time_point_type current_lap_start_time
Definition: timer.h:324
unsigned int min_index
Definition: mpi.h:391
Utilities::MPI::MinMaxAvg accumulated_wall_time_data
Definition: timer.h:400
~TimerOutput()
Definition: timer.cc:350
unsigned int max_index
Definition: mpi.h:401
ClockMeasurements< cpu_clock_type > cpu_times
Definition: timer.h:367
void leave_subsection(const std::string &section_name=std::string())
Definition: timer.cc:442