From: Daniel Arndt Date: Fri, 23 Jun 2017 22:55:34 +0000 (+0200) Subject: Fixup timer X-Git-Tag: v9.0.0-rc1~1462^2~3 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=9ef67983fdda5da0e5c1f0d8eb5cb5c619d67f3c;p=dealii.git Fixup timer --- diff --git a/doc/news/changes/minor/20170623DanielArndt b/doc/news/changes/minor/20170623DanielArndt new file mode 100644 index 0000000000..e50b1a6bd2 --- /dev/null +++ b/doc/news/changes/minor/20170623DanielArndt @@ -0,0 +1,5 @@ +New: The class Timer has the new functions get_total_data(), +print_total_data, last_wall_time, cpu_time and last_cpu_time() +for accessing information about the total run and the last lap now. +
+(Daniel Arndt, 2017/06/23) diff --git a/include/deal.II/base/timer.h b/include/deal.II/base/timer.h index 44a9886ec8..64dff3e921 100644 --- a/include/deal.II/base/timer.h +++ b/include/deal.II/base/timer.h @@ -95,24 +95,37 @@ public: * time for the MPI operations are not included in the timing but may slow * down your program. * - * This constructor is only available if the deal.II compiler is an MPI - * compiler. + * This constructor is only available if deal.II is compiled with MPI support. */ Timer (MPI_Comm mpi_communicator, const bool sync_wall_time = false); /** - * Return a reference to the data structure with global timing information. - * Filled after calling stop(). + * Return a reference to the data structure with global timing information + * for the last lap. Filled after calling stop(). */ const Utilities::MPI::MinMaxAvg &get_data() const; /** - * Prints the data to the given stream. + * Return a reference to the data structure with global timing information + * for the total run. + * Filled after calling stop(). + */ + const Utilities::MPI::MinMaxAvg &get_total_data() const; + + /** + * Prints the data returned by get_data(), i.e. for the last lap, + * to the given stream. */ template void print_data(StreamType &stream) const; + /** + * Prints the data returned by get_total_data(), i.e. for the total run, + * to the given stream. + */ + template + void print_total_data(StreamType &stream) const; #endif @@ -140,20 +153,45 @@ public: void restart(); /** - * Access to the current CPU time without disturbing time measurement. The - * elapsed time is returned in units of seconds. + * Access to the current CPU time without stopping the timer. + * The elapsed time is returned in units of seconds. + * + * @deprecated Use cpu_time() instead. */ - double operator() () const; + double operator() () const DEAL_II_DEPRECATED; /** - * Access to the current wall time without disturbing time measurement. The - * elapsed time is returned in units of seconds. + * Access to the current total wall time without stopping the timer. + * The elapsed time is returned in units of seconds. */ double wall_time () const; /** - * Return the last lap time; the time taken between the last start()/stop() + * Access to the wall time for the currently running lap without stopping + * the timer. + * The elapsed time is returned in units of seconds. + */ + double last_wall_time() const; + + /** + * Access to the current total wall time without stopping the timer. + * The elapsed time is returned in units of seconds. + */ + double cpu_time() const; + + /** + * Access to the wall time for the currently running lap without stopping + * the timer. + * The elapsed time is returned in units of seconds. + */ + double last_cpu_time() const; + + /** + * Return the last lap time; the time taken between the last + * start()/stop() * call. + * + * @deprecated Use last_wall_time() instead. */ double get_lap_time () const; @@ -197,10 +235,15 @@ private: double cumulative_wall_time; /** - * Stores the last lap time; the time between the last start()/stop() cycle. + * Stores the last lap wall time; the time between the last start()/stop() cycle. */ double last_lap_time; + /** + * Stores the last lap cpu time; the time between the last start()/stop() cycle. + */ + double last_lap_cpu_time; + /** * Store whether the timer is presently running. */ @@ -221,9 +264,17 @@ private: * A structure for parallel wall time measurement that includes the minimum * time recorded among all processes, the maximum time as well as the * average time defined as the sum of all individual times divided by the - * number of MPI processes in the MPI_Comm. + * number of MPI processes in the MPI_Comm for the last lap. */ Utilities::MPI::MinMaxAvg mpi_data; + + /** + * A structure for parallel wall time measurement that includes the minimum + * time recorded among all processes, the maximum time as well as the + * average time defined as the sum of all individual times divided by the + * number of MPI processes in the MPI_Comm for the total run time. + */ + Utilities::MPI::MinMaxAvg mpi_total_data; #endif }; @@ -732,18 +783,37 @@ Timer::get_data() const +inline +const Utilities::MPI::MinMaxAvg & +Timer::get_total_data() const +{ + return mpi_total_data; +} + + + template inline void Timer::print_data(StreamType &stream) const { - unsigned int my_id = dealii::Utilities::MPI::this_mpi_process(mpi_communicator); - if (my_id==0) - stream << mpi_data.max << " wall," - << " max @" << mpi_data.max_index - << ", min=" << mpi_data.min << " @" << mpi_data.min_index - << ", avg=" << mpi_data.avg - << std::endl; + const Utilities::MPI::MinMaxAvg statistic = get_data(); + stream << statistic.max << " wall," + << " max @" << statistic.max_index << ", min=" << statistic.min << " @" + << statistic.min_index << ", avg=" << statistic.avg << std::endl; +} + + + +template +inline +void +Timer::print_total_data(StreamType &stream) const +{ + const Utilities::MPI::MinMaxAvg statistic = get_total_data(); + stream << statistic.max << " wall," + << " max @" << statistic.max_index << ", min=" << statistic.min << " @" + << statistic.min_index << ", avg=" << statistic.avg << std::endl; } #endif diff --git a/source/base/timer.cc b/source/base/timer.cc index 14bae115b5..c7441509a6 100644 --- a/source/base/timer.cc +++ b/source/base/timer.cc @@ -48,6 +48,7 @@ Timer::Timer() cumulative_time (0.), cumulative_wall_time (0.), last_lap_time (0.), + last_lap_cpu_time(0.), running (false) #ifdef DEAL_II_WITH_MPI , @@ -56,8 +57,10 @@ Timer::Timer() #endif { #ifdef DEAL_II_WITH_MPI - mpi_data.sum = mpi_data.min = mpi_data.max = mpi_data.avg = numbers::signaling_nan(); - mpi_data.min_index = mpi_data.max_index = numbers::invalid_unsigned_int; + mpi_data.sum = mpi_data.min = mpi_data.max = mpi_data.avg = 0.; + mpi_data.min_index = mpi_data.max_index = 0; + mpi_total_data.sum = mpi_total_data.min = mpi_total_data.max = mpi_total_data.avg = 0.; + mpi_total_data.min_index = mpi_total_data.max_index = 0; #endif start(); @@ -77,13 +80,16 @@ Timer::Timer(MPI_Comm mpi_communicator, cumulative_time (0.), cumulative_wall_time (0.), last_lap_time (0.), + last_lap_cpu_time (0.), running (false), mpi_communicator (mpi_communicator), sync_wall_time(sync_wall_time_) { #ifdef DEAL_II_WITH_MPI - mpi_data.sum = mpi_data.min = mpi_data.max = mpi_data.avg = numbers::signaling_nan(); - mpi_data.min_index = mpi_data.max_index = numbers::invalid_unsigned_int; + mpi_data.sum = mpi_data.min = mpi_data.max = mpi_data.avg = 0.; + mpi_data.min_index = mpi_data.max_index = 0; + mpi_total_data.sum = mpi_total_data.min = mpi_total_data.max = mpi_total_data.avg = 0.; + mpi_total_data.min_index = mpi_total_data.max_index = 0; #endif start(); @@ -175,13 +181,13 @@ double Timer::stop () rusage usage; getrusage (RUSAGE_SELF, &usage); const double dtime = usage.ru_utime.tv_sec + 1.e-6 * usage.ru_utime.tv_usec; - cumulative_time += dtime - start_time; + last_lap_cpu_time = dtime - start_time; rusage usage_children; getrusage (RUSAGE_CHILDREN, &usage_children); const double dtime_children = usage_children.ru_utime.tv_sec + 1.e-6 * usage_children.ru_utime.tv_usec; - cumulative_time += dtime_children - start_time_children; + last_lap_cpu_time += dtime_children - start_time_children; struct timeval wall_timer; gettimeofday(&wall_timer, nullptr); @@ -189,37 +195,34 @@ double Timer::stop () - start_wall_time; #elif defined(DEAL_II_MSVC) last_lap_time = windows::wall_clock() - start_wall_time; - cumulative_time += windows::cpu_clock() - start_time; + last_lap_cpu_time = windows::cpu_clock() - start_time; #else # error Unsupported platform. Porting not finished. #endif #ifdef DEAL_II_WITH_MPI + this->mpi_data = Utilities::MPI::min_max_avg (last_lap_time, + mpi_communicator); if (sync_wall_time && Utilities::MPI::job_supports_mpi()) { - this->mpi_data - = Utilities::MPI::min_max_avg (last_lap_time, mpi_communicator); last_lap_time = this->mpi_data.max; - cumulative_wall_time += last_lap_time; + last_lap_cpu_time = Utilities::MPI::min_max_avg (last_lap_cpu_time, + mpi_communicator).max; } - else + cumulative_wall_time += last_lap_time; + cumulative_time += last_lap_cpu_time; + this->mpi_total_data = Utilities::MPI::min_max_avg (cumulative_wall_time, + mpi_communicator); #endif - cumulative_wall_time += last_lap_time; + cumulative_wall_time += last_lap_time; + cumulative_time += last_lap_cpu_time; } return cumulative_time; } -double Timer::get_lap_time() const -{ - // time already has the difference between the last start()/stop() cycle. - return Utilities::MPI::max (last_lap_time, mpi_communicator); -} - - - -double Timer::operator() () const +double Timer::cpu_time() const { if (running) { @@ -257,6 +260,27 @@ double Timer::operator() () const +double Timer::last_cpu_time() const +{ + return last_lap_cpu_time; +} + + + +double Timer::get_lap_time() const +{ + return last_lap_time; +} + + + +double Timer::operator() () const +{ + return cpu_time(); +} + + + double Timer::wall_time () const { if (running) @@ -279,12 +303,26 @@ double Timer::wall_time () const +double Timer::last_wall_time () const +{ + return last_lap_time; +} + + + void Timer::reset () { last_lap_time = 0.; + last_lap_cpu_time = 0.; cumulative_time = 0.; cumulative_wall_time = 0.; running = false; +#ifdef DEAL_II_WITH_MPI + mpi_data.sum = mpi_data.min = mpi_data.max = mpi_data.avg = 0.; + mpi_data.min_index = mpi_data.max_index = 0; + mpi_total_data.sum = mpi_total_data.min = mpi_total_data.max = mpi_total_data.avg = 0.; + mpi_total_data.min_index = mpi_total_data.max_index = 0; +#endif } diff --git a/tests/base/timer_04.cc b/tests/base/timer_04.cc new file mode 100644 index 0000000000..8fc498bfe7 --- /dev/null +++ b/tests/base/timer_04.cc @@ -0,0 +1,72 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2017 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + +// the same as timer.cc but also test the new functions last_wall_time(), +// cpu_time() and last_cpu_time(). + +#include "../tests.h" +#include +#include +#include +#include +#include + +// burn computer time + +double s = 0.; +void burn (unsigned int n) +{ + for (unsigned int i=0 ; i 0., ExcInternalError()); + const double old_cpu_time = t.wall_time(); + AssertThrow(old_cpu_time > 0., ExcInternalError()); + AssertThrow(t.last_wall_time() == 0., ExcInternalError()); + AssertThrow(t.last_cpu_time() == 0, ExcInternalError()); + + burn(50); + AssertThrow(t.stop() > 0., ExcInternalError()); + + AssertThrow(t.wall_time() > old_wall_time, ExcInternalError()); + AssertThrow(t.cpu_time() > old_cpu_time, ExcInternalError()); + AssertThrow(t.last_wall_time() > 0., ExcInternalError()); + AssertThrow(t.last_cpu_time() > 0, ExcInternalError()); + + t.reset(); + AssertThrow(t.wall_time() == 0., ExcInternalError()); + AssertThrow(t.cpu_time()== 0., ExcInternalError()); + AssertThrow(t.last_wall_time() == 0., ExcInternalError()); + AssertThrow(t.last_cpu_time() == 0, ExcInternalError()); + + deallog << "OK" << std::endl; +} + diff --git a/tests/base/timer_04.output b/tests/base/timer_04.output new file mode 100644 index 0000000000..0fd8fc12f0 --- /dev/null +++ b/tests/base/timer_04.output @@ -0,0 +1,2 @@ + +DEAL::OK diff --git a/tests/base/timer_05.cc b/tests/base/timer_05.cc new file mode 100644 index 0000000000..31a6e4a6e2 --- /dev/null +++ b/tests/base/timer_05.cc @@ -0,0 +1,90 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2017 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + +// the same as timer_04.cc but this time with MPI. + +#include "../tests.h" +#include +#include +#include +#include +#include + +// burn computer time + +double s = 0.; +void burn (unsigned int n) +{ + for (unsigned int i=0 ; i 0., ExcInternalError()); + const double old_cpu_time = t.wall_time(); + AssertThrow(old_cpu_time > 0., ExcInternalError()); + AssertThrow(t.last_wall_time() == 0., ExcInternalError()); + AssertThrow(t.last_cpu_time() == 0., ExcInternalError()); + AssertThrow(t.get_data().min == 0., ExcInternalError()); + AssertThrow(t.get_total_data().min == 0., ExcInternalError()); + + burn(50); + AssertThrow(t.stop() > 0., ExcInternalError()); + + AssertThrow(t.wall_time() > old_wall_time, ExcInternalError()); + AssertThrow(t.cpu_time() > old_cpu_time, ExcInternalError()); + AssertThrow(t.last_wall_time() > 0., ExcInternalError()); + AssertThrow(t.last_cpu_time() > 0., ExcInternalError()); + AssertThrow(t.get_data().min > 0., ExcInternalError()); + AssertThrow(t.get_total_data().min > 0., ExcInternalError()); + + t.reset(); + AssertThrow(t.wall_time() == 0., ExcInternalError()); + AssertThrow(t.cpu_time()== 0., ExcInternalError()); + AssertThrow(t.last_wall_time() == 0., ExcInternalError()); + AssertThrow(t.last_cpu_time() == 0., ExcInternalError()); + AssertThrow(t.get_data().min == 0., ExcInternalError()); + AssertThrow(t.get_total_data().min == 0, ExcInternalError()); + + deallog << "OK" << std::endl; +} + + + +int main (int argc, char **argv) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); + mpi_initlog(); + deallog.threshold_double(1.e-10); + + Timer t1(MPI_COMM_WORLD, false); + Timer t2(MPI_COMM_WORLD, true); + + test_timer(t1); + test_timer(t2); +} + diff --git a/tests/base/timer_05.with_mpi=true.mpirun=1.output b/tests/base/timer_05.with_mpi=true.mpirun=1.output new file mode 100644 index 0000000000..8b3b075900 --- /dev/null +++ b/tests/base/timer_05.with_mpi=true.mpirun=1.output @@ -0,0 +1,3 @@ + +DEAL::OK +DEAL::OK