Merge branch release-5-1
[gromacs.git] / src / gromacs / timing / walltime_accounting.cpp
blob41317eb5a1a381b4aff925f0462d3304e97952e8
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013, The GROMACS development team.
5 * Copyright (c) 2013,2014,2015, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
36 #include "gmxpre.h"
38 #include "walltime_accounting.h"
40 #include "config.h"
42 #include <ctime>
44 #ifdef HAVE_UNISTD_H
45 #include <unistd.h>
46 #endif
47 #ifdef HAVE_SYS_TIME_H
48 #include <sys/time.h>
49 #endif
51 #include "gromacs/utility/basedefinitions.h"
52 #include "gromacs/utility/smalloc.h"
54 /* TODO in future: convert gmx_walltime_accounting to a class,
55 * resolve who should have responsibility for recording the number of
56 * steps done, consider whether parts of finish_time, print_perf,
57 * wallcycle_print belong in this module.
59 * If/when any kind of task parallelism is implemented (even OpenMP
60 * regions simultaneously assigned to different tasks), consider
61 * whether this data structure (and/or cycle counters) should be
62 * maintained on a per-OpenMP-thread basis. */
64 /*! \brief Manages caching wall-clock time measurements for
65 * simulations */
66 typedef struct gmx_walltime_accounting {
67 //! Seconds since the epoch recorded at the start of the simulation
68 double start_time_stamp;
69 //! Seconds since the epoch recorded at the start of the simulation for this thread
70 double start_time_stamp_per_thread;
71 //! Total seconds elapsed over the simulation
72 double elapsed_time;
73 //! Total seconds elapsed over the simulation running this thread
74 double elapsed_time_over_all_threads;
75 /*! \brief Number of OpenMP threads that will be launched by this
76 * MPI rank.
78 * This is used to scale elapsed_time_over_all_threads so
79 * that any combination of real MPI, thread MPI and OpenMP (even
80 * mdrun -ntomp_pme) processes/threads would (when run at maximum
81 * efficiency) return values such that the sum of
82 * elapsed_time_over_all_threads over all threads was constant
83 * with respect to parallelism implementation. */
84 int numOpenMPThreads;
85 //! Set by integrators to report the amount of work they did
86 gmx_int64_t nsteps_done;
87 } t_gmx_walltime_accounting;
89 /*! \brief Calls system timing routines (e.g. clock_gettime) to get
90 * the (fractional) number of seconds elapsed since the epoch when
91 * this thread was executing.
93 * This can be used to measure system load. This can be unreliable if
94 * threads migrate between sockets. If thread-specific timers are not
95 * supported by the OS (e.g. if the OS is not POSIX-compliant), this
96 * function is implemented by gmx_gettime. */
97 static double gmx_gettime_per_thread();
99 // TODO In principle, all this should get protected by checks that
100 // walltime_accounting is not null. In practice, that NULL condition
101 // does not happen, and future refactoring will likely enforce it by
102 // having the gmx_walltime_accounting_t object be owned by the runner
103 // object. When these become member functions, existence will be
104 // guaranteed.
106 gmx_walltime_accounting_t
107 walltime_accounting_init(int numOpenMPThreads)
109 gmx_walltime_accounting_t walltime_accounting;
111 snew(walltime_accounting, 1);
112 walltime_accounting->start_time_stamp = 0;
113 walltime_accounting->start_time_stamp_per_thread = 0;
114 walltime_accounting->elapsed_time = 0;
115 walltime_accounting->nsteps_done = 0;
116 walltime_accounting->numOpenMPThreads = numOpenMPThreads;
118 return walltime_accounting;
121 void
122 walltime_accounting_destroy(gmx_walltime_accounting_t walltime_accounting)
124 sfree(walltime_accounting);
127 void
128 walltime_accounting_start(gmx_walltime_accounting_t walltime_accounting)
130 walltime_accounting->start_time_stamp = gmx_gettime();
131 walltime_accounting->start_time_stamp_per_thread = gmx_gettime_per_thread();
132 walltime_accounting->elapsed_time = 0;
133 walltime_accounting->nsteps_done = 0;
136 void
137 walltime_accounting_end(gmx_walltime_accounting_t walltime_accounting)
139 double now, now_per_thread;
141 now = gmx_gettime();
142 now_per_thread = gmx_gettime_per_thread();
144 walltime_accounting->elapsed_time = now - walltime_accounting->start_time_stamp;
145 walltime_accounting->elapsed_time_over_all_threads = now_per_thread - walltime_accounting->start_time_stamp_per_thread;
146 /* For thread-MPI, the per-thread CPU timer makes this just
147 * work. For OpenMP threads, the per-thread CPU timer measurement
148 * needs to be multiplied by the number of OpenMP threads used,
149 * under the current assumption that all regions ever opened
150 * within a process are of the same size, and each thread should
151 * keep one core busy.
153 walltime_accounting->elapsed_time_over_all_threads *= walltime_accounting->numOpenMPThreads;
156 double
157 walltime_accounting_get_current_elapsed_time(gmx_walltime_accounting_t walltime_accounting)
159 return gmx_gettime() - walltime_accounting->start_time_stamp;
162 double
163 walltime_accounting_get_elapsed_time(gmx_walltime_accounting_t walltime_accounting)
165 return walltime_accounting->elapsed_time;
168 double
169 walltime_accounting_get_elapsed_time_over_all_threads(gmx_walltime_accounting_t walltime_accounting)
171 return walltime_accounting->elapsed_time_over_all_threads;
174 double
175 walltime_accounting_get_start_time_stamp(gmx_walltime_accounting_t walltime_accounting)
177 return walltime_accounting->start_time_stamp;
180 gmx_int64_t
181 walltime_accounting_get_nsteps_done(gmx_walltime_accounting_t walltime_accounting)
183 return walltime_accounting->nsteps_done;
186 void
187 walltime_accounting_set_nsteps_done(gmx_walltime_accounting_t walltime_accounting,
188 gmx_int64_t nsteps_done)
190 walltime_accounting->nsteps_done = nsteps_done;
193 double
194 gmx_gettime()
196 #if defined HAVE_CLOCK_GETTIME && _POSIX_TIMERS >= 0 && !(defined __bgq__ && defined __clang__)
197 /* Mac and Windows do not support this. For added fun, Windows
198 * defines _POSIX_TIMERS without actually providing the
199 * implementation. The BlueGene/Q CNK only supports gettimeofday,
200 * and bgclang doesn't provide a fully functional implementation
201 * for clock_gettime (unlike xlc). */
202 struct timespec t;
203 double seconds;
205 clock_gettime(CLOCK_REALTIME, &t);
206 seconds = static_cast<double>(t.tv_sec) + 1e-9*t.tv_nsec;
208 return seconds;
209 #elif defined HAVE_GETTIMEOFDAY
210 // Note that gettimeofday() is deprecated by POSIX, but since Mac
211 // and Windows do not yet support POSIX, we are still stuck.
212 // Also, this is the only supported API call on Bluegene/Q.
213 struct timeval t;
214 double seconds;
216 gettimeofday(&t, NULL);
217 seconds = static_cast<double>(t.tv_sec) + 1e-6*t.tv_usec;
219 return seconds;
220 #else
221 double seconds;
223 seconds = time(NULL);
225 return seconds;
226 #endif
229 static double
230 gmx_gettime_per_thread()
232 #if defined HAVE_CLOCK_GETTIME && _POSIX_THREAD_CPUTIME >= 0
233 struct timespec t;
234 double seconds;
236 clock_gettime(CLOCK_THREAD_CPUTIME_ID, &t);
237 seconds = static_cast<double>(t.tv_sec) + 1e-9*t.tv_nsec;
239 return seconds;
240 #else
241 return gmx_gettime();
242 #endif