Reference documentation for deal.II version GIT relicensing-218-g1f873f3929 2024-03-28 15:00:02+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
timer.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1998 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
16#include <deal.II/base/mpi.h>
18#include <deal.II/base/timer.h>
20
21#include <boost/io/ios_state.hpp>
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
43
44namespace 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
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
138template <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
147template <typename clock_type_>
148void
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
164Timer::Timer(const 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
175void
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
192double
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<
212 decltype(wall_times)::duration_type>(last_lap_wall_time_data.max);
213 cpu_times.last_lap_time = internal::TimerImplementation::from_seconds<
214 decltype(cpu_times)::duration_type>(
216 internal::TimerImplementation::to_seconds(
217 cpu_times.last_lap_time),
219 .max);
220 }
222 cpu_times.accumulated_time += cpu_times.last_lap_time;
224 Utilities::MPI::min_max_avg(internal::TimerImplementation::to_seconds(
227 }
228 return internal::TimerImplementation::to_seconds(cpu_times.accumulated_time);
229}
230
231
232
233double
235{
236 if (running)
237 {
238 const double running_time = internal::TimerImplementation::to_seconds(
239 cpu_clock_type::now() - cpu_times.current_lap_start_time +
240 cpu_times.accumulated_time);
241 return Utilities::MPI::sum(running_time, mpi_communicator);
242 }
243 else
244 {
245 return Utilities::MPI::sum(internal::TimerImplementation::to_seconds(
246 cpu_times.accumulated_time),
248 }
249}
250
251
252
253double
255{
256 return internal::TimerImplementation::to_seconds(cpu_times.last_lap_time);
257}
258
259
260
261double
263{
264 wall_clock_type::duration current_elapsed_wall_time;
265 if (running)
266 current_elapsed_wall_time = wall_clock_type::now() -
269 else
270 current_elapsed_wall_time = wall_times.accumulated_time;
271
272 return internal::TimerImplementation::to_seconds(current_elapsed_wall_time);
273}
274
275
276
277double
279{
280 return internal::TimerImplementation::to_seconds(wall_times.last_lap_time);
281}
282
283
284
285void
287{
289 cpu_times.reset();
290 running = false;
291 internal::TimerImplementation::clear_timing_data(last_lap_wall_time_data);
292 internal::TimerImplementation::clear_timing_data(accumulated_wall_time_data);
293}
294
296
297/* ---------------------------- TimerOutput -------------------------- */
298
299TimerOutput::TimerOutput(std::ostream &stream,
300 const OutputFrequency output_frequency,
301 const OutputType output_type)
302 : output_frequency(output_frequency)
303 , output_type(output_type)
304 , out_stream(stream, true)
305 , output_is_enabled(true)
306 , mpi_communicator(MPI_COMM_SELF)
307{}
308
309
310
312 const OutputFrequency output_frequency,
313 const OutputType output_type)
314 : output_frequency(output_frequency)
315 , output_type(output_type)
316 , out_stream(stream)
317 , output_is_enabled(true)
318 , mpi_communicator(MPI_COMM_SELF)
319{}
320
321
322
323TimerOutput::TimerOutput(const MPI_Comm mpi_communicator,
324 std::ostream &stream,
325 const OutputFrequency output_frequency,
326 const OutputType output_type)
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
336TimerOutput::TimerOutput(const MPI_Comm mpi_communicator,
337 ConditionalOStream &stream,
338 const OutputFrequency output_frequency,
339 const OutputType output_type)
340 : output_frequency(output_frequency)
341 , output_type(output_type)
342 , out_stream(stream)
343 , output_is_enabled(true)
344 , mpi_communicator(mpi_communicator)
345{}
346
347
348
350{
351 auto do_exit = [this]() {
352 try
353 {
354 while (active_sections.size() > 0)
356 // don't print unless we leave all subsections
357 if ((output_frequency == summary ||
359 output_is_enabled == true)
361 }
362 catch (...)
363 {}
364 };
365
366 // avoid communicating with other processes if there is an uncaught
367 // exception
368#ifdef DEAL_II_WITH_MPI
369# if __cpp_lib_uncaught_exceptions >= 201411
370 // std::uncaught_exception() is deprecated in c++17
371 if (std::uncaught_exceptions() > 0 && mpi_communicator != MPI_COMM_SELF)
372# else
373 if (std::uncaught_exception() == true && mpi_communicator != MPI_COMM_SELF)
374# endif
375 {
376 const unsigned int myid =
378 if (myid == 0)
379 std::cerr
380 << "---------------------------------------------------------\n"
381 << "TimerOutput objects finalize timed values printed to the\n"
382 << "screen by communicating over MPI in their destructors.\n"
383 << "Since an exception is currently uncaught, this\n"
384 << "synchronization (and subsequent output) will be skipped\n"
385 << "to avoid a possible deadlock.\n"
386 << "---------------------------------------------------------"
387 << std::endl;
388 }
389 else
390 {
391 do_exit();
392 }
393#else
394 do_exit();
395#endif
396}
397
398
399
400void
401TimerOutput::enter_subsection(const std::string &section_name)
402{
403 std::lock_guard<std::mutex> lock(mutex);
404
405 Assert(section_name.empty() == false, ExcMessage("Section string is empty."));
406
407 Assert(std::find(active_sections.begin(),
408 active_sections.end(),
409 section_name) == active_sections.end(),
410 ExcMessage(std::string("Cannot enter the already active section <") +
411 section_name + ">."));
412
413 if (sections.find(section_name) == sections.end())
414 {
415 if (mpi_communicator != MPI_COMM_SELF)
416 {
417 // create a new timer for this section. the second argument
418 // will ensure that we have an MPI barrier before starting
419 // and stopping a timer, and this ensures that we get the
420 // maximum run time for this section over all processors.
421 // The mpi_communicator from TimerOutput is passed to the
422 // Timer here, so this Timer will collect timing information
423 // among all processes inside mpi_communicator.
424 sections[section_name].timer = Timer(mpi_communicator, true);
425 }
426
427
428 sections[section_name].total_cpu_time = 0;
429 sections[section_name].total_wall_time = 0;
430 sections[section_name].n_calls = 0;
431 }
432
433 sections[section_name].timer.reset();
434 sections[section_name].timer.start();
435 ++sections[section_name].n_calls;
436
437 active_sections.push_back(section_name);
438}
439
440
441
442void
443TimerOutput::leave_subsection(const std::string &section_name)
444{
445 Assert(!active_sections.empty(),
446 ExcMessage("Cannot exit any section because none has been entered!"));
447
448 std::lock_guard<std::mutex> lock(mutex);
449
450 if (!section_name.empty())
451 {
452 Assert(sections.find(section_name) != sections.end(),
453 ExcMessage("Cannot delete a section that was never created."));
454 Assert(std::find(active_sections.begin(),
455 active_sections.end(),
456 section_name) != active_sections.end(),
457 ExcMessage("Cannot delete a section that has not been entered."));
458 }
459
460 // if no string is given, exit the last
461 // active section.
462 const std::string actual_section_name =
463 (section_name.empty() ? active_sections.back() : section_name);
464
465 sections[actual_section_name].timer.stop();
466 sections[actual_section_name].total_wall_time +=
467 sections[actual_section_name].timer.last_wall_time();
468
469 // Get cpu time. On MPI systems, if constructed with an mpi_communicator
470 // like MPI_COMM_WORLD, then the Timer will sum up the CPU time between
471 // processors among the provided mpi_communicator. Therefore, no
472 // communication is needed here.
473 const double cpu_time = sections[actual_section_name].timer.last_cpu_time();
474 sections[actual_section_name].total_cpu_time += cpu_time;
475
476 // in case we have to print out something, do that here...
479 output_is_enabled == true)
480 {
481 std::string output_time;
482 std::ostringstream cpu;
483 cpu << cpu_time << "s";
484 std::ostringstream wall;
485 wall << sections[actual_section_name].timer.last_wall_time() << "s";
486 if (output_type == cpu_times)
487 output_time = ", CPU time: " + cpu.str();
488 else if (output_type == wall_times)
489 output_time = ", wall time: " + wall.str() + ".";
490 else
491 output_time =
492 ", CPU/wall time: " + cpu.str() + " / " + wall.str() + ".";
493
494 out_stream << actual_section_name << output_time << std::endl;
495 }
496
497 // delete the index from the list of
498 // active ones
499 active_sections.erase(std::find(active_sections.begin(),
500 active_sections.end(),
501 actual_section_name));
502}
503
504
505
506std::map<std::string, double>
508{
509 std::map<std::string, double> output;
510 for (const auto &section : sections)
511 {
512 switch (kind)
513 {
515 output[section.first] = section.second.total_cpu_time;
516 break;
518 output[section.first] = section.second.total_wall_time;
519 break;
521 output[section.first] = section.second.n_calls;
522 break;
523 default:
525 }
526 }
527 return output;
528}
529
530
531
532void
534{
535 // we are going to change the precision and width of output below. store the
536 // old values so the get restored when exiting this function
537 const boost::io::ios_base_all_saver restore_stream(out_stream.get_stream());
538
539 // get the maximum width among all sections
540 unsigned int max_width = 0;
541 for (const auto &i : sections)
542 max_width = std::max(max_width, static_cast<unsigned int>(i.first.size()));
543
544 // 32 is the default width until | character
545 max_width = std::max(max_width + 1, static_cast<unsigned int>(32));
546 const std::string extra_dash = std::string(max_width - 32, '-');
547 const std::string extra_space = std::string(max_width - 32, ' ');
548
550 {
551 // in case we want to write CPU times
552 if (output_type != wall_times)
553 {
554 double total_cpu_time =
556
557 // check that the sum of all times is less or equal than the total
558 // time. otherwise, we might have generated a lot of overhead in this
559 // function.
560 double check_time = 0.;
561 for (const auto &i : sections)
562 check_time += i.second.total_cpu_time;
563
564 const double time_gap = check_time - total_cpu_time;
565 if (time_gap > 0.0)
566 total_cpu_time = check_time;
567
568 // generate a nice table
569 out_stream << "\n\n"
570 << "+---------------------------------------------"
571 << extra_dash << "+------------"
572 << "+------------+\n"
573 << "| Total CPU time elapsed since start "
574 << extra_space << "|";
575 out_stream << std::setw(10) << std::setprecision(3) << std::right;
576 out_stream << total_cpu_time << "s | |\n";
577 out_stream << "| "
578 << extra_space << "| "
579 << "| |\n";
580 out_stream << "| Section " << extra_space
581 << "| no. calls |";
582 out_stream << std::setw(10);
583 out_stream << std::setprecision(3);
584 out_stream << " CPU time "
585 << " | % of total |\n";
586 out_stream << "+---------------------------------" << extra_dash
587 << "+-----------+------------"
588 << "+------------+";
589 for (const auto &i : sections)
590 {
591 std::string name_out = i.first;
592
593 // resize the array so that it is always of the same size
594 unsigned int pos_non_space = name_out.find_first_not_of(' ');
595 name_out.erase(0, pos_non_space);
596 name_out.resize(max_width, ' ');
597 out_stream << std::endl;
598 out_stream << "| " << name_out;
599 out_stream << "| ";
600 out_stream << std::setw(9);
601 out_stream << i.second.n_calls << " |";
602 out_stream << std::setw(10);
603 out_stream << std::setprecision(3);
604 out_stream << i.second.total_cpu_time << "s |";
605 out_stream << std::setw(10);
606 if (total_cpu_time != 0)
607 {
608 // if run time was less than 0.1%, just print a zero to avoid
609 // printing silly things such as "2.45e-6%". otherwise print
610 // the actual percentage
611 const double fraction =
612 i.second.total_cpu_time / total_cpu_time;
613 if (fraction > 0.001)
614 {
615 out_stream << std::setprecision(2);
616 out_stream << fraction * 100;
617 }
618 else
619 out_stream << 0.0;
620
621 out_stream << "% |";
622 }
623 else
624 out_stream << 0.0 << "% |";
625 }
626 out_stream << std::endl
627 << "+---------------------------------" << extra_dash
628 << "+-----------+"
629 << "------------+------------+\n"
630 << std::endl;
631
632 if (time_gap > 0.0)
634 << std::endl
635 << "Note: The sum of counted times is " << time_gap
636 << " seconds larger than the total time.\n"
637 << "(Timer function may have introduced too much overhead, or different\n"
638 << "section timers may have run at the same time.)" << std::endl;
639 }
640
641 // in case we want to write out wallclock times
642 if (output_type != cpu_times)
643 {
645
646 // now generate a nice table
647 out_stream << "\n\n"
648 << "+---------------------------------------------"
649 << extra_dash << "+------------"
650 << "+------------+\n"
651 << "| Total wallclock time elapsed since start "
652 << extra_space << "|";
653 out_stream << std::setw(10) << std::setprecision(3) << std::right;
654 out_stream << total_wall_time << "s | |\n";
655 out_stream << "| "
656 << extra_space << "| "
657 << "| |\n";
658 out_stream << "| Section " << extra_space
659 << "| no. calls |";
660 out_stream << std::setw(10);
661 out_stream << std::setprecision(3);
662 out_stream << " wall time | % of total |\n";
663 out_stream << "+---------------------------------" << extra_dash
664 << "+-----------+------------"
665 << "+------------+";
666 for (const auto &i : sections)
667 {
668 std::string name_out = i.first;
669
670 // resize the array so that it is always of the same size
671 unsigned int pos_non_space = name_out.find_first_not_of(' ');
672 name_out.erase(0, pos_non_space);
673 name_out.resize(max_width, ' ');
674 out_stream << std::endl;
675 out_stream << "| " << name_out;
676 out_stream << "| ";
677 out_stream << std::setw(9);
678 out_stream << i.second.n_calls << " |";
679 out_stream << std::setw(10);
680 out_stream << std::setprecision(3);
681 out_stream << i.second.total_wall_time << "s |";
682 out_stream << std::setw(10);
683
684 if (total_wall_time != 0)
685 {
686 // if run time was less than 0.1%, just print a zero to avoid
687 // printing silly things such as "2.45e-6%". otherwise print
688 // the actual percentage
689 const double fraction =
690 i.second.total_wall_time / total_wall_time;
691 if (fraction > 0.001)
692 {
693 out_stream << std::setprecision(2);
694 out_stream << fraction * 100;
695 }
696 else
697 out_stream << 0.0;
698
699 out_stream << "% |";
700 }
701 else
702 out_stream << 0.0 << "% |";
703 }
704 out_stream << std::endl
705 << "+---------------------------------" << extra_dash
706 << "+-----------+"
707 << "------------+------------+\n"
708 << std::endl;
709 }
710 }
711 else
712 // output_type == cpu_and_wall_times_grouped
713 {
714 const double total_wall_time = timer_all.wall_time();
715 double total_cpu_time =
717
718 // check that the sum of all times is less or equal than the total time.
719 // otherwise, we might have generated a lot of overhead in this function.
720 double check_time = 0.;
721
722 for (const auto &i : sections)
723 check_time += i.second.total_cpu_time;
724
725 const double time_gap = check_time - total_cpu_time;
726 if (time_gap > 0.0)
727 total_cpu_time = check_time;
728
729 // generate a nice table
730 out_stream << "\n\n+---------------------------------------------"
731 << extra_dash << "+"
732 << "------------+------------+"
733 << "------------+------------+" << '\n'
734 << "| Total CPU/wall time elapsed since start "
735 << extra_space << "|" << std::setw(10) << std::setprecision(3)
736 << std::right << total_cpu_time << "s | |"
737 << total_wall_time << "s | |"
738 << "\n| "
739 << extra_space << "|"
740 << " | |"
741 << " | |"
742 << "\n| Section " << extra_space
743 << "| no. calls |"
744 << " CPU time | % of total |"
745 << " wall time | % of total |"
746 << "\n+---------------------------------" << extra_dash
747 << "+-----------+"
748 << "------------+------------+"
749 << "------------+------------+" << std::endl;
750
751 for (const auto &i : sections)
752 {
753 std::string name_out = i.first;
754
755 // resize the array so that it is always of the same size
756 unsigned int pos_non_space = name_out.find_first_not_of(' ');
757 name_out.erase(0, pos_non_space);
758 name_out.resize(max_width, ' ');
759 out_stream << "| " << name_out << "| ";
760
761 out_stream << std::setw(9);
762 out_stream << i.second.n_calls << " |";
763
764 if (output_type != wall_times)
765 {
766 out_stream << std::setw(10);
767 out_stream << std::setprecision(3);
768 out_stream << i.second.total_cpu_time << "s |";
769 out_stream << std::setw(10);
770 if (total_cpu_time != 0)
771 {
772 // if run time was less than 0.1%, just print a zero to avoid
773 // printing silly things such as "2.45e-6%". otherwise print
774 // the actual percentage
775 const double fraction =
776 i.second.total_cpu_time / total_cpu_time;
777 if (fraction > 0.001)
778 {
779 out_stream << std::setprecision(2);
780 out_stream << fraction * 100;
781 }
782 else
783 out_stream << 0.0;
784
785 out_stream << "% |";
786 }
787 else
788 out_stream << 0.0 << "% |";
789 }
790
791 if (output_type != cpu_times)
792 {
793 out_stream << std::setw(10);
794 out_stream << std::setprecision(3);
795 out_stream << i.second.total_wall_time << "s |";
796 out_stream << std::setw(10);
797
798 if (total_wall_time != 0)
799 {
800 // if run time was less than 0.1%, just print a zero to avoid
801 // printing silly things such as "2.45e-6%". otherwise print
802 // the actual percentage
803 const double fraction =
804 i.second.total_wall_time / total_wall_time;
805 if (fraction > 0.001)
806 {
807 out_stream << std::setprecision(2);
808 out_stream << fraction * 100;
809 }
810 else
811 out_stream << 0.0;
812
813 out_stream << "% |";
814 }
815 else
816 out_stream << 0.0 << "% |";
817 }
818 out_stream << std::endl;
819 }
820
821 out_stream << "+---------------------------------" << extra_dash
822 << "+-----------+"
823 << "------------+------------+"
824 << "------------+------------+" << std::endl
825 << std::endl;
826
827 if (output_type != wall_times && time_gap > 0.0)
829 << std::endl
830 << "Note: The sum of counted times is " << time_gap
831 << " seconds larger than the total time.\n"
832 << "(Timer function may have introduced too much overhead, or different\n"
833 << "section timers may have run at the same time.)" << std::endl;
834 }
835}
836
837
838
839void
841 const double quantile) const
842{
843 // we are going to change the precision and width of output below. store the
844 // old values so the get restored when exiting this function
845 const boost::io::ios_base_all_saver restore_stream(out_stream.get_stream());
846
848 Utilities::MPI::max(sections.size(), mpi_comm));
849 Assert(quantile >= 0. && quantile <= 0.5,
850 ExcMessage("The quantile must be between 0 and 0.5"));
851
852 // get the maximum width among all sections
853 unsigned int max_width = 0;
854 for (const auto &i : sections)
855 max_width = std::max(max_width, static_cast<unsigned int>(i.first.size()));
856
857 // 17 is the default width until | character
858 max_width = std::max(max_width + 1, static_cast<unsigned int>(17));
859 const std::string extra_dash = std::string(max_width - 17, '-');
860 const std::string extra_space = std::string(max_width - 17, ' ');
861
862 // function to print data in a nice table
863 const auto print_statistics = [&](const double given_time) {
864 const unsigned int n_ranks = Utilities::MPI::n_mpi_processes(mpi_comm);
865 if (n_ranks == 1 || quantile == 0.)
866 {
868 Utilities::MPI::min_max_avg(given_time, mpi_comm);
869
870 out_stream << std::setw(10) << std::setprecision(4) << std::right;
871 out_stream << data.min << "s ";
872 out_stream << std::setw(5) << std::right;
873 out_stream << data.min_index << (n_ranks > 99999 ? "" : " ") << "|";
874 out_stream << std::setw(10) << std::setprecision(4) << std::right;
875 out_stream << data.avg << "s |";
876 out_stream << std::setw(10) << std::setprecision(4) << std::right;
877 out_stream << data.max << "s ";
878 out_stream << std::setw(5) << std::right;
879 out_stream << data.max_index << (n_ranks > 99999 ? "" : " ") << "|\n";
880 }
881 else
882 {
883 const unsigned int my_rank = Utilities::MPI::this_mpi_process(mpi_comm);
884 std::vector<double> receive_data(my_rank == 0 ? n_ranks : 0);
885 std::vector<double> result(9);
886#ifdef DEAL_II_WITH_MPI
887 int ierr = MPI_Gather(&given_time,
888 1,
889 MPI_DOUBLE,
890 receive_data.data(),
891 1,
892 MPI_DOUBLE,
893 0,
894 mpi_comm);
895 AssertThrowMPI(ierr);
896 if (my_rank == 0)
897 {
898 // fill the received data in a pair and sort; on the way, also
899 // compute the average
900 std::vector<std::pair<double, unsigned int>> data_rank;
901 data_rank.reserve(n_ranks);
902 for (unsigned int i = 0; i < n_ranks; ++i)
903 {
904 data_rank.emplace_back(receive_data[i], i);
905 result[4] += receive_data[i];
906 }
907 result[4] /= n_ranks;
908 std::sort(data_rank.begin(), data_rank.end());
909
910 const unsigned int quantile_index =
911 static_cast<unsigned int>(std::round(quantile * n_ranks));
912 AssertIndexRange(quantile_index, data_rank.size());
913 result[0] = data_rank[0].first;
914 result[1] = data_rank[0].second;
915 result[2] = data_rank[quantile_index].first;
916 result[3] = data_rank[quantile_index].second;
917 result[5] = data_rank[n_ranks - 1 - quantile_index].first;
918 result[6] = data_rank[n_ranks - 1 - quantile_index].second;
919 result[7] = data_rank[n_ranks - 1].first;
920 result[8] = data_rank[n_ranks - 1].second;
921 }
922 ierr = MPI_Bcast(result.data(), 9, MPI_DOUBLE, 0, mpi_comm);
923 AssertThrowMPI(ierr);
924#endif
925 out_stream << std::setw(10) << std::setprecision(4) << std::right;
926 out_stream << result[0] << "s ";
927 out_stream << std::setw(5) << std::right;
928 out_stream << static_cast<unsigned int>(result[1])
929 << (n_ranks > 99999 ? "" : " ") << "|";
930 out_stream << std::setw(10) << std::setprecision(4) << std::right;
931 out_stream << result[2] << "s ";
932 out_stream << std::setw(5) << std::right;
933 out_stream << static_cast<unsigned int>(result[3])
934 << (n_ranks > 99999 ? "" : " ") << "|";
935 out_stream << std::setw(10) << std::setprecision(4) << std::right;
936 out_stream << result[4] << "s |";
937 out_stream << std::setw(10) << std::setprecision(4) << std::right;
938 out_stream << result[5] << "s ";
939 out_stream << std::setw(5) << std::right;
940 out_stream << static_cast<unsigned int>(result[6])
941 << (n_ranks > 99999 ? "" : " ") << "|";
942 out_stream << std::setw(10) << std::setprecision(4) << std::right;
943 out_stream << result[7] << "s ";
944 out_stream << std::setw(5) << std::right;
945 out_stream << static_cast<unsigned int>(result[8])
946 << (n_ranks > 99999 ? "" : " ") << "|\n";
947 }
948 };
949
950 // in case we want to write out wallclock times
951 {
952 const unsigned int n_ranks = Utilities::MPI::n_mpi_processes(mpi_comm);
953
954 const std::string time_rank_column = "------------------+";
955 const std::string time_rank_space = " |";
956
957 // now generate a nice table
958 out_stream << '\n'
959 << "+------------------------------" << extra_dash << "+"
960 << time_rank_column
961 << (n_ranks > 1 && quantile > 0. ? time_rank_column : "")
962 << "------------+"
963 << (n_ranks > 1 && quantile > 0. ? time_rank_column : "")
964 << time_rank_column << '\n'
965 << "| Total wallclock time elapsed " << extra_space << "|";
966
967 print_statistics(timer_all.wall_time());
968
969 out_stream << "| " << extra_space << "|"
970 << time_rank_space
971 << (n_ranks > 1 && quantile > 0. ? time_rank_space : "")
972 << " "
973 << (n_ranks > 1 && quantile > 0. ? time_rank_space : "")
974 << time_rank_space << '\n';
975 out_stream << "| Section " << extra_space << "| no. calls "
976 << "| min time rank |";
977 if (n_ranks > 1 && quantile > 0.)
978 out_stream << " " << std::setw(5) << std::setprecision(2) << std::right
979 << quantile << "-tile rank |";
980 out_stream << " avg time |";
981 if (n_ranks > 1 && quantile > 0.)
982 out_stream << " " << std::setw(5) << std::setprecision(2) << std::right
983 << 1. - quantile << "-tile rank |";
984 out_stream << " max time rank |\n";
985 out_stream << "+------------------------------" << extra_dash << "+"
986 << time_rank_column
987 << (n_ranks > 1 && quantile > 0. ? time_rank_column : "")
988 << "------------+"
989 << (n_ranks > 1 && quantile > 0. ? time_rank_column : "")
990 << time_rank_column << '\n';
991 for (const auto &i : sections)
992 {
993 std::string name_out = i.first;
994
995 // resize the array so that it is always of the same size
996 unsigned int pos_non_space = name_out.find_first_not_of(' ');
997 name_out.erase(0, pos_non_space);
998 name_out.resize(max_width, ' ');
999 out_stream << "| " << name_out;
1000 out_stream << "| ";
1001 out_stream << std::setw(9);
1002 out_stream << i.second.n_calls << " |";
1003
1004 print_statistics(i.second.total_wall_time);
1005 }
1006 out_stream << "+------------------------------" << extra_dash << "+"
1007 << time_rank_column
1008 << (n_ranks > 1 && quantile > 0. ? time_rank_column : "")
1009 << "------------+"
1010 << (n_ranks > 1 && quantile > 0. ? time_rank_column : "")
1011 << time_rank_column << '\n';
1012 }
1013}
1014
1015
1016
1017void
1022
1023
1024
1025void
1030
1031void
1033{
1034 std::lock_guard<std::mutex> lock(mutex);
1035 sections.clear();
1036 active_sections.clear();
1038}
1039
1041{
1042 try
1043 {
1044 stop();
1045 }
1046 catch (...)
1047 {}
1048}
1049
1050
std::ostream & get_stream() const
void reset()
Definition timer.cc:1032
MPI_Comm mpi_communicator
Definition timer.h:881
void print_summary() const
Definition timer.cc:533
@ cpu_and_wall_times_grouped
Definition timer.h:659
@ cpu_times
Definition timer.h:647
@ wall_times
Definition timer.h:651
OutputFrequency output_frequency
Definition timer.h:828
void disable_output()
Definition timer.cc:1018
std::map< std::string, Section > sections
Definition timer.h:857
void enable_output()
Definition timer.cc:1026
void leave_subsection(const std::string &section_name="")
Definition timer.cc:443
OutputType output_type
Definition timer.h:833
OutputFrequency
Definition timer.h:599
@ every_call
Definition timer.h:603
@ every_call_and_summary
Definition timer.h:611
void print_wall_time_statistics(const MPI_Comm mpi_comm, const double print_quantile=0.) const
Definition timer.cc:840
std::list< std::string > active_sections
Definition timer.h:876
Threads::Mutex mutex
Definition timer.h:887
ConditionalOStream out_stream
Definition timer.h:862
@ total_wall_time
Definition timer.h:631
@ total_cpu_time
Definition timer.h:627
std::map< std::string, double > get_summary_data(const OutputData kind) const
Definition timer.cc:507
bool output_is_enabled
Definition timer.h:868
~TimerOutput()
Definition timer.cc:349
Timer timer_all
Definition timer.h:840
void enter_subsection(const std::string &section_name)
Definition timer.cc:401
TimerOutput(std::ostream &stream, const OutputFrequency output_frequency, const OutputType output_type)
Definition timer.cc:299
Definition timer.h:117
double last_cpu_time() const
Definition timer.cc:254
bool sync_lap_times
Definition timer.h:334
void start()
Definition timer.cc:176
bool running
Definition timer.h:321
Utilities::MPI::MinMaxAvg accumulated_wall_time_data
Definition timer.h:349
double cpu_time() const
Definition timer.cc:234
Timer()
Definition timer.cc:158
double wall_time() const
Definition timer.cc:262
MPI_Comm mpi_communicator
Definition timer.h:328
ClockMeasurements< wall_clock_type > wall_times
Definition timer.h:311
ClockMeasurements< cpu_clock_type > cpu_times
Definition timer.h:316
void reset()
Definition timer.cc:286
double stop()
Definition timer.cc:193
Utilities::MPI::MinMaxAvg last_lap_wall_time_data
Definition timer.h:341
void restart()
Definition timer.h:896
double last_wall_time() const
Definition timer.cc:278
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_WARNING(desc)
Definition config.h:591
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define DEAL_II_NOT_IMPLEMENTED()
T sum(const T &t, const MPI_Comm mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
Definition mpi.cc:128
T max(const T &t, const MPI_Comm mpi_communicator)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
Definition mpi.cc:143
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm mpi_communicator)
Definition mpi.cc:102
static const unsigned int invalid_unsigned_int
Definition types.h:220
STL namespace.
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
static time_point now() noexcept
Definition timer.cc:110
std::chrono::time_point< CPUClock, duration > time_point
Definition timer.h:58
clock_type_ clock_type
Definition timer.h:256
time_point_type current_lap_start_time
Definition timer.h:272
typename clock_type::duration duration_type
Definition timer.h:266
duration_type accumulated_time
Definition timer.h:277
duration_type last_lap_time
Definition timer.h:282
unsigned int max_index
Definition mpi.h:979
unsigned int min_index
Definition mpi.h:969