Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
timer.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2020 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 
17 #include <deal.II/base/mpi.h>
19 #include <deal.II/base/timer.h>
20 #include <deal.II/base/utilities.h>
21 
22 #include <boost/io/ios_state.hpp>
23 
24 #include <algorithm>
25 #include <chrono>
26 #include <iomanip>
27 #include <iostream>
28 #include <map>
29 #include <sstream>
30 #include <string>
31 #include <type_traits>
32 
33 #ifdef DEAL_II_HAVE_SYS_RESOURCE_H
34 # include <sys/resource.h>
35 #endif
36 
37 #ifdef DEAL_II_MSVC
38 # include <windows.h>
39 #endif
40 
41 
42 
44 
45 namespace internal
46 {
47  namespace TimerImplementation
48  {
49  namespace
50  {
55  template <typename T>
56  struct is_duration : std::false_type
57  {};
58 
62  template <typename Rep, typename Period>
63  struct is_duration<std::chrono::duration<Rep, Period>> : std::true_type
64  {};
65 
71  template <typename T>
72  T
73  from_seconds(const double time)
74  {
75  static_assert(is_duration<T>::value,
76  "The template type should be a duration type.");
77  return T(std::lround(T::period::den * (time / T::period::num)));
78  }
79 
84  template <typename Rep, typename Period>
85  double
86  to_seconds(const std::chrono::duration<Rep, Period> duration)
87  {
88  return Period::num * double(duration.count()) / Period::den;
89  }
90 
94  void
95  clear_timing_data(Utilities::MPI::MinMaxAvg &data)
96  {
97  data.sum = numbers::signaling_nan<double>();
98  data.min = numbers::signaling_nan<double>();
99  data.max = numbers::signaling_nan<double>();
100  data.avg = numbers::signaling_nan<double>();
103  }
104  } // namespace
105  } // namespace TimerImplementation
106 } // namespace internal
107 
108 
109 
111 CPUClock::now() noexcept
112 {
113  double system_cpu_duration = 0.0;
114 #ifdef DEAL_II_MSVC
115  FILETIME cpuTime, sysTime, createTime, exitTime;
116  const auto succeeded = GetProcessTimes(
117  GetCurrentProcess(), &createTime, &exitTime, &sysTime, &cpuTime);
118  if (succeeded)
119  {
120  system_cpu_duration =
121  (double)(((unsigned long long)cpuTime.dwHighDateTime << 32) |
122  cpuTime.dwLowDateTime) /
123  1e7;
124  }
125  // keep the zero value if GetProcessTimes didn't work
126 #elif defined(DEAL_II_HAVE_SYS_RESOURCE_H)
127  rusage usage;
128  getrusage(RUSAGE_SELF, &usage);
129  system_cpu_duration = usage.ru_utime.tv_sec + 1.e-6 * usage.ru_utime.tv_usec;
130 #else
131  DEAL_II_WARNING("Unsupported platform. Porting not finished.")
132 #endif
133  return time_point(
134  internal::TimerImplementation::from_seconds<duration>(system_cpu_duration));
135 }
136 
137 
138 
139 template <typename clock_type_>
141  : current_lap_start_time(clock_type::now())
142  , accumulated_time(duration_type::zero())
143  , last_lap_time(duration_type::zero())
144 {}
145 
146 
147 
148 template <typename clock_type_>
149 void
151 {
152  current_lap_start_time = clock_type::now();
153  accumulated_time = duration_type::zero();
154  last_lap_time = duration_type::zero();
155 }
156 
157 
158 
160  : Timer(MPI_COMM_SELF, /*sync_lap_times=*/false)
161 {}
162 
163 
164 
165 Timer::Timer(MPI_Comm mpi_communicator, const bool sync_lap_times_)
166  : running(false)
167  , mpi_communicator(mpi_communicator)
168  , sync_lap_times(sync_lap_times_)
169 {
170  reset();
171  start();
172 }
173 
174 
175 
176 void
178 {
179  running = true;
180 #ifdef DEAL_II_WITH_MPI
181  if (sync_lap_times)
182  {
183  const int ierr = MPI_Barrier(mpi_communicator);
184  AssertThrowMPI(ierr);
185  }
186 #endif
187  wall_times.current_lap_start_time = wall_clock_type::now();
188  cpu_times.current_lap_start_time = cpu_clock_type::now();
189 }
190 
191 
192 
193 double
195 {
196  if (running)
197  {
198  running = false;
199 
201  wall_clock_type::now() - wall_times.current_lap_start_time;
202  cpu_times.last_lap_time =
203  cpu_clock_type::now() - cpu_times.current_lap_start_time;
204 
206  Utilities::MPI::min_max_avg(internal::TimerImplementation::to_seconds(
209  if (sync_lap_times)
210  {
212  internal::TimerImplementation::from_seconds<decltype(
213  wall_times)::duration_type>(last_lap_wall_time_data.max);
214  cpu_times.last_lap_time =
215  internal::TimerImplementation::from_seconds<decltype(
216  cpu_times)::duration_type>(
218  internal::TimerImplementation::to_seconds(
219  cpu_times.last_lap_time),
221  .max);
222  }
224  cpu_times.accumulated_time += cpu_times.last_lap_time;
226  Utilities::MPI::min_max_avg(internal::TimerImplementation::to_seconds(
229  }
230  return internal::TimerImplementation::to_seconds(cpu_times.accumulated_time);
231 }
232 
233 
234 
235 double
237 {
238  if (running)
239  {
240  const double running_time = internal::TimerImplementation::to_seconds(
241  cpu_clock_type::now() - cpu_times.current_lap_start_time +
242  cpu_times.accumulated_time);
243  return Utilities::MPI::sum(running_time, mpi_communicator);
244  }
245  else
246  {
247  return Utilities::MPI::sum(internal::TimerImplementation::to_seconds(
248  cpu_times.accumulated_time),
250  }
251 }
252 
253 
254 
255 double
257 {
258  return internal::TimerImplementation::to_seconds(cpu_times.last_lap_time);
259 }
260 
261 
262 
263 double
265 {
266  wall_clock_type::duration current_elapsed_wall_time;
267  if (running)
268  current_elapsed_wall_time = wall_clock_type::now() -
271  else
272  current_elapsed_wall_time = wall_times.accumulated_time;
273 
274  return internal::TimerImplementation::to_seconds(current_elapsed_wall_time);
275 }
276 
277 
278 
279 double
281 {
282  return internal::TimerImplementation::to_seconds(wall_times.last_lap_time);
283 }
284 
285 
286 
287 void
289 {
290  wall_times.reset();
291  cpu_times.reset();
292  running = false;
293  internal::TimerImplementation::clear_timing_data(last_lap_wall_time_data);
294  internal::TimerImplementation::clear_timing_data(accumulated_wall_time_data);
295 }
296 
297 
298 
299 /* ---------------------------- TimerOutput -------------------------- */
300 
301 TimerOutput::TimerOutput(std::ostream & stream,
302  const OutputFrequency output_frequency,
303  const OutputType output_type)
304  : output_frequency(output_frequency)
305  , output_type(output_type)
306  , out_stream(stream, true)
307  , output_is_enabled(true)
308  , mpi_communicator(MPI_COMM_SELF)
309 {}
310 
311 
312 
314  const OutputFrequency output_frequency,
315  const OutputType output_type)
316  : output_frequency(output_frequency)
317  , output_type(output_type)
318  , out_stream(stream)
319  , output_is_enabled(true)
320  , mpi_communicator(MPI_COMM_SELF)
321 {}
322 
323 
324 
326  std::ostream & stream,
327  const OutputFrequency output_frequency,
328  const OutputType output_type)
329  : output_frequency(output_frequency)
330  , output_type(output_type)
331  , out_stream(stream, true)
332  , output_is_enabled(true)
333  , mpi_communicator(mpi_communicator)
334 {}
335 
336 
337 
339  ConditionalOStream & stream,
340  const OutputFrequency output_frequency,
341  const OutputType output_type)
342  : output_frequency(output_frequency)
343  , output_type(output_type)
344  , out_stream(stream)
345  , output_is_enabled(true)
346  , mpi_communicator(mpi_communicator)
347 {}
348 
349 
350 
352 {
353  auto do_exit = [this]() {
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 __cpp_lib_uncaught_exceptions >= 201411
372  // std::uncaught_exception() is deprecated in c++17
373  if (std::uncaught_exceptions() > 0 && mpi_communicator != MPI_COMM_SELF)
374 # else
375  if (std::uncaught_exception() == true && mpi_communicator != MPI_COMM_SELF)
376 # endif
377  {
378  const unsigned int myid =
380  if (myid == 0)
381  std::cerr
382  << "---------------------------------------------------------\n"
383  << "TimerOutput objects finalize timed values printed to the\n"
384  << "screen by communicating over MPI in their destructors.\n"
385  << "Since an exception is currently uncaught, this\n"
386  << "synchronization (and subsequent output) will be skipped\n"
387  << "to avoid a possible deadlock.\n"
388  << "---------------------------------------------------------"
389  << std::endl;
390  }
391  else
392  {
393  do_exit();
394  }
395 #else
396  do_exit();
397 #endif
398 }
399 
400 
401 
402 void
403 TimerOutput::enter_subsection(const std::string &section_name)
404 {
405  std::lock_guard<std::mutex> lock(mutex);
406 
407  Assert(section_name.empty() == false, ExcMessage("Section string is empty."));
408 
409  Assert(std::find(active_sections.begin(),
410  active_sections.end(),
411  section_name) == active_sections.end(),
412  ExcMessage(std::string("Cannot enter the already active section <") +
413  section_name + ">."));
414 
415  if (sections.find(section_name) == sections.end())
416  {
417  if (mpi_communicator != MPI_COMM_SELF)
418  {
419  // create a new timer for this section. the second argument
420  // will ensure that we have an MPI barrier before starting
421  // and stopping a timer, and this ensures that we get the
422  // maximum run time for this section over all processors.
423  // The mpi_communicator from TimerOutput is passed to the
424  // Timer here, so this Timer will collect timing information
425  // among all processes inside mpi_communicator.
426  sections[section_name].timer = Timer(mpi_communicator, true);
427  }
428 
429 
430  sections[section_name].total_cpu_time = 0;
431  sections[section_name].total_wall_time = 0;
432  sections[section_name].n_calls = 0;
433  }
434 
435  sections[section_name].timer.reset();
436  sections[section_name].timer.start();
437  sections[section_name].n_calls++;
438 
439  active_sections.push_back(section_name);
440 }
441 
442 
443 
444 void
445 TimerOutput::leave_subsection(const std::string &section_name)
446 {
447  Assert(!active_sections.empty(),
448  ExcMessage("Cannot exit any section because none has been entered!"));
449 
450  std::lock_guard<std::mutex> lock(mutex);
451 
452  if (!section_name.empty())
453  {
454  Assert(sections.find(section_name) != sections.end(),
455  ExcMessage("Cannot delete a section that was never created."));
456  Assert(std::find(active_sections.begin(),
457  active_sections.end(),
458  section_name) != active_sections.end(),
459  ExcMessage("Cannot delete a section that has not been entered."));
460  }
461 
462  // if no string is given, exit the last
463  // active section.
464  const std::string actual_section_name =
465  (section_name.empty() ? active_sections.back() : section_name);
466 
467  sections[actual_section_name].timer.stop();
468  sections[actual_section_name].total_wall_time +=
469  sections[actual_section_name].timer.last_wall_time();
470 
471  // Get cpu time. On MPI systems, if constructed with an mpi_communicator
472  // like MPI_COMM_WORLD, then the Timer will sum up the CPU time between
473  // processors among the provided mpi_communicator. Therefore, no
474  // communication is needed here.
475  const double cpu_time = sections[actual_section_name].timer.last_cpu_time();
476  sections[actual_section_name].total_cpu_time += cpu_time;
477 
478  // in case we have to print out something, do that here...
479  if ((output_frequency == every_call ||
481  output_is_enabled == true)
482  {
483  std::string output_time;
484  std::ostringstream cpu;
485  cpu << cpu_time << "s";
486  std::ostringstream wall;
487  wall << sections[actual_section_name].timer.last_wall_time() << "s";
488  if (output_type == cpu_times)
489  output_time = ", CPU time: " + cpu.str();
490  else if (output_type == wall_times)
491  output_time = ", wall time: " + wall.str() + ".";
492  else
493  output_time =
494  ", CPU/wall time: " + cpu.str() + " / " + wall.str() + ".";
495 
496  out_stream << actual_section_name << output_time << std::endl;
497  }
498 
499  // delete the index from the list of
500  // active ones
501  active_sections.erase(std::find(active_sections.begin(),
502  active_sections.end(),
503  actual_section_name));
504 }
505 
506 
507 
508 std::map<std::string, double>
510 {
511  std::map<std::string, double> output;
512  for (const auto &section : sections)
513  {
514  switch (kind)
515  {
516  case TimerOutput::OutputData::total_cpu_time:
517  output[section.first] = section.second.total_cpu_time;
518  break;
519  case TimerOutput::OutputData::total_wall_time:
520  output[section.first] = section.second.total_wall_time;
521  break;
522  case TimerOutput::OutputData::n_calls:
523  output[section.first] = section.second.n_calls;
524  break;
525  default:
526  Assert(false, ExcNotImplemented());
527  }
528  }
529  return output;
530 }
531 
532 
533 
534 void
536 {
537  // we are going to change the precision and width of output below. store the
538  // old values so the get restored when exiting this function
539  const boost::io::ios_base_all_saver restore_stream(out_stream.get_stream());
540 
541  // get the maximum width among all sections
542  unsigned int max_width = 0;
543  for (const auto &i : sections)
544  max_width =
545  std::max(max_width, static_cast<unsigned int>(i.first.length()));
546 
547  // 32 is the default width until | character
548  max_width = std::max(max_width + 1, static_cast<unsigned int>(32));
549  const std::string extra_dash = std::string(max_width - 32, '-');
550  const std::string extra_space = std::string(max_width - 32, ' ');
551 
553  {
554  // in case we want to write CPU times
555  if (output_type != wall_times)
556  {
557  double total_cpu_time =
559 
560  // check that the sum of all times is less or equal than the total
561  // time. otherwise, we might have generated a lot of overhead in this
562  // function.
563  double check_time = 0.;
564  for (const auto &i : sections)
565  check_time += i.second.total_cpu_time;
566 
567  const double time_gap = check_time - total_cpu_time;
568  if (time_gap > 0.0)
569  total_cpu_time = check_time;
570 
571  // generate a nice table
572  out_stream << "\n\n"
573  << "+---------------------------------------------"
574  << extra_dash << "+------------"
575  << "+------------+\n"
576  << "| Total CPU time elapsed since start "
577  << extra_space << "|";
578  out_stream << std::setw(10) << std::setprecision(3) << std::right;
579  out_stream << total_cpu_time << "s | |\n";
580  out_stream << "| "
581  << extra_space << "| "
582  << "| |\n";
583  out_stream << "| Section " << extra_space
584  << "| no. calls |";
585  out_stream << std::setw(10);
586  out_stream << std::setprecision(3);
587  out_stream << " CPU time "
588  << " | % of total |\n";
589  out_stream << "+---------------------------------" << extra_dash
590  << "+-----------+------------"
591  << "+------------+";
592  for (const auto &i : sections)
593  {
594  std::string name_out = i.first;
595 
596  // resize the array so that it is always of the same size
597  unsigned int pos_non_space = name_out.find_first_not_of(' ');
598  name_out.erase(0, pos_non_space);
599  name_out.resize(max_width, ' ');
600  out_stream << std::endl;
601  out_stream << "| " << name_out;
602  out_stream << "| ";
603  out_stream << std::setw(9);
604  out_stream << i.second.n_calls << " |";
605  out_stream << std::setw(10);
606  out_stream << std::setprecision(3);
607  out_stream << i.second.total_cpu_time << "s |";
608  out_stream << std::setw(10);
609  if (total_cpu_time != 0)
610  {
611  // if run time was less than 0.1%, just print a zero to avoid
612  // printing silly things such as "2.45e-6%". otherwise print
613  // the actual percentage
614  const double fraction =
615  i.second.total_cpu_time / total_cpu_time;
616  if (fraction > 0.001)
617  {
618  out_stream << std::setprecision(2);
619  out_stream << fraction * 100;
620  }
621  else
622  out_stream << 0.0;
623 
624  out_stream << "% |";
625  }
626  else
627  out_stream << 0.0 << "% |";
628  }
629  out_stream << std::endl
630  << "+---------------------------------" << extra_dash
631  << "+-----------+"
632  << "------------+------------+\n"
633  << std::endl;
634 
635  if (time_gap > 0.0)
636  out_stream
637  << std::endl
638  << "Note: The sum of counted times is " << time_gap
639  << " seconds larger than the total time.\n"
640  << "(Timer function may have introduced too much overhead, or different\n"
641  << "section timers may have run at the same time.)" << std::endl;
642  }
643 
644  // in case we want to write out wallclock times
645  if (output_type != cpu_times)
646  {
648 
649  // now generate a nice table
650  out_stream << "\n\n"
651  << "+---------------------------------------------"
652  << extra_dash << "+------------"
653  << "+------------+\n"
654  << "| Total wallclock time elapsed since start "
655  << extra_space << "|";
656  out_stream << std::setw(10) << std::setprecision(3) << std::right;
657  out_stream << total_wall_time << "s | |\n";
658  out_stream << "| "
659  << extra_space << "| "
660  << "| |\n";
661  out_stream << "| Section " << extra_space
662  << "| no. calls |";
663  out_stream << std::setw(10);
664  out_stream << std::setprecision(3);
665  out_stream << " wall time | % of total |\n";
666  out_stream << "+---------------------------------" << extra_dash
667  << "+-----------+------------"
668  << "+------------+";
669  for (const auto &i : sections)
670  {
671  std::string name_out = i.first;
672 
673  // resize the array so that it is always of the same size
674  unsigned int pos_non_space = name_out.find_first_not_of(' ');
675  name_out.erase(0, pos_non_space);
676  name_out.resize(max_width, ' ');
677  out_stream << std::endl;
678  out_stream << "| " << name_out;
679  out_stream << "| ";
680  out_stream << std::setw(9);
681  out_stream << i.second.n_calls << " |";
682  out_stream << std::setw(10);
683  out_stream << std::setprecision(3);
684  out_stream << i.second.total_wall_time << "s |";
685  out_stream << std::setw(10);
686 
687  if (total_wall_time != 0)
688  {
689  // if run time was less than 0.1%, just print a zero to avoid
690  // printing silly things such as "2.45e-6%". otherwise print
691  // the actual percentage
692  const double fraction =
693  i.second.total_wall_time / total_wall_time;
694  if (fraction > 0.001)
695  {
696  out_stream << std::setprecision(2);
697  out_stream << fraction * 100;
698  }
699  else
700  out_stream << 0.0;
701 
702  out_stream << "% |";
703  }
704  else
705  out_stream << 0.0 << "% |";
706  }
707  out_stream << std::endl
708  << "+---------------------------------" << extra_dash
709  << "+-----------+"
710  << "------------+------------+\n"
711  << std::endl;
712  }
713  }
714  else
715  // output_type == cpu_and_wall_times_grouped
716  {
717  const double total_wall_time = timer_all.wall_time();
718  double total_cpu_time =
720 
721  // check that the sum of all times is less or equal than the total time.
722  // otherwise, we might have generated a lot of overhead in this function.
723  double check_time = 0.;
724 
725  for (const auto &i : sections)
726  check_time += i.second.total_cpu_time;
727 
728  const double time_gap = check_time - total_cpu_time;
729  if (time_gap > 0.0)
730  total_cpu_time = check_time;
731 
732  // generate a nice table
733  out_stream << "\n\n+---------------------------------------------"
734  << extra_dash << "+"
735  << "------------+------------+"
736  << "------------+------------+"
737  << "\n"
738  << "| Total CPU/wall time elapsed since start "
739  << extra_space << "|" << std::setw(10) << std::setprecision(3)
740  << std::right << total_cpu_time << "s | |"
741  << total_wall_time << "s | |"
742  << "\n| "
743  << extra_space << "|"
744  << " | |"
745  << " | |"
746  << "\n| Section " << extra_space
747  << "| no. calls |"
748  << " CPU time | % of total |"
749  << " wall time | % of total |"
750  << "\n+---------------------------------" << extra_dash
751  << "+-----------+"
752  << "------------+------------+"
753  << "------------+------------+" << std::endl;
754 
755  for (const auto &i : sections)
756  {
757  std::string name_out = i.first;
758 
759  // resize the array so that it is always of the same size
760  unsigned int pos_non_space = name_out.find_first_not_of(' ');
761  name_out.erase(0, pos_non_space);
762  name_out.resize(max_width, ' ');
763  out_stream << "| " << name_out << "| ";
764 
765  out_stream << std::setw(9);
766  out_stream << i.second.n_calls << " |";
767 
768  if (output_type != wall_times)
769  {
770  out_stream << std::setw(10);
771  out_stream << std::setprecision(3);
772  out_stream << i.second.total_cpu_time << "s |";
773  out_stream << std::setw(10);
774  if (total_cpu_time != 0)
775  {
776  // if run time was less than 0.1%, just print a zero to avoid
777  // printing silly things such as "2.45e-6%". otherwise print
778  // the actual percentage
779  const double fraction =
780  i.second.total_cpu_time / total_cpu_time;
781  if (fraction > 0.001)
782  {
783  out_stream << std::setprecision(2);
784  out_stream << fraction * 100;
785  }
786  else
787  out_stream << 0.0;
788 
789  out_stream << "% |";
790  }
791  else
792  out_stream << 0.0 << "% |";
793  }
794 
795  if (output_type != cpu_times)
796  {
797  out_stream << std::setw(10);
798  out_stream << std::setprecision(3);
799  out_stream << i.second.total_wall_time << "s |";
800  out_stream << std::setw(10);
801 
802  if (total_wall_time != 0)
803  {
804  // if run time was less than 0.1%, just print a zero to avoid
805  // printing silly things such as "2.45e-6%". otherwise print
806  // the actual percentage
807  const double fraction =
808  i.second.total_wall_time / total_wall_time;
809  if (fraction > 0.001)
810  {
811  out_stream << std::setprecision(2);
812  out_stream << fraction * 100;
813  }
814  else
815  out_stream << 0.0;
816 
817  out_stream << "% |";
818  }
819  else
820  out_stream << 0.0 << "% |";
821  }
822  out_stream << std::endl;
823  }
824 
825  out_stream << "+---------------------------------" << extra_dash
826  << "+-----------+"
827  << "------------+------------+"
828  << "------------+------------+" << std::endl
829  << std::endl;
830 
831  if (output_type != wall_times && time_gap > 0.0)
832  out_stream
833  << std::endl
834  << "Note: The sum of counted times is " << time_gap
835  << " seconds larger than the total time.\n"
836  << "(Timer function may have introduced too much overhead, or different\n"
837  << "section timers may have run at the same time.)" << std::endl;
838  }
839 }
840 
841 
842 
843 void
845  const double quantile) const
846 {
847  // we are going to change the precision and width of output below. store the
848  // old values so the get restored when exiting this function
849  const boost::io::ios_base_all_saver restore_stream(out_stream.get_stream());
850 
851  AssertDimension(sections.size(),
852  Utilities::MPI::max(sections.size(), mpi_comm));
853  Assert(quantile >= 0. && quantile <= 0.5,
854  ExcMessage("The quantile must be between 0 and 0.5"));
855 
856  // get the maximum width among all sections
857  unsigned int max_width = 0;
858  for (const auto &i : sections)
859  max_width =
860  std::max(max_width, static_cast<unsigned int>(i.first.length()));
861 
862  // 17 is the default width until | character
863  max_width = std::max(max_width + 1, static_cast<unsigned int>(17));
864  const std::string extra_dash = std::string(max_width - 17, '-');
865  const std::string extra_space = std::string(max_width - 17, ' ');
866 
867  // function to print data in a nice table
868  const auto print_statistics = [&](const double given_time) {
869  const unsigned int n_ranks = Utilities::MPI::n_mpi_processes(mpi_comm);
870  if (n_ranks == 1 || quantile == 0.)
871  {
873  Utilities::MPI::min_max_avg(given_time, mpi_comm);
874 
875  out_stream << std::setw(10) << std::setprecision(4) << std::right;
876  out_stream << data.min << "s ";
877  out_stream << std::setw(5) << std::right;
878  out_stream << data.min_index << (n_ranks > 99999 ? "" : " ") << "|";
879  out_stream << std::setw(10) << std::setprecision(4) << std::right;
880  out_stream << data.avg << "s |";
881  out_stream << std::setw(10) << std::setprecision(4) << std::right;
882  out_stream << data.max << "s ";
883  out_stream << std::setw(5) << std::right;
884  out_stream << data.max_index << (n_ranks > 99999 ? "" : " ") << "|\n";
885  }
886  else
887  {
888  const unsigned int my_rank = Utilities::MPI::this_mpi_process(mpi_comm);
889  std::vector<double> receive_data(my_rank == 0 ? n_ranks : 0);
890  std::vector<double> result(9);
891 #ifdef DEAL_II_WITH_MPI
892  MPI_Gather(&given_time,
893  1,
894  MPI_DOUBLE,
895  receive_data.data(),
896  1,
897  MPI_DOUBLE,
898  0,
899  mpi_comm);
900  if (my_rank == 0)
901  {
902  // fill the received data in a pair and sort; on the way, also
903  // compute the average
904  std::vector<std::pair<double, unsigned int>> data_rank;
905  data_rank.reserve(n_ranks);
906  for (unsigned int i = 0; i < n_ranks; ++i)
907  {
908  data_rank.emplace_back(receive_data[i], i);
909  result[4] += receive_data[i];
910  }
911  result[4] /= n_ranks;
912  std::sort(data_rank.begin(), data_rank.end());
913 
914  const unsigned int quantile_index =
915  static_cast<unsigned int>(std::round(quantile * n_ranks));
916  AssertIndexRange(quantile_index, data_rank.size());
917  result[0] = data_rank[0].first;
918  result[1] = data_rank[0].second;
919  result[2] = data_rank[quantile_index].first;
920  result[3] = data_rank[quantile_index].second;
921  result[5] = data_rank[n_ranks - 1 - quantile_index].first;
922  result[6] = data_rank[n_ranks - 1 - quantile_index].second;
923  result[7] = data_rank[n_ranks - 1].first;
924  result[8] = data_rank[n_ranks - 1].second;
925  }
926  MPI_Bcast(result.data(), 9, MPI_DOUBLE, 0, mpi_comm);
927 #endif
928  out_stream << std::setw(10) << std::setprecision(4) << std::right;
929  out_stream << result[0] << "s ";
930  out_stream << std::setw(5) << std::right;
931  out_stream << static_cast<unsigned int>(result[1])
932  << (n_ranks > 99999 ? "" : " ") << "|";
933  out_stream << std::setw(10) << std::setprecision(4) << std::right;
934  out_stream << result[2] << "s ";
935  out_stream << std::setw(5) << std::right;
936  out_stream << static_cast<unsigned int>(result[3])
937  << (n_ranks > 99999 ? "" : " ") << "|";
938  out_stream << std::setw(10) << std::setprecision(4) << std::right;
939  out_stream << result[4] << "s |";
940  out_stream << std::setw(10) << std::setprecision(4) << std::right;
941  out_stream << result[5] << "s ";
942  out_stream << std::setw(5) << std::right;
943  out_stream << static_cast<unsigned int>(result[6])
944  << (n_ranks > 99999 ? "" : " ") << "|";
945  out_stream << std::setw(10) << std::setprecision(4) << std::right;
946  out_stream << result[7] << "s ";
947  out_stream << std::setw(5) << std::right;
948  out_stream << static_cast<unsigned int>(result[8])
949  << (n_ranks > 99999 ? "" : " ") << "|\n";
950  }
951  };
952 
953  // in case we want to write out wallclock times
954  {
955  const unsigned int n_ranks = Utilities::MPI::n_mpi_processes(mpi_comm);
956 
957  const std::string time_rank_column = "------------------+";
958  const std::string time_rank_space = " |";
959 
960  // now generate a nice table
961  out_stream << "\n"
962  << "+------------------------------" << extra_dash << "+"
963  << time_rank_column
964  << (n_ranks > 1 && quantile > 0. ? time_rank_column : "")
965  << "------------+"
966  << (n_ranks > 1 && quantile > 0. ? time_rank_column : "")
967  << time_rank_column << "\n"
968  << "| Total wallclock time elapsed " << extra_space << "|";
969 
970  print_statistics(timer_all.wall_time());
971 
972  out_stream << "| " << extra_space << "|"
973  << time_rank_space
974  << (n_ranks > 1 && quantile > 0. ? time_rank_space : "")
975  << " "
976  << (n_ranks > 1 && quantile > 0. ? time_rank_space : "")
977  << time_rank_space << "\n";
978  out_stream << "| Section " << extra_space << "| no. calls "
979  << "| min time rank |";
980  if (n_ranks > 1 && quantile > 0.)
981  out_stream << " " << std::setw(5) << std::setprecision(2) << std::right
982  << quantile << "-tile rank |";
983  out_stream << " avg time |";
984  if (n_ranks > 1 && quantile > 0.)
985  out_stream << " " << std::setw(5) << std::setprecision(2) << std::right
986  << 1. - quantile << "-tile rank |";
987  out_stream << " max time rank |\n";
988  out_stream << "+------------------------------" << extra_dash << "+"
989  << time_rank_column
990  << (n_ranks > 1 && quantile > 0. ? time_rank_column : "")
991  << "------------+"
992  << (n_ranks > 1 && quantile > 0. ? time_rank_column : "")
993  << time_rank_column << "\n";
994  for (const auto &i : sections)
995  {
996  std::string name_out = i.first;
997 
998  // resize the array so that it is always of the same size
999  unsigned int pos_non_space = name_out.find_first_not_of(' ');
1000  name_out.erase(0, pos_non_space);
1001  name_out.resize(max_width, ' ');
1002  out_stream << "| " << name_out;
1003  out_stream << "| ";
1004  out_stream << std::setw(9);
1005  out_stream << i.second.n_calls << " |";
1006 
1007  print_statistics(i.second.total_wall_time);
1008  }
1009  out_stream << "+------------------------------" << extra_dash << "+"
1010  << time_rank_column
1011  << (n_ranks > 1 && quantile > 0. ? time_rank_column : "")
1012  << "------------+"
1013  << (n_ranks > 1 && quantile > 0. ? time_rank_column : "")
1014  << time_rank_column << "\n";
1015  }
1016 }
1017 
1018 
1019 
1020 void
1022 {
1023  output_is_enabled = false;
1024 }
1025 
1026 
1027 
1028 void
1030 {
1031  output_is_enabled = true;
1032 }
1033 
1034 void
1036 {
1037  std::lock_guard<std::mutex> lock(mutex);
1038  sections.clear();
1039  active_sections.clear();
1040  timer_all.restart();
1041 }
1042 
1044 {
1045  try
1046  {
1047  stop();
1048  }
1049  catch (...)
1050  {}
1051 }
1052 
1053 
TimerOutput::OutputType
OutputType
Definition: timer.h:640
TimerOutput::enter_subsection
void enter_subsection(const std::string &section_name)
Definition: timer.cc:403
TimerOutput::cpu_and_wall_times_grouped
@ cpu_and_wall_times_grouped
Definition: timer.h:657
TimerOutput::print_summary
void print_summary() const
Definition: timer.cc:535
TimerOutput::Scope::~Scope
~Scope()
Definition: timer.cc:1043
TimerOutput::mpi_communicator
MPI_Comm mpi_communicator
Definition: timer.h:895
TimerOutput::output_frequency
OutputFrequency output_frequency
Definition: timer.h:842
Utilities::MPI::MinMaxAvg::sum
double sum
Definition: mpi.h:693
Timer::cpu_times
ClockMeasurements< cpu_clock_type > cpu_times
Definition: timer.h:319
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
Utilities::MPI::MinMaxAvg::min_index
unsigned int min_index
Definition: mpi.h:715
TimerOutput::cpu_times
@ cpu_times
Definition: timer.h:645
CPUClock
Definition: timer.h:38
utilities.h
TimerOutput::reset
void reset()
Definition: timer.cc:1035
Timer::stop
double stop()
Definition: timer.cc:194
Timer::ClockMeasurements< CPUClock >::duration_type
typename clock_type::duration duration_type
Definition: timer.h:269
TimerOutput::output_type
OutputType output_type
Definition: timer.h:847
TimerOutput::summary
@ summary
Definition: timer.h:605
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
Timer::mpi_communicator
MPI_Comm mpi_communicator
Definition: timer.h:331
Timer::ClockMeasurements::ClockMeasurements
ClockMeasurements()
Definition: timer.cc:140
TimerOutput::output_is_enabled
bool output_is_enabled
Definition: timer.h:882
Timer::last_lap_wall_time_data
Utilities::MPI::MinMaxAvg last_lap_wall_time_data
Definition: timer.h:344
MPI_Comm
TimerOutput::wall_times
@ wall_times
Definition: timer.h:649
Utilities::MPI::MinMaxAvg::max
double max
Definition: mpi.h:705
ConditionalOStream::get_stream
std::ostream & get_stream() const
Definition: conditional_ostream.h:169
TimerOutput::Scope::stop
void stop()
Definition: timer.h:984
CPUClock::now
static time_point now() noexcept
Definition: timer.cc:111
TimerOutput::OutputData
OutputData
Definition: timer.h:620
Timer
Definition: timer.h:119
Utilities::MPI::MinMaxAvg
Definition: mpi.h:687
Timer::ClockMeasurements::current_lap_start_time
time_point_type current_lap_start_time
Definition: timer.h:275
TimerOutput::every_call_and_summary
@ every_call_and_summary
Definition: timer.h:609
TimerOutput::mutex
Threads::Mutex mutex
Definition: timer.h:901
Timer::ClockMeasurements::accumulated_time
duration_type accumulated_time
Definition: timer.h:280
Timer::wall_time
double wall_time() const
Definition: timer.cc:264
LAPACKSupport::T
static const char T
Definition: lapack_support.h:163
Timer::wall_times
ClockMeasurements< wall_clock_type > wall_times
Definition: timer.h:314
Timer::ClockMeasurements::reset
void reset()
Definition: timer.cc:150
Timer::sync_lap_times
bool sync_lap_times
Definition: timer.h:337
TimerOutput::total_cpu_time
@ total_cpu_time
Definition: timer.h:625
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
double
TimerOutput::active_sections
std::list< std::string > active_sections
Definition: timer.h:890
Timer::start
void start()
Definition: timer.cc:177
TimerOutput::timer_all
Timer timer_all
Definition: timer.h:854
mpi.h
TimerOutput::leave_subsection
void leave_subsection(const std::string &section_name="")
Definition: timer.cc:445
timer.h
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
TimerOutput::disable_output
void disable_output()
Definition: timer.cc:1021
TimerOutput::get_summary_data
std::map< std::string, double > get_summary_data(const OutputData kind) const
Definition: timer.cc:509
Utilities::MPI::sum
T sum(const T &t, const MPI_Comm &mpi_communicator)
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
TimerOutput::TimerOutput
TimerOutput(std::ostream &stream, const OutputFrequency output_frequency, const OutputType output_type)
Definition: timer.cc:301
AssertThrowMPI
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1707
exceptions.h
Timer::restart
void restart()
Definition: timer.h:910
TimerOutput::print_wall_time_statistics
void print_wall_time_statistics(const MPI_Comm mpi_comm, const double print_quantile=0.) const
Definition: timer.cc:844
value
static const bool value
Definition: dof_tools_constraints.cc:433
TimerOutput::total_wall_time
@ total_wall_time
Definition: timer.h:629
Timer::Timer
Timer()
Definition: timer.cc:159
DEAL_II_WARNING
#define DEAL_II_WARNING(desc)
Definition: config.h:432
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
TimerOutput::OutputFrequency
OutputFrequency
Definition: timer.h:596
Timer::cpu_time
double cpu_time() const
Definition: timer.cc:236
Utilities::MPI::this_mpi_process
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:128
Utilities::MPI::MinMaxAvg::avg
double avg
Definition: mpi.h:731
Utilities::MPI::MinMaxAvg::max_index
unsigned int max_index
Definition: mpi.h:725
CPUClock::time_point
std::chrono::time_point< CPUClock, duration > time_point
Definition: timer.h:60
TimerOutput::every_call
@ every_call
Definition: timer.h:601
Utilities::MPI::n_mpi_processes
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:117
Timer::running
bool running
Definition: timer.h:324
ConditionalOStream
Definition: conditional_ostream.h:81
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
TimerOutput::enable_output
void enable_output()
Definition: timer.cc:1029
LAPACKSupport::zero
static const types::blas_int zero
Definition: lapack_support.h:179
Utilities::MPI::min_max_avg
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm &mpi_communicator)
Definition: mpi.cc:91
signaling_nan.h
Timer::accumulated_wall_time_data
Utilities::MPI::MinMaxAvg accumulated_wall_time_data
Definition: timer.h:352
internal
Definition: aligned_vector.h:369
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Timer::last_cpu_time
double last_cpu_time() const
Definition: timer.cc:256
Timer::last_wall_time
double last_wall_time() const
Definition: timer.cc:280
Utilities::MPI::MinMaxAvg::min
double min
Definition: mpi.h:699
Timer::ClockMeasurements::last_lap_time
duration_type last_lap_time
Definition: timer.h:285
TimerOutput::~TimerOutput
~TimerOutput()
Definition: timer.cc:351
TimerOutput::sections
std::map< std::string, Section > sections
Definition: timer.h:871
TimerOutput::out_stream
ConditionalOStream out_stream
Definition: timer.h:876
Utilities::MPI::max
T max(const T &t, const MPI_Comm &mpi_communicator)
Timer::reset
void reset()
Definition: timer.cc:288