]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Clean up some variable names in the Timer class. 5032/head
authorDavid Wells <wellsd2@rpi.edu>
Tue, 5 Sep 2017 21:36:08 +0000 (17:36 -0400)
committerDavid Wells <wellsd2@rpi.edu>
Tue, 5 Sep 2017 22:29:58 +0000 (18:29 -0400)
This is compatible with 8.5 but not with some previous patches to this class
made since then.

include/deal.II/base/timer.h
source/base/timer.cc
tests/base/timer_05.cc

index 7696d284b35da780747ce46eb37c5760e01280f5..3ded786c7f3b361fb7e83abcfdcab9e6c05b25e7 100644 (file)
@@ -121,7 +121,7 @@ public:
 
   /**
    * Constructor specifying that CPU times should be summed over the given
-   * communicator. If @p sync_wall_time is <code>true</code> then the Timer
+   * communicator. If @p sync_lap_times is <code>true</code> then the Timer
    * will set the elapsed wall and CPU times over the last lap to their
    * maximum values across the provided communicator. This synchronization is
    * only performed if Timer::stop() is called before the timer is queried for
@@ -134,15 +134,15 @@ public:
    * measured.
    */
   Timer (MPI_Comm mpi_communicator,
-         const bool sync_wall_time = false);
+         const bool sync_lap_times = false);
 
   /**
    * Return a reference to the data structure with global timing information
    * for the last lap. This structure does not contain meaningful values until
    * Timer::stop() has been called.
    *
-   * @deprecated Use Timer::get_last_lap_data() instead, which returns a
-   * reference to the same structure.
+   * @deprecated Use Timer::get_last_lap_wall_time_data() instead, which
+   * returns a reference to the same structure.
    */
   const Utilities::MPI::MinMaxAvg &get_data() const DEAL_II_DEPRECATED;
 
@@ -152,7 +152,7 @@ public:
    * communicator. This structure does not contain meaningful values until
    * Timer::stop() has been called.
    */
-  const Utilities::MPI::MinMaxAvg &get_last_lap_data() const;
+  const Utilities::MPI::MinMaxAvg &get_last_lap_wall_time_data() const;
 
   /**
    * Return a reference to the data structure containing basic statistics on
@@ -177,18 +177,18 @@ public:
    * Prints the data returned by get_data(), i.e. for the last lap,
    * to the given stream.
    *
-   * @deprecated Use Timer::print_last_lap_data() instead, which prints the
-   * same information.
+   * @deprecated Use Timer::print_last_lap_wall_time_data() instead, which
+   * prints the same information.
    */
   template <class StreamType>
   void print_data(StreamType &stream) const DEAL_II_DEPRECATED;
 
   /**
-   * Print the data returned by Timer::get_last_lap_data() to the given
-   * stream.
+   * Print the data returned by Timer::get_last_lap_wall_time_data() to the
+   * given stream.
    */
   template <class StreamType>
-  void print_last_lap_data(StreamType &stream) const;
+  void print_last_lap_wall_time_data(StreamType &stream) const;
 
   /**
    * Prints the data returned by get_total_data(), i.e. for the total run,
@@ -208,14 +208,18 @@ public:
   void print_accumulated_wall_time_data(StreamType &stream) const;
 
   /**
-   * Begin measuring a new lap. If <code>sync_wall_time</code> is
+   * Begin measuring a new lap. If <code>sync_lap_times</code> is
    * <code>true</code> then an MPI barrier is used to ensure that all
    * processes begin the lap at the same wall time.
    */
   void start ();
 
   /**
-   * Stop the timer. This updates the lap times and accumulated times.
+   * Stop the timer. This updates the lap times and accumulated times. If
+   * <code>sync_lap_times</code> is <code>true</code> then the lap times are
+   * synchronized over all processors in the communicator (i.e., the lap times
+   * are set to the maximum lap time).
+   *
    * Returns the accumulated CPU time in seconds.
    */
   double stop ();
@@ -232,8 +236,8 @@ public:
   void restart();
 
   /**
-   * Access to the current CPU time without stopping the timer.
-   * 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.
    */
@@ -321,8 +325,8 @@ private:
     duration_type last_lap_time;
 
     /**
-     * Constructor. Sets <code>current_lap_start_time</code> to the
-     * current clock time and the durations to zero.
+     * Constructor. Sets <code>current_lap_start_time</code> to the current
+     * clock time and the durations to zero.
      */
     ClockMeasurements();
 
@@ -356,27 +360,27 @@ private:
   /**
    * Whether or not the timer is presently running.
    */
-  bool                running;
+  bool running;
 
   /**
    * The communicator over which various time values are synchronized and
    * combined: see the documentation of the relevant constructor for
    * additional information.
    */
-  MPI_Comm            mpi_communicator;
+  MPI_Comm mpi_communicator;
 
   /**
-   * Store whether or not the wall time and CPU time are synchronized in
-   * Timer::start() and Timer::stop().
+   * Store whether or not the wall time and CPU time are synchronized across
+   * the communicator in Timer::start() and Timer::stop().
    */
-  bool sync_wall_time;
+  bool sync_lap_times;
 
   /**
    * A structure for parallel wall time measurement that includes the minimum,
    * maximum, and average over all processors known to the MPI communicator of
    * the last lap time.
    */
-  Utilities::MPI::MinMaxAvg last_lap_data;
+  Utilities::MPI::MinMaxAvg last_lap_wall_time_data;
 
   /**
    * A structure for parallel wall time measurement that includes the minimum
@@ -879,16 +883,16 @@ inline
 const Utilities::MPI::MinMaxAvg &
 Timer::get_data() const
 {
-  return last_lap_data;
+  return last_lap_wall_time_data;
 }
 
 
 
 inline
 const Utilities::MPI::MinMaxAvg &
-Timer::get_last_lap_data() const
+Timer::get_last_lap_wall_time_data() const
 {
-  return last_lap_data;
+  return last_lap_wall_time_data;
 }
 
 
@@ -916,7 +920,7 @@ inline
 void
 Timer::print_data(StreamType &stream) const
 {
-  print_last_lap_data(stream);
+  print_last_lap_wall_time_data(stream);
 }
 
 
@@ -924,9 +928,9 @@ Timer::print_data(StreamType &stream) const
 template <class StreamType>
 inline
 void
-Timer::print_last_lap_data(StreamType &stream) const
+Timer::print_last_lap_wall_time_data(StreamType &stream) const
 {
-  const Utilities::MPI::MinMaxAvg statistic = get_last_lap_data();
+  const Utilities::MPI::MinMaxAvg statistic = get_last_lap_wall_time_data();
   stream << statistic.max << " wall,"
          << " max @" << statistic.max_index << ", min=" << statistic.min << " @"
          << statistic.min_index << ", avg=" << statistic.avg << std::endl;
index ab8825498d3382362d46f3af6be04974fa6cf79c..dea34f34920512f45b0e8b9122069417a6429592 100644 (file)
@@ -152,17 +152,17 @@ Timer::ClockMeasurements<clock_type_>::reset()
 
 Timer::Timer()
   :
-  Timer(MPI_COMM_SELF, /*sync_wall_time=*/false)
+  Timer(MPI_COMM_SELF, /*sync_lap_times=*/false)
 {}
 
 
 
 Timer::Timer(MPI_Comm mpi_communicator,
-             const bool sync_wall_time_)
+             const bool sync_lap_times_)
   :
   running (false),
   mpi_communicator (mpi_communicator),
-  sync_wall_time(sync_wall_time_)
+  sync_lap_times(sync_lap_times_)
 {
   reset();
   start();
@@ -174,7 +174,7 @@ void Timer::start ()
 {
   running = true;
 #ifdef DEAL_II_WITH_MPI
-  if (sync_wall_time)
+  if (sync_lap_times)
     {
       const int ierr = MPI_Barrier(mpi_communicator);
       AssertThrowMPI(ierr);
@@ -195,13 +195,13 @@ double Timer::stop ()
       wall_times.last_lap_time = wall_clock_type::now() - wall_times.current_lap_start_time;
       cpu_times.last_lap_time = cpu_clock_type::now() - cpu_times.current_lap_start_time;
 
-      last_lap_data = Utilities::MPI::min_max_avg
-                      (internal::Timer::to_seconds(wall_times.last_lap_time),
-                       mpi_communicator);
-      if (sync_wall_time)
+      last_lap_wall_time_data = Utilities::MPI::min_max_avg
+                                (internal::Timer::to_seconds(wall_times.last_lap_time),
+                                 mpi_communicator);
+      if (sync_lap_times)
         {
           wall_times.last_lap_time = internal::Timer::from_seconds<decltype(wall_times)::duration_type>
-                                     (last_lap_data.max);
+                                     (last_lap_wall_time_data.max);
           cpu_times.last_lap_time = internal::Timer::from_seconds<decltype(cpu_times)::duration_type>
                                     (Utilities::MPI::min_max_avg
                                      (internal::Timer::to_seconds(cpu_times.last_lap_time),
@@ -285,7 +285,7 @@ void Timer::reset ()
   wall_times.reset();
   cpu_times.reset();
   running = false;
-  internal::Timer::clear_timing_data(last_lap_data);
+  internal::Timer::clear_timing_data(last_lap_wall_time_data);
   internal::Timer::clear_timing_data(accumulated_wall_time_data);
 }
 
index 493640bd2c811346a1287d523e3a49eec2ab7b03..9b50ae0ab9bc30f96adc2bcdb8f087638c69bcf1 100644 (file)
@@ -72,14 +72,14 @@ test_timer(Timer &t)
   AssertThrow(old_cpu_time > 0., ExcInternalError());
   assert_min_max_avg_invalid(t.get_data());
   assert_min_max_avg_invalid(t.get_total_data());
-  assert_min_max_avg_invalid(t.get_last_lap_data());
+  assert_min_max_avg_invalid(t.get_last_lap_wall_time_data());
   assert_min_max_avg_invalid(t.get_accumulated_wall_time_data());
 
   burn(50);
   AssertThrow(t.stop() > 0., ExcInternalError());
   assert_min_max_avg_valid(t.get_data());
   assert_min_max_avg_valid(t.get_total_data());
-  assert_min_max_avg_valid(t.get_last_lap_data());
+  assert_min_max_avg_valid(t.get_last_lap_wall_time_data());
   assert_min_max_avg_valid(t.get_accumulated_wall_time_data());
 
   AssertThrow(t.wall_time() > old_wall_time, ExcInternalError());
@@ -92,7 +92,7 @@ test_timer(Timer &t)
   AssertThrow(t.cpu_time()== 0., ExcInternalError());
   assert_min_max_avg_invalid(t.get_data());
   assert_min_max_avg_invalid(t.get_total_data());
-  assert_min_max_avg_invalid(t.get_last_lap_data());
+  assert_min_max_avg_invalid(t.get_last_lap_wall_time_data());
   assert_min_max_avg_invalid(t.get_accumulated_wall_time_data());
 
   deallog << "OK" << std::endl;

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.