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> 30 #include <type_traits> 32 #ifdef DEAL_II_HAVE_SYS_RESOURCE_H 33 # include <sys/resource.h> 42 DEAL_II_NAMESPACE_OPEN
46 namespace TimerImplementation
54 struct is_duration : std::false_type {};
59 template <
typename Rep,
typename Period>
60 struct is_duration<
std::chrono::duration<Rep, Period>> : std::true_type {};
69 from_seconds(
const double time)
71 static_assert(is_duration<T>::value,
72 "The template type should be a duration type.");
73 return T(std::lround(T::period::den*(time/T::period::num)));
80 template <
typename Rep,
typename Period>
82 to_seconds(
const std::chrono::duration<Rep, Period> duration)
84 return Period::num*double(duration.count())/Period::den;
93 data.
sum = numbers::signaling_nan<double>();
94 data.
min = numbers::signaling_nan<double>();
95 data.
max = numbers::signaling_nan<double>();
96 data.
avg = numbers::signaling_nan<double>();
108 double system_cpu_duration = 0.0;
110 FILETIME cpuTime, sysTime, createTime, exitTime;
111 const auto succeeded = GetProcessTimes
112 (GetCurrentProcess(), &createTime, &exitTime, &sysTime, &cpuTime);
115 system_cpu_duration = (double)
116 (((
unsigned long long)cpuTime.dwHighDateTime << 32)
117 | cpuTime.dwLowDateTime) / 1e7;
120 #elif defined(DEAL_II_HAVE_SYS_RESOURCE_H) 122 getrusage (RUSAGE_SELF, &usage);
123 system_cpu_duration = usage.ru_utime.tv_sec + 1.e-6 * usage.ru_utime.tv_usec;
125 # warning "Unsupported platform. Porting not finished." 127 return time_point(internal::TimerImplementation::from_seconds<duration>(system_cpu_duration));
132 template <
typename clock_type_>
142 template <
typename clock_type_>
146 current_lap_start_time = clock_type::now();
147 accumulated_time = duration_type::zero();
148 last_lap_time = duration_type::zero();
155 Timer(MPI_COMM_SELF, false)
161 const bool sync_lap_times_)
164 mpi_communicator (mpi_communicator),
165 sync_lap_times(sync_lap_times_)
176 #ifdef DEAL_II_WITH_MPI 205 cpu_times.last_lap_time = internal::TimerImplementation::from_seconds<decltype(cpu_times)::duration_type>
207 (internal::TimerImplementation::to_seconds(
cpu_times.last_lap_time),
216 return internal::TimerImplementation::to_seconds(
cpu_times.accumulated_time);
225 const double running_time = internal::TimerImplementation::to_seconds
242 return internal::TimerImplementation::to_seconds(
cpu_times.last_lap_time);
263 wall_clock_type::duration current_elapsed_wall_time;
265 current_elapsed_wall_time = wall_clock_type::now()
271 return internal::TimerImplementation::to_seconds(current_elapsed_wall_time);
300 output_frequency (output_frequency),
301 output_type (output_type),
302 out_stream (stream, true),
303 output_is_enabled (true),
304 mpi_communicator (MPI_COMM_SELF)
313 output_frequency (output_frequency),
314 output_type (output_type),
316 output_is_enabled (true),
317 mpi_communicator (MPI_COMM_SELF)
323 std::ostream &stream,
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)
341 output_frequency (output_frequency),
342 output_type (output_type),
344 output_is_enabled (true),
345 mpi_communicator (mpi_communicator)
352 auto do_exit = [
this]()
370 #ifdef DEAL_II_WITH_MPI 373 std::cerr <<
"---------------------------------------------------------" 375 <<
"TimerOutput objects finalize timed values printed to the" 377 <<
"screen by communicating over MPI in their destructors." 379 <<
"Since an exception is currently uncaught, this" 381 <<
"synchronization (and subsequent output) will be skipped to" 383 <<
"avoid a possible deadlock." 385 <<
"---------------------------------------------------------" 404 Assert (section_name.empty() ==
false,
409 ExcMessage (std::string(
"Cannot enter the already active section <")
410 + section_name +
">."));
427 sections[section_name].total_cpu_time = 0;
428 sections[section_name].total_wall_time = 0;
432 sections[section_name].timer.reset();
433 sections[section_name].timer.start();
445 ExcMessage(
"Cannot exit any section because none has been entered!"));
449 if (section_name !=
"")
452 ExcMessage (
"Cannot delete a section that was never created."));
455 ExcMessage (
"Cannot delete a section that has not been entered."));
460 const std::string actual_section_name = (section_name ==
"" ?
464 sections[actual_section_name].timer.stop();
465 sections[actual_section_name].total_wall_time
466 +=
sections[actual_section_name].timer.last_wall_time();
472 const double cpu_time =
sections[actual_section_name].timer.last_cpu_time();
473 sections[actual_section_name].total_cpu_time += cpu_time;
479 std::string output_time;
480 std::ostringstream cpu;
481 cpu << cpu_time <<
"s";
482 std::ostringstream wall;
483 wall <<
sections[actual_section_name].timer.last_wall_time() <<
"s";
485 output_time =
", CPU time: " + cpu.str();
487 output_time =
", wall time: " + wall.str() +
".";
489 output_time =
", CPU/wall time: " + cpu.str() +
" / " + wall.str() +
".";
491 out_stream << actual_section_name << output_time
498 actual_section_name));
503 std::map<std::string, double>
506 std::map<std::string, double> output;
507 for (
const auto §ion :
sections)
511 case TimerOutput::OutputData::total_cpu_time:
512 output[section.first] = section.second.total_cpu_time;
514 case TimerOutput::OutputData::total_wall_time:
515 output[section.first] = section.second.total_wall_time;
517 case TimerOutput::OutputData::n_calls:
518 output[section.first] = section.second.n_calls;
550 double check_time = 0.;
551 for (std::map<std::string, Section>::const_iterator
553 check_time += i->second.total_cpu_time;
561 <<
"+---------------------------------------------+------------" 562 <<
"+------------+\n" 563 <<
"| Total CPU time elapsed since start |";
564 out_stream << std::setw(10) << std::setprecision(3) << std::right;
571 out_stream <<
" CPU time " <<
" | % of total |\n";
572 out_stream <<
"+---------------------------------+-----------+------------" 574 for (std::map<std::string, Section>::const_iterator
577 std::string name_out = i->first;
581 unsigned int pos_non_space = name_out.find_first_not_of (
' ');
582 name_out.erase(0, pos_non_space);
583 name_out.resize (32,
' ');
591 out_stream << i->second.total_cpu_time <<
"s |";
599 if (fraction > 0.001)
613 <<
"+---------------------------------+-----------+" 614 <<
"------------+------------+\n" 619 <<
"Note: The sum of counted times is " << time_gap
620 <<
" seconds larger than the total time.\n" 621 <<
"(Timer function may have introduced too much overhead, or different\n" 622 <<
"section timers may have run at the same time.)" << std::endl;
632 <<
"+---------------------------------------------+------------" 633 <<
"+------------+\n" 634 <<
"| Total wallclock time elapsed since start |";
635 out_stream << std::setw(10) << std::setprecision(3) << std::right;
643 out_stream <<
"+---------------------------------+-----------+------------" 645 for (std::map<std::string, Section>::const_iterator
648 std::string name_out = i->first;
652 unsigned int pos_non_space = name_out.find_first_not_of (
' ');
653 name_out.erase(0, pos_non_space);
654 name_out.resize (32,
' ');
662 out_stream << i->second.total_wall_time <<
"s |";
671 if (fraction > 0.001)
685 <<
"+---------------------------------+-----------+" 686 <<
"------------+------------+\n" 732 DEAL_II_NAMESPACE_CLOSE
void print_summary() const
static const unsigned int invalid_unsigned_int
double last_wall_time() const
double get_lap_time() const
MPI_Comm mpi_communicator
MPI_Comm mpi_communicator
static time_point now() noexcept
clock_type::duration duration_type
Utilities::MPI::MinMaxAvg last_lap_wall_time_data
static ::ExceptionBase & ExcMessage(std::string arg1)
T sum(const T &t, const MPI_Comm &mpi_communicator)
#define Assert(cond, exc)
double operator()() const
duration_type accumulated_time
std::chrono::time_point< CPUClock, duration > time_point
ClockMeasurements< wall_clock_type > wall_times
duration_type last_lap_time
std::list< std::string > active_sections
std::map< std::string, double > get_summary_data(const OutputData kind) const
#define AssertThrowMPI(error_code)
double last_cpu_time() const
TimerOutput(std::ostream &stream, const OutputFrequency output_frequency, const OutputType output_type)
std::map< std::string, Section > sections
ConditionalOStream out_stream
std::ostream & get_stream() const
static ::ExceptionBase & ExcNotImplemented()
void enter_subsection(const std::string §ion_name)
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm &mpi_communicator)
OutputFrequency output_frequency
time_point_type current_lap_start_time
Utilities::MPI::MinMaxAvg accumulated_wall_time_data
ClockMeasurements< cpu_clock_type > cpu_times
void leave_subsection(const std::string §ion_name=std::string())