Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
timer.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2019 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.md at
12 // the top level directory of deal.II.
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 <iomanip>
26 #include <iostream>
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  {
54  template <typename T>
55  struct is_duration : std::false_type
56  {};
57 
61  template <typename Rep, typename Period>
62  struct is_duration<std::chrono::duration<Rep, Period>> : std::true_type
63  {};
64 
70  template <typename T>
71  T
72  from_seconds(const double time)
73  {
74  static_assert(is_duration<T>::value,
75  "The template type should be a duration type.");
76  return T(std::lround(T::period::den * (time / T::period::num)));
77  }
78 
83  template <typename Rep, typename Period>
84  double
85  to_seconds(const std::chrono::duration<Rep, Period> duration)
86  {
87  return Period::num * double(duration.count()) / Period::den;
88  }
89 
93  void
94  clear_timing_data(Utilities::MPI::MinMaxAvg &data)
95  {
96  data.sum = numbers::signaling_nan<double>();
97  data.min = numbers::signaling_nan<double>();
98  data.max = numbers::signaling_nan<double>();
99  data.avg = numbers::signaling_nan<double>();
102  }
103  } // namespace
104  } // namespace TimerImplementation
105 } // namespace internal
106 
107 
108 
110 CPUClock::now() noexcept
111 {
112  double system_cpu_duration = 0.0;
113 #ifdef DEAL_II_MSVC
114  FILETIME cpuTime, sysTime, createTime, exitTime;
115  const auto succeeded = GetProcessTimes(
116  GetCurrentProcess(), &createTime, &exitTime, &sysTime, &cpuTime);
117  if (succeeded)
118  {
119  system_cpu_duration =
120  (double)(((unsigned long long)cpuTime.dwHighDateTime << 32) |
121  cpuTime.dwLowDateTime) /
122  1e7;
123  }
124  // keep the zero value if GetProcessTimes didn't work
125 #elif defined(DEAL_II_HAVE_SYS_RESOURCE_H)
126  rusage usage;
127  getrusage(RUSAGE_SELF, &usage);
128  system_cpu_duration = usage.ru_utime.tv_sec + 1.e-6 * usage.ru_utime.tv_usec;
129 #else
130  DEAL_II_WARNING("Unsupported platform. Porting not finished.")
131 #endif
132  return time_point(
133  internal::TimerImplementation::from_seconds<duration>(system_cpu_duration));
134 }
135 
136 
137 
138 template <typename clock_type_>
140  : current_lap_start_time(clock_type::now())
141  , accumulated_time(duration_type::zero())
142  , last_lap_time(duration_type::zero())
143 {}
144 
145 
146 
147 template <typename clock_type_>
148 void
150 {
151  current_lap_start_time = clock_type::now();
152  accumulated_time = duration_type::zero();
153  last_lap_time = duration_type::zero();
154 }
155 
156 
157 
159  : Timer(MPI_COMM_SELF, /*sync_lap_times=*/false)
160 {}
161 
162 
163 
164 Timer::Timer(MPI_Comm mpi_communicator, const bool sync_lap_times_)
165  : running(false)
166  , mpi_communicator(mpi_communicator)
167  , sync_lap_times(sync_lap_times_)
168 {
169  reset();
170  start();
171 }
172 
173 
174 
175 void
177 {
178  running = true;
179 #ifdef DEAL_II_WITH_MPI
180  if (sync_lap_times)
181  {
182  const int ierr = MPI_Barrier(mpi_communicator);
183  AssertThrowMPI(ierr);
184  }
185 #endif
186  wall_times.current_lap_start_time = wall_clock_type::now();
187  cpu_times.current_lap_start_time = cpu_clock_type::now();
188 }
189 
190 
191 
192 double
194 {
195  if (running)
196  {
197  running = false;
198 
200  wall_clock_type::now() - wall_times.current_lap_start_time;
201  cpu_times.last_lap_time =
202  cpu_clock_type::now() - cpu_times.current_lap_start_time;
203 
205  Utilities::MPI::min_max_avg(internal::TimerImplementation::to_seconds(
208  if (sync_lap_times)
209  {
211  internal::TimerImplementation::from_seconds<decltype(
212  wall_times)::duration_type>(last_lap_wall_time_data.max);
213  cpu_times.last_lap_time =
214  internal::TimerImplementation::from_seconds<decltype(
215  cpu_times)::duration_type>(
217  internal::TimerImplementation::to_seconds(
218  cpu_times.last_lap_time),
220  .max);
221  }
223  cpu_times.accumulated_time += cpu_times.last_lap_time;
225  Utilities::MPI::min_max_avg(internal::TimerImplementation::to_seconds(
228  }
229  return internal::TimerImplementation::to_seconds(cpu_times.accumulated_time);
230 }
231 
232 
233 
234 double
236 {
237  if (running)
238  {
239  const double running_time = internal::TimerImplementation::to_seconds(
240  cpu_clock_type::now() - cpu_times.current_lap_start_time +
241  cpu_times.accumulated_time);
242  return Utilities::MPI::sum(running_time, mpi_communicator);
243  }
244  else
245  {
246  return Utilities::MPI::sum(internal::TimerImplementation::to_seconds(
247  cpu_times.accumulated_time),
249  }
250 }
251 
252 
253 
254 double
256 {
257  return internal::TimerImplementation::to_seconds(cpu_times.last_lap_time);
258 }
259 
260 
261 
262 double
264 {
265  return internal::TimerImplementation::to_seconds(wall_times.last_lap_time);
266 }
267 
268 
269 
270 double
272 {
273  return cpu_time();
274 }
275 
276 
277 
278 double
280 {
281  wall_clock_type::duration current_elapsed_wall_time;
282  if (running)
283  current_elapsed_wall_time = wall_clock_type::now() -
286  else
287  current_elapsed_wall_time = wall_times.accumulated_time;
288 
289  return internal::TimerImplementation::to_seconds(current_elapsed_wall_time);
290 }
291 
292 
293 
294 double
296 {
297  return internal::TimerImplementation::to_seconds(wall_times.last_lap_time);
298 }
299 
300 
301 
302 void
304 {
305  wall_times.reset();
306  cpu_times.reset();
307  running = false;
308  internal::TimerImplementation::clear_timing_data(last_lap_wall_time_data);
309  internal::TimerImplementation::clear_timing_data(accumulated_wall_time_data);
310 }
311 
312 
313 
314 /* ---------------------------- TimerOutput -------------------------- */
315 
316 TimerOutput::TimerOutput(std::ostream & stream,
317  const OutputFrequency output_frequency,
318  const OutputType output_type)
319  : output_frequency(output_frequency)
320  , output_type(output_type)
321  , out_stream(stream, true)
322  , output_is_enabled(true)
323  , mpi_communicator(MPI_COMM_SELF)
324 {}
325 
326 
327 
329  const OutputFrequency output_frequency,
330  const OutputType output_type)
331  : output_frequency(output_frequency)
332  , output_type(output_type)
333  , out_stream(stream)
334  , output_is_enabled(true)
335  , mpi_communicator(MPI_COMM_SELF)
336 {}
337 
338 
339 
340 TimerOutput::TimerOutput(MPI_Comm mpi_communicator,
341  std::ostream & stream,
342  const OutputFrequency output_frequency,
343  const OutputType output_type)
344  : output_frequency(output_frequency)
345  , output_type(output_type)
346  , out_stream(stream, true)
347  , output_is_enabled(true)
348  , mpi_communicator(mpi_communicator)
349 {}
350 
351 
352 
353 TimerOutput::TimerOutput(MPI_Comm mpi_communicator,
354  ConditionalOStream & stream,
355  const OutputFrequency output_frequency,
356  const OutputType output_type)
357  : output_frequency(output_frequency)
358  , output_type(output_type)
359  , out_stream(stream)
360  , output_is_enabled(true)
361  , mpi_communicator(mpi_communicator)
362 {}
363 
364 
365 
367 {
368  auto do_exit = [this]() {
369  try
370  {
371  while (active_sections.size() > 0)
373  // don't print unless we leave all subsections
374  if ((output_frequency == summary ||
376  output_is_enabled == true)
377  print_summary();
378  }
379  catch (...)
380  {}
381  };
382 
383  // avoid communicating with other processes if there is an uncaught
384  // exception
385 #ifdef DEAL_II_WITH_MPI
386  if (std::uncaught_exception() && mpi_communicator != MPI_COMM_SELF)
387  {
388  std::cerr << "---------------------------------------------------------"
389  << std::endl
390  << "TimerOutput objects finalize timed values printed to the"
391  << std::endl
392  << "screen by communicating over MPI in their destructors."
393  << std::endl
394  << "Since an exception is currently uncaught, this" << std::endl
395  << "synchronization (and subsequent output) will be skipped to"
396  << std::endl
397  << "avoid a possible deadlock." << std::endl
398  << "---------------------------------------------------------"
399  << std::endl;
400  }
401  else
402  {
403  do_exit();
404  }
405 #else
406  do_exit();
407 #endif
408 }
409 
410 
411 
412 void
413 TimerOutput::enter_subsection(const std::string &section_name)
414 {
415  std::lock_guard<std::mutex> lock(mutex);
416 
417  Assert(section_name.empty() == false, ExcMessage("Section string is empty."));
418 
419  Assert(std::find(active_sections.begin(),
420  active_sections.end(),
421  section_name) == active_sections.end(),
422  ExcMessage(std::string("Cannot enter the already active section <") +
423  section_name + ">."));
424 
425  if (sections.find(section_name) == sections.end())
426  {
427  if (mpi_communicator != MPI_COMM_SELF)
428  {
429  // create a new timer for this section. the second argument
430  // will ensure that we have an MPI barrier before starting
431  // and stopping a timer, and this ensures that we get the
432  // maximum run time for this section over all processors.
433  // The mpi_communicator from TimerOutput is passed to the
434  // Timer here, so this Timer will collect timing information
435  // among all processes inside mpi_communicator.
436  sections[section_name].timer = Timer(mpi_communicator, true);
437  }
438 
439 
440  sections[section_name].total_cpu_time = 0;
441  sections[section_name].total_wall_time = 0;
442  sections[section_name].n_calls = 0;
443  }
444 
445  sections[section_name].timer.reset();
446  sections[section_name].timer.start();
447  sections[section_name].n_calls++;
448 
449  active_sections.push_back(section_name);
450 }
451 
452 
453 
454 void
455 TimerOutput::leave_subsection(const std::string &section_name)
456 {
457  Assert(!active_sections.empty(),
458  ExcMessage("Cannot exit any section because none has been entered!"));
459 
460  std::lock_guard<std::mutex> lock(mutex);
461 
462  if (!section_name.empty())
463  {
464  Assert(sections.find(section_name) != sections.end(),
465  ExcMessage("Cannot delete a section that was never created."));
466  Assert(std::find(active_sections.begin(),
467  active_sections.end(),
468  section_name) != active_sections.end(),
469  ExcMessage("Cannot delete a section that has not been entered."));
470  }
471 
472  // if no string is given, exit the last
473  // active section.
474  const std::string actual_section_name =
475  (section_name.empty() ? active_sections.back() : section_name);
476 
477  sections[actual_section_name].timer.stop();
478  sections[actual_section_name].total_wall_time +=
479  sections[actual_section_name].timer.last_wall_time();
480 
481  // Get cpu time. On MPI systems, if constructed with an mpi_communicator
482  // like MPI_COMM_WORLD, then the Timer will sum up the CPU time between
483  // processors among the provided mpi_communicator. Therefore, no
484  // communication is needed here.
485  const double cpu_time = sections[actual_section_name].timer.last_cpu_time();
486  sections[actual_section_name].total_cpu_time += cpu_time;
487 
488  // in case we have to print out something, do that here...
489  if ((output_frequency == every_call ||
491  output_is_enabled == true)
492  {
493  std::string output_time;
494  std::ostringstream cpu;
495  cpu << cpu_time << "s";
496  std::ostringstream wall;
497  wall << sections[actual_section_name].timer.last_wall_time() << "s";
498  if (output_type == cpu_times)
499  output_time = ", CPU time: " + cpu.str();
500  else if (output_type == wall_times)
501  output_time = ", wall time: " + wall.str() + ".";
502  else
503  output_time =
504  ", CPU/wall time: " + cpu.str() + " / " + wall.str() + ".";
505 
506  out_stream << actual_section_name << output_time << std::endl;
507  }
508 
509  // delete the index from the list of
510  // active ones
511  active_sections.erase(std::find(active_sections.begin(),
512  active_sections.end(),
513  actual_section_name));
514 }
515 
516 
517 
518 std::map<std::string, double>
520 {
521  std::map<std::string, double> output;
522  for (const auto &section : sections)
523  {
524  switch (kind)
525  {
526  case TimerOutput::OutputData::total_cpu_time:
527  output[section.first] = section.second.total_cpu_time;
528  break;
529  case TimerOutput::OutputData::total_wall_time:
530  output[section.first] = section.second.total_wall_time;
531  break;
532  case TimerOutput::OutputData::n_calls:
533  output[section.first] = section.second.n_calls;
534  break;
535  default:
536  Assert(false, ExcNotImplemented());
537  }
538  }
539  return output;
540 }
541 
542 
543 
544 void
546 {
547  // we are going to change the precision and width of output below. store the
548  // old values so we can restore it later on
549  const std::istream::fmtflags old_flags = out_stream.get_stream().flags();
550  const std::streamsize old_precision = out_stream.get_stream().precision();
551  const std::streamsize old_width = out_stream.get_stream().width();
552 
553  // get the maximum width among all sections
554  unsigned int max_width = 0;
555  for (const auto &i : sections)
556  max_width =
557  std::max(max_width, static_cast<unsigned int>(i.first.length()));
558 
559  // 32 is the default width until | character
560  max_width = std::max(max_width + 1, static_cast<unsigned int>(32));
561  const std::string extra_dash = std::string(max_width - 32, '-');
562  const std::string extra_space = std::string(max_width - 32, ' ');
563 
565  {
566  // in case we want to write CPU times
567  if (output_type != wall_times)
568  {
569  double total_cpu_time =
571 
572  // check that the sum of all times is less or equal than the total
573  // time. otherwise, we might have generated a lot of overhead in this
574  // function.
575  double check_time = 0.;
576  for (const auto &i : sections)
577  check_time += i.second.total_cpu_time;
578 
579  const double time_gap = check_time - total_cpu_time;
580  if (time_gap > 0.0)
581  total_cpu_time = check_time;
582 
583  // generate a nice table
584  out_stream << "\n\n"
585  << "+---------------------------------------------"
586  << extra_dash << "+------------"
587  << "+------------+\n"
588  << "| Total CPU time elapsed since start "
589  << extra_space << "|";
590  out_stream << std::setw(10) << std::setprecision(3) << std::right;
591  out_stream << total_cpu_time << "s | |\n";
592  out_stream << "| "
593  << extra_space << "| "
594  << "| |\n";
595  out_stream << "| Section " << extra_space
596  << "| no. calls |";
597  out_stream << std::setw(10);
598  out_stream << std::setprecision(3);
599  out_stream << " CPU time "
600  << " | % of total |\n";
601  out_stream << "+---------------------------------" << extra_dash
602  << "+-----------+------------"
603  << "+------------+";
604  for (const auto &i : sections)
605  {
606  std::string name_out = i.first;
607 
608  // resize the array so that it is always of the same size
609  unsigned int pos_non_space = name_out.find_first_not_of(' ');
610  name_out.erase(0, pos_non_space);
611  name_out.resize(max_width, ' ');
612  out_stream << std::endl;
613  out_stream << "| " << name_out;
614  out_stream << "| ";
615  out_stream << std::setw(9);
616  out_stream << i.second.n_calls << " |";
617  out_stream << std::setw(10);
618  out_stream << std::setprecision(3);
619  out_stream << i.second.total_cpu_time << "s |";
620  out_stream << std::setw(10);
621  if (total_cpu_time != 0)
622  {
623  // if run time was less than 0.1%, just print a zero to avoid
624  // printing silly things such as "2.45e-6%". otherwise print
625  // the actual percentage
626  const double fraction =
627  i.second.total_cpu_time / total_cpu_time;
628  if (fraction > 0.001)
629  {
630  out_stream << std::setprecision(2);
631  out_stream << fraction * 100;
632  }
633  else
634  out_stream << 0.0;
635 
636  out_stream << "% |";
637  }
638  else
639  out_stream << 0.0 << "% |";
640  }
641  out_stream << std::endl
642  << "+---------------------------------" << extra_dash
643  << "+-----------+"
644  << "------------+------------+\n"
645  << std::endl;
646 
647  if (time_gap > 0.0)
648  out_stream
649  << std::endl
650  << "Note: The sum of counted times is " << time_gap
651  << " seconds larger than the total time.\n"
652  << "(Timer function may have introduced too much overhead, or different\n"
653  << "section timers may have run at the same time.)" << std::endl;
654  }
655 
656  // in case we want to write out wallclock times
657  if (output_type != cpu_times)
658  {
660 
661  // now generate a nice table
662  out_stream << "\n\n"
663  << "+---------------------------------------------"
664  << extra_dash << "+------------"
665  << "+------------+\n"
666  << "| Total wallclock time elapsed since start "
667  << extra_space << "|";
668  out_stream << std::setw(10) << std::setprecision(3) << std::right;
669  out_stream << total_wall_time << "s | |\n";
670  out_stream << "| "
671  << extra_space << "| "
672  << "| |\n";
673  out_stream << "| Section " << extra_space
674  << "| no. calls |";
675  out_stream << std::setw(10);
676  out_stream << std::setprecision(3);
677  out_stream << " wall time | % of total |\n";
678  out_stream << "+---------------------------------" << extra_dash
679  << "+-----------+------------"
680  << "+------------+";
681  for (const auto &i : sections)
682  {
683  std::string name_out = i.first;
684 
685  // resize the array so that it is always of the same size
686  unsigned int pos_non_space = name_out.find_first_not_of(' ');
687  name_out.erase(0, pos_non_space);
688  name_out.resize(max_width, ' ');
689  out_stream << std::endl;
690  out_stream << "| " << name_out;
691  out_stream << "| ";
692  out_stream << std::setw(9);
693  out_stream << i.second.n_calls << " |";
694  out_stream << std::setw(10);
695  out_stream << std::setprecision(3);
696  out_stream << i.second.total_wall_time << "s |";
697  out_stream << std::setw(10);
698 
699  if (total_wall_time != 0)
700  {
701  // if run time was less than 0.1%, just print a zero to avoid
702  // printing silly things such as "2.45e-6%". otherwise print
703  // the actual percentage
704  const double fraction =
705  i.second.total_wall_time / total_wall_time;
706  if (fraction > 0.001)
707  {
708  out_stream << std::setprecision(2);
709  out_stream << fraction * 100;
710  }
711  else
712  out_stream << 0.0;
713 
714  out_stream << "% |";
715  }
716  else
717  out_stream << 0.0 << "% |";
718  }
719  out_stream << std::endl
720  << "+---------------------------------" << extra_dash
721  << "+-----------+"
722  << "------------+------------+\n"
723  << std::endl;
724  }
725  }
726  else
727  // output_type == cpu_and_wall_times_grouped
728  {
729  const double total_wall_time = timer_all.wall_time();
730  double total_cpu_time =
732 
733  // check that the sum of all times is less or equal than the total time.
734  // otherwise, we might have generated a lot of overhead in this function.
735  double check_time = 0.;
736 
737  for (const auto &i : sections)
738  check_time += i.second.total_cpu_time;
739 
740  const double time_gap = check_time - total_cpu_time;
741  if (time_gap > 0.0)
742  total_cpu_time = check_time;
743 
744  // generate a nice table
745  out_stream << "\n\n+---------------------------------------------"
746  << extra_dash << "+"
747  << "------------+------------+"
748  << "------------+------------+"
749  << "\n"
750  << "| Total CPU/wall time elapsed since start "
751  << extra_space << "|" << std::setw(10) << std::setprecision(3)
752  << std::right << total_cpu_time << "s | |"
753  << total_wall_time << "s | |"
754  << "\n| "
755  << extra_space << "|"
756  << " | |"
757  << " | |"
758  << "\n| Section " << extra_space
759  << "| no. calls |"
760  << " CPU time | % of total |"
761  << " wall time | % of total |"
762  << "\n+---------------------------------" << extra_dash
763  << "+-----------+"
764  << "------------+------------+"
765  << "------------+------------+" << std::endl;
766 
767  for (const auto &i : sections)
768  {
769  std::string name_out = i.first;
770 
771  // resize the array so that it is always of the same size
772  unsigned int pos_non_space = name_out.find_first_not_of(' ');
773  name_out.erase(0, pos_non_space);
774  name_out.resize(max_width, ' ');
775  out_stream << "| " << name_out << "| ";
776 
777  out_stream << std::setw(9);
778  out_stream << i.second.n_calls << " |";
779 
780  if (output_type != wall_times)
781  {
782  out_stream << std::setw(10);
783  out_stream << std::setprecision(3);
784  out_stream << i.second.total_cpu_time << "s |";
785  out_stream << std::setw(10);
786  if (total_cpu_time != 0)
787  {
788  // if run time was less than 0.1%, just print a zero to avoid
789  // printing silly things such as "2.45e-6%". otherwise print
790  // the actual percentage
791  const double fraction =
792  i.second.total_cpu_time / total_cpu_time;
793  if (fraction > 0.001)
794  {
795  out_stream << std::setprecision(2);
796  out_stream << fraction * 100;
797  }
798  else
799  out_stream << 0.0;
800 
801  out_stream << "% |";
802  }
803  else
804  out_stream << 0.0 << "% |";
805  }
806 
807  if (output_type != cpu_times)
808  {
809  out_stream << std::setw(10);
810  out_stream << std::setprecision(3);
811  out_stream << i.second.total_wall_time << "s |";
812  out_stream << std::setw(10);
813 
814  if (total_wall_time != 0)
815  {
816  // if run time was less than 0.1%, just print a zero to avoid
817  // printing silly things such as "2.45e-6%". otherwise print
818  // the actual percentage
819  const double fraction =
820  i.second.total_wall_time / total_wall_time;
821  if (fraction > 0.001)
822  {
823  out_stream << std::setprecision(2);
824  out_stream << fraction * 100;
825  }
826  else
827  out_stream << 0.0;
828 
829  out_stream << "% |";
830  }
831  else
832  out_stream << 0.0 << "% |";
833  }
834  out_stream << std::endl;
835  }
836 
837  out_stream << "+---------------------------------" << extra_dash
838  << "+-----------+"
839  << "------------+------------+"
840  << "------------+------------+" << std::endl
841  << std::endl;
842 
843  if (output_type != wall_times && time_gap > 0.0)
844  out_stream
845  << std::endl
846  << "Note: The sum of counted times is " << time_gap
847  << " seconds larger than the total time.\n"
848  << "(Timer function may have introduced too much overhead, or different\n"
849  << "section timers may have run at the same time.)" << std::endl;
850  }
851 
852  // restore previous precision and width
853  out_stream.get_stream().precision(old_precision);
854  out_stream.get_stream().width(old_width);
855  out_stream.get_stream().flags(old_flags);
856 }
857 
858 
859 
860 void
862 {
863  output_is_enabled = false;
864 }
865 
866 
867 
868 void
870 {
871  output_is_enabled = true;
872 }
873 
874 void
876 {
877  std::lock_guard<std::mutex> lock(mutex);
878  sections.clear();
879  active_sections.clear();
880  timer_all.restart();
881 }
882 
884 {
885  try
886  {
887  stop();
888  }
889  catch (...)
890  {}
891 }
892 
893 
894 DEAL_II_NAMESPACE_CLOSE
void print_summary() const
Definition: timer.cc:545
void reset()
Definition: timer.cc:875
static const unsigned int invalid_unsigned_int
Definition: types.h:173
double last_wall_time() const
Definition: timer.cc:295
double get_lap_time() const
Definition: timer.cc:263
Definition: timer.h:38
MPI_Comm mpi_communicator
Definition: timer.h:907
void reset()
Definition: timer.cc:303
MPI_Comm mpi_communicator
Definition: timer.h:397
STL namespace.
static time_point now() noexcept
Definition: timer.cc:110
double stop()
Definition: timer.cc:193
Timer timer_all
Definition: timer.h:866
Threads::Mutex mutex
Definition: timer.h:913
OutputType output_type
Definition: timer.h:859
bool output_is_enabled
Definition: timer.h:894
Utilities::MPI::MinMaxAvg last_lap_wall_time_data
Definition: timer.h:410
bool sync_lap_times
Definition: timer.h:403
static ::ExceptionBase & ExcMessage(std::string arg1)
T sum(const T &t, const MPI_Comm &mpi_communicator)
double wall_time() const
Definition: timer.cc:279
#define Assert(cond, exc)
Definition: exceptions.h:1407
double operator()() const
Definition: timer.cc:271
duration_type accumulated_time
Definition: timer.h:346
Timer()
Definition: timer.cc:158
typename clock_type::duration duration_type
Definition: timer.h:335
OutputFrequency
Definition: timer.h:630
ClockMeasurements< wall_clock_type > wall_times
Definition: timer.h:380
void disable_output()
Definition: timer.cc:861
bool running
Definition: timer.h:390
duration_type last_lap_time
Definition: timer.h:351
std::list< std::string > active_sections
Definition: timer.h:902
void start()
Definition: timer.cc:176
void leave_subsection(const std::string &section_name="")
Definition: timer.cc:455
double cpu_time() const
Definition: timer.cc:235
std::map< std::string, double > get_summary_data(const OutputData kind) const
Definition: timer.cc:519
void restart()
Definition: timer.h:922
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1695
double last_cpu_time() const
Definition: timer.cc:255
TimerOutput(std::ostream &stream, const OutputFrequency output_frequency, const OutputType output_type)
Definition: timer.cc:316
std::map< std::string, Section > sections
Definition: timer.h:883
ConditionalOStream out_stream
Definition: timer.h:888
std::ostream & get_stream() const
void enable_output()
Definition: timer.cc:869
std::chrono::time_point< CPUClock, duration > time_point
Definition: timer.h:60
static ::ExceptionBase & ExcNotImplemented()
Definition: timer.h:119
void enter_subsection(const std::string &section_name)
Definition: timer.cc:413
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm &mpi_communicator)
Definition: mpi.cc:429
OutputFrequency output_frequency
Definition: timer.h:854
time_point_type current_lap_start_time
Definition: timer.h:341
unsigned int min_index
Definition: mpi.h:495
Utilities::MPI::MinMaxAvg accumulated_wall_time_data
Definition: timer.h:418
~TimerOutput()
Definition: timer.cc:366
unsigned int max_index
Definition: mpi.h:505
ClockMeasurements< cpu_clock_type > cpu_times
Definition: timer.h:385