Line data Source code
1 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3 : !!!! !!!!
4 : !!!! ParaMonte: Parallel Monte Carlo and Machine Learning Library. !!!!
5 : !!!! !!!!
6 : !!!! Copyright (C) 2012-present, The Computational Data Science Lab !!!!
7 : !!!! !!!!
8 : !!!! This file is part of the ParaMonte library. !!!!
9 : !!!! !!!!
10 : !!!! LICENSE !!!!
11 : !!!! !!!!
12 : !!!! https://github.com/cdslaborg/paramonte/blob/main/LICENSE.md !!!!
13 : !!!! !!!!
14 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
15 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
16 :
17 : !> \brief
18 : !> This module contains the timer procedures and derived types to facilitate timing applications at runtime.<br>
19 : !>
20 : !> \details
21 : !> Currently, five types of timers are available in this module in addition to a generic timer class.<br>
22 : !> <ol>
23 : !> <li> [timer_type](@ref pm_timer::timer_type)
24 : !> <ol>
25 : !> <li> This is a generic timer class which defaults to the most appropriate timer type based on the ParaMonte library build.<br>
26 : !> <ol>
27 : !> <li> If the ParaMonte library is built for serial applications,
28 : !> then [timerSYS_type](@ref pm_timer::timerSYS_type) is used.<br>
29 : !> <li> If the ParaMonte library is built for MPI-parallel applications,
30 : !> then [timerMPI_type](@ref pm_timer::timerMPI_type) is used.<br>
31 : !> <li> If the ParaMonte library is built for OpenMP-parallel applications,
32 : !> then [timerOMP_type](@ref pm_timer::timerOMP_type) is used.<br>
33 : !> </ol>
34 : !> <li> Use this timer class if you are not sure which timer is suitable for the task.<br>
35 : !> </ol>
36 : !> <li> [timerCPU_type](@ref pm_timer::timerCPU_type)
37 : !> <ol>
38 : !> <li> This timer relies on the Fortran intrinsic `cpu_time()`.<br>
39 : !> <li> This timer is **not ideal** for timing parallel or multi-threaded applications
40 : !> because `cpu_time()` returns the sum of the total time spent by all processes.<br>
41 : !> <li> The timing precision of the Fortran intrinsic `cpu_time()` is processor-dependent.<br>
42 : !> <li> There is currently no direct way of measuring the time resolution of `cpu_time()`.<br>
43 : !> <li> This module uses a simple cycling method to compute the resolution of the time returned by `cpu_time()`.<br>
44 : !> It computes the resolution as the least noticeable non-zero time difference between any two time points returned by `cpu_time()`.<br>
45 : !> <li> On processors with no clock, the time returned by `cpu_time()` is negative to indicate the lack of a processor clock.<br>
46 : !> The clock existence is verified only if the library is built with the preprocessor macro `CHECK_ENABLED=1`.<br>
47 : !> </ol>
48 : !> <li> [timerDAT_type](@ref pm_timer::timerDAT_type)
49 : !> <ol>
50 : !> <li> This timer relies on the Fortran intrinsic `date_and_time()`.<br>
51 : !> <li> This timer is **not ideal** for high resolution time measurements
52 : !> because the resolution is limited to milliseconds.<br>
53 : !> </ol>
54 : !> <li> [timerMPI_type](@ref pm_timer::timerMPI_type)
55 : !> <ol>
56 : !> <li> This timer relies on the MPI library intrinsic `mpi_wtime()`.<br>
57 : !> <li> The clock time resolution is computed by calling `mpi_wtick()`.<br>
58 : !> <li> This timer is **ideal** for high resolution (e.g., nanseconds) time measurements in distributed parallel applications.<br>
59 : !> <li> This timer is available only when the ParaMonte library is built with MPI support enabled.<br>
60 : !> <li> This timer requires the MPI library to have been already initialized via `mpi_init()`
61 : !> and not have been finalized via `mpi_finalize()`.<br>
62 : !> \vericons
63 : !> </ol>
64 : !> <li> [timerOMP_type](@ref pm_timer::timerOMP_type)
65 : !> <ol>
66 : !> <li> This timer relies on the OpenMP library intrinsic `omp_get_wtime()`.<br>
67 : !> <li> The clock time resolution is computed by calling `omp_get_wtick()`.<br>
68 : !> <li> This timer is **ideal** for high resolution time measurements in multithreaded parallel applications.<br>
69 : !> <li> This timer is available only when the ParaMonte library is built with OpenMP support enabled.<br>
70 : !> <li> This timer measures elapsed wall clock time in seconds.<br>
71 : !> <li> The time is measured per thread.<br>
72 : !> <li> No guarantee can be made that two distinct threads measure the same time.<br>
73 : !> <li> Time is measured from **some time in the past**, which is an arbitrary time
74 : !> guaranteed not to change during the execution of the program.<br>
75 : !> </ol>
76 : !> <li> [timerSYS_type](@ref pm_timer::timerSYS_type)
77 : !> <ol>
78 : !> <li> This timer relies on the Fortran intrinsic `system_clock()`.<br>
79 : !> <li> This timer can be used in most parallel or serial timing scenarios.<br>
80 : !> <li> The resolution of this timer is processor-dependent.<br>
81 : !> <li> Time is measured from **some time in the past**, which is an arbitrary
82 : !> time guaranteed not to change during the execution of the program.<br>
83 : !> <li> The time is measured per processor.<br>
84 : !> <li> On processors with no clock, `sysctem_clock()` returns a negative clock `count` to indicate the lack of a processor clock.<br>
85 : !> The clock existence is verified only if the library is built with the preprocessor macro `CHECK_ENABLED=1`.<br>
86 : !> </ol>
87 : !> </ol>
88 : !>
89 : !> \test
90 : !> [test_pm_timer](@ref test_pm_timer)
91 : !>
92 : !> \finmain
93 : !>
94 : !> \author
95 : !> \AmirShahmoradi, March 22, 2012, 2:21 PM, National Institute for Fusion Studies, The University of Texas at Austin
96 :
97 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
98 :
99 : module pm_timer
100 :
101 : #if CHECK_ENABLED
102 : #define CHECK_ASSERTION(LINE,ASSERTION,MSG) \
103 : block; \
104 : use pm_err, only: getFine; \
105 : use pm_val2str, only: getStr; \
106 : use pm_err, only: setAsserted; \
107 : call setAsserted(ASSERTION,getFine(__FILE__,LINE)//MODULE_NAME//MSG); \
108 : end block;
109 : #else
110 : #define CHECK_ASSERTION(LINE,ASSERTION,MSG) continue;
111 : #endif
112 :
113 : use pm_kind, only: SK, LK, IKD, RKD
114 :
115 : implicit none
116 :
117 : character(*, SK), parameter :: MODULE_NAME = "@pm_timer"
118 :
119 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
120 :
121 : !> \brief
122 : !> This is the `abstract interface` of the [setTime](@ref pm_timer::timer_type) static type-bound procedure
123 : !> component of [timer_type](@ref pm_timer::timer_type) `abstract` type that performs the timing since a user-specified or a processor-dependent origin.
124 : !>
125 : !> \finmain{getTime_proc}
126 : !>
127 : !> \author
128 : !> \AmirShahmoradi, March 23, 2012, 2:21 AM, National Institute for Fusion Studies, The University of Texas at Austin
129 : abstract interface
130 : function getTime_proc(since) result(timeInSec)
131 : use pm_kind, only: RKD
132 : real(RKD), intent(in), optional :: since
133 : real(RKD) :: timeInSec
134 : end function
135 : end interface
136 :
137 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
138 :
139 : !> \brief
140 : !> This is the `abstract interface` of the [wait](@ref pm_timer::timer_type) static type-bound procedure component of
141 : !> [timer_type](@ref pm_timer::timer_type) `abstract` type that keeps the system waiting for the specified amount of seconds.<br>
142 : !>
143 : !> \details
144 : !> This procedure is functionally equivalent to the popular non-standard `sleep()` procedure offered by some compiler vendors.<br>
145 : !>
146 : !> \note
147 : !> Note that even though the system sleeps during the idle time, the timer clock is still ticking.<br>
148 : !>
149 : !> \finmain{setIdle_proc}
150 : !>
151 : !> \author
152 : !> \AmirShahmoradi, March 23, 2012, 2:21 AM, National Institute for Fusion Studies, The University of Texas at Austin
153 : abstract interface
154 : subroutine setIdle_proc(seconds)
155 : use pm_kind, only: RKD
156 : real(RKD), intent(in) :: seconds
157 : end subroutine
158 : end interface
159 :
160 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
161 :
162 : !> \brief
163 : !> This is the `abstract` base derived type that serves as a simple
164 : !> container template for other timer classes in the ParaMonte library.
165 : !>
166 : !> \details
167 : !> Timer classes derived from based on this `abstract` type minimally contain components to store,
168 : !> <ol>
169 : !> <li> the **start** time of a timer in seconds (since the processor-dependent epoch).
170 : !> <li> the **clock** in seconds (optionally computed since the start of the timer).
171 : !> <li> the **delta** elapsed time (optionally computed since the last call to the timer in seconds).
172 : !> <li> the **resolution** of the timer in seconds.
173 : !> </ol>
174 : !> See the documentation of [constructTimer](@ref pm_timer::constructTimer)
175 : !> for the non-default class constructor interface for this type.<br>
176 : !> See also the documentation details of [pm_timer](@ref pm_timer).
177 : !>
178 : !> \interface{timer_type}
179 : !> \code{.F90}
180 : !>
181 : !> use pm_timer, only: timer_type
182 : !> class(timer_type), allocatable :: timer
183 : !>
184 : !> timer = timer_type()
185 : !>
186 : !> \endcode
187 : !>
188 : !> \see
189 : !> [timer_type](@ref pm_timer::timer_type)<br>
190 : !> [timerDAT_type](@ref pm_timer::timerDAT_type)<br>
191 : !> [timerMPI_type](@ref pm_timer::timerMPI_type)<br>
192 : !> [timerOMP_type](@ref pm_timer::timerOMP_type)<br>
193 : !> [timerSYS_type](@ref pm_timer::timerSYS_type)<br>
194 : !>
195 : !> \example{timer_type}
196 : !> \include{lineno} example/pm_timer/timer_type/main.F90
197 : !> \compilef{timer_type}
198 : !> \output{timer_type}
199 : !> \include{lineno} example/pm_timer/timer_type/main.out.F90
200 : !>
201 : !> \test
202 : !> [test_pm_timer](@ref test_pm_timer)
203 : !>
204 : !> \todo
205 : !> \pvhigh
206 : !> The `start` component of this derived type must become `protected` as soon as Fortran 2023 is supported by the compilers.<br>
207 : !>
208 : !> \finmain{timer_type}
209 : !>
210 : !> \author
211 : !> \AmirShahmoradi, March 22, 2012, 2:21 PM, National Institute for Fusion Studies, The University of Texas at Austin
212 : type, abstract :: timer_type
213 : real(RKD) :: start !< \public The `protected` scalar `real` of kind `double precision` \RKD provided as a convenience for the user to contain the start time of the timer since the processor epoch (a processor-dependent past time) in seconds.<br>
214 : !! \warning The object of class [timer_type](@ref pm_timer::timer_type) must be initialized with the constructor (e.g., `timer = timer_type()`) to appropriately set the object `start` component to the processor epoch.<br>
215 : real(RKD) :: clock !< \public The scalar `real` of kind `double precision` \RKD provided as a convenience for the user to contain the total time in seconds elapsed since the start of the timer.<br>
216 : !! The clock time can be readily can be computed as `timer%%clock = timer%%time(since = timer%%start)`.<br>
217 : real(RKD) :: delta !< \public The scalar `real` of kind `double precision` \RKD provided as a convenience for the user to contain the delta time in seconds since the last timing (last timer call).<br>
218 : !! The delta time can be readily can be computed as `timer%%delta = timer%%time(since = timer%%start) - timer%%clock` where `timer%%clock` is the last clock read as `timer%%clock = timer%%time(since = timer%%start)`.<br>
219 : real(RKD) :: resol !< \public The `protected` scalar `real` of kind `double precision` \RKD provided as a convenience for the user to contain the time in seconds between the timer clock tics.<br>
220 : !! \warning The object of class [timer_type](@ref pm_timer::timer_type) must be initialized with the constructor (e.g., `timer = timer_type()`) to appropriately set the object `resol` component to the clock time resolution.<br>
221 : contains
222 : procedure(getTime_proc), nopass, deferred :: time !< \public See [getTime_proc](@ref pm_timer::getTime_proc).
223 : procedure(setIdle_proc), nopass, deferred :: wait !< \public See [setIdle_proc](@ref pm_timer::setIdle_proc).
224 : end type
225 :
226 : interface timer_type
227 : module procedure :: constructTimer
228 : end interface
229 :
230 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
231 :
232 : !> \brief
233 : !> This is the `timerCPU_type` class, containing attributes and static methods for
234 : !> setting up a timer based on the CPU-clock, using the Fortran intrinsic `cpu_time()`.
235 : !>
236 : !> \details
237 : !> See the documentation of [constructTimerCPU](@ref pm_timer::constructTimerCPU)
238 : !> for the non-default constructor interface of this type.<br>
239 : !> See also the documentation details of [pm_timer](@ref pm_timer).
240 : !>
241 : !> \interface{timerCPU_type}
242 : !> \code{.F90}
243 : !>
244 : !> use pm_timer, only: timerCPU_type
245 : !> type(timerCPU_type) :: timer
246 : !>
247 : !> timer = timerCPU_type()
248 : !>
249 : !> \endcode
250 : !>
251 : !> \see
252 : !> [timer_type](@ref pm_timer::timer_type)<br>
253 : !> [timerDAT_type](@ref pm_timer::timerDAT_type)<br>
254 : !> [timerMPI_type](@ref pm_timer::timerMPI_type)<br>
255 : !> [timerOMP_type](@ref pm_timer::timerOMP_type)<br>
256 : !> [timerSYS_type](@ref pm_timer::timerSYS_type)<br>
257 : !>
258 : !> \example{timerCPU_type}
259 : !> \include{lineno} example/pm_timer/timerCPU_type/main.F90
260 : !> \compilef{timerCPU_type}
261 : !> \output{timerCPU_type}
262 : !> \include{lineno} example/pm_timer/timerCPU_type/main.out.F90
263 : !>
264 : !> \test
265 : !> [test_pm_timer](@ref test_pm_timer)
266 : !>
267 : !> \finmain{timerCPU_type}
268 : !>
269 : !> \author
270 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
271 : type, extends(timer_type) :: timerCPU_type
272 : contains
273 : procedure, nopass :: time => getTimeCPU !< \public See [getTime_proc](@ref pm_timer::getTime_proc).
274 : procedure, nopass :: wait => setIdleCPU !< \public See [setIdle_proc](@ref pm_timer::setIdle_proc).
275 : end type
276 :
277 : interface timerCPU_type
278 : module procedure :: constructTimerCPU
279 : end interface
280 :
281 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
282 :
283 : !> \brief
284 : !> This is the `timerDAT_type` class, containing attributes and static methods for
285 : !> setting up a timer based on the system-clock, using the Fortran intrinsic `date_and_time()`.
286 : !>
287 : !> \details
288 : !> See the documentation of [constructTimerDAT](@ref pm_timer::constructTimerDAT)
289 : !> for the non-default constructor interface of this type.<br>
290 : !> See also the documentation details of [pm_timer](@ref pm_timer).
291 : !>
292 : !> \interface{timerDAT_type}
293 : !> \code{.F90}
294 : !>
295 : !> use pm_timer, only: timerDAT_type
296 : !> type(timerDAT_type) :: timer
297 : !>
298 : !> timer = timerDAT_type()
299 : !>
300 : !> \endcode
301 : !>
302 : !> \see
303 : !> [timer_type](@ref pm_timer::timer_type)<br>
304 : !> [timerDAT_type](@ref pm_timer::timerDAT_type)<br>
305 : !> [timerMPI_type](@ref pm_timer::timerMPI_type)<br>
306 : !> [timerOMP_type](@ref pm_timer::timerOMP_type)<br>
307 : !> [timerSYS_type](@ref pm_timer::timerSYS_type)<br>
308 : !>
309 : !> \example{timerDAT_type}
310 : !> \include{lineno} example/pm_timer/timerDAT_type/main.F90
311 : !> \compilef{timerDAT_type}
312 : !> \output{timerDAT_type}
313 : !> \include{lineno} example/pm_timer/timerDAT_type/main.out.F90
314 : !>
315 : !> \test
316 : !> [test_pm_timer](@ref test_pm_timer)
317 : !>
318 : !> \finmain{timerDAT_type}
319 : !>
320 : !> \author
321 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
322 : type, extends(timer_type) :: timerDAT_type
323 : contains
324 : procedure, nopass :: time => getTimeDAT !< \public See [getTime_proc](@ref pm_timer::getTime_proc).
325 : procedure, nopass :: wait => setIdleDAT !< \public See [setIdle_proc](@ref pm_timer::setIdle_proc).
326 : end type
327 :
328 : interface timerDAT_type
329 : module procedure :: constructTimerDAT
330 : end interface
331 :
332 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
333 :
334 : !> \brief
335 : !> This is the `timerMPI_type` class, containing attributes and static
336 : !> methods for setting up a timer based on the MPI intrinsic `MPI_Wtime()`.
337 : !>
338 : !> \details
339 : !> See the documentation of [constructTimerMPI](@ref pm_timer::constructTimerMPI)
340 : !> for the non-default constructor interface of this type.<br>
341 : !> See also the documentation details of [pm_timer](@ref pm_timer).
342 : !>
343 : !> \interface{timerMPI_type}
344 : !> \code{.F90}
345 : !>
346 : !> use pm_timer, only: timerMPI_type
347 : !> type(timerMPI_type) :: timer
348 : !>
349 : !> timer = timerMPI_type()
350 : !>
351 : !> \endcode
352 : !>
353 : !> \note
354 : !> This type is available only if the library is built with the preprocessor macro `MPI_ENABLED=1`.
355 : !>
356 : !> \see
357 : !> [timer_type](@ref pm_timer::timer_type)<br>
358 : !> [timerDAT_type](@ref pm_timer::timerDAT_type)<br>
359 : !> [timerMPI_type](@ref pm_timer::timerMPI_type)<br>
360 : !> [timerOMP_type](@ref pm_timer::timerOMP_type)<br>
361 : !> [timerSYS_type](@ref pm_timer::timerSYS_type)<br>
362 : !>
363 : !> \example{timerMPI_type}
364 : !> \include{lineno} example/pm_timer/timerMPI_type/main.F90
365 : !> \compilef{timerMPI_type}
366 : !> \output{timerMPI_type}
367 : !> \include{lineno} example/pm_timer/timerMPI_type/main.out.F90
368 : !>
369 : !> \test
370 : !> [test_pm_timer](@ref test_pm_timer)
371 : !>
372 : !> \finmain{timerMPI_type}
373 : !>
374 : !> \author
375 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
376 : #if MPI_ENABLED
377 : type, extends(timer_type) :: timerMPI_type
378 : contains
379 : procedure, nopass :: time => getTimeMPI !< \public See [getTime_proc](@ref pm_timer::getTime_proc).
380 : procedure, nopass :: wait => setIdleMPI !< \public See [setIdle_proc](@ref pm_timer::setIdle_proc).
381 : end type
382 :
383 : interface timerMPI_type
384 : module procedure :: constructTimerMPI
385 : end interface
386 : #endif
387 :
388 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
389 :
390 : !> \brief
391 : !> This is the `timerMPI_type` class, containing attributes and static
392 : !> methods for setting up a timer based on the OpenMP intrinsic `omp_get_wtime()`.
393 : !>
394 : !> \details
395 : !> See the documentation of [constructTimerOMP](@ref pm_timer::constructTimerOMP)
396 : !> for the non-default constructor interface of this type.<br>
397 : !> See also the documentation details of [pm_timer](@ref pm_timer).
398 : !>
399 : !> \interface{timerOMP_type}
400 : !> \code{.F90}
401 : !>
402 : !> use pm_timer, only: timerOMP_type
403 : !> type(timerOMP_type) :: timer
404 : !>
405 : !> timer = timerOMP_type()
406 : !>
407 : !> \endcode
408 : !>
409 : !> \note
410 : !> This type is available only if the library is built with the preprocessor macro `OMP_ENABLED=1`.
411 : !>
412 : !> \see
413 : !> [timer_type](@ref pm_timer::timer_type)<br>
414 : !> [timerDAT_type](@ref pm_timer::timerDAT_type)<br>
415 : !> [timerMPI_type](@ref pm_timer::timerMPI_type)<br>
416 : !> [timerOMP_type](@ref pm_timer::timerOMP_type)<br>
417 : !> [timerSYS_type](@ref pm_timer::timerSYS_type)<br>
418 : !>
419 : !> \example{timerOMP_type}
420 : !> \include{lineno} example/pm_timer/timerOMP_type/main.F90
421 : !> \compilef{timerOMP_type}
422 : !> \output{timerOMP_type}
423 : !> \include{lineno} example/pm_timer/timerOMP_type/main.out.F90
424 : !>
425 : !> \test
426 : !> [test_pm_timer](@ref test_pm_timer)
427 : !>
428 : !> \finmain{timerOMP_type}
429 : !>
430 : !> \author
431 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
432 : #if OMP_ENABLED
433 : type, extends(timer_type) :: timerOMP_type
434 : contains
435 : procedure, nopass :: time => getTimeOMP !< \public See [getTime_proc](@ref pm_timer::getTime_proc).
436 : procedure, nopass :: wait => setIdleOMP !< \public See [setIdle_proc](@ref pm_timer::setIdle_proc).
437 : end type
438 :
439 : interface timerOMP_type
440 : module procedure :: constructTimerOMP
441 : end interface
442 : #endif
443 :
444 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
445 :
446 : !!> \brief
447 : !!> This s the derived type for generating objects containing information about the system clock count.
448 : !!>
449 : !!> \finmain{Count_type}
450 : !!>
451 : !!> \author
452 : !!> Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
453 : !type :: SysClockCount_type
454 : ! integer(IKD) :: start = 0_IKD !< The scalar `integer` of kind double precision \IKD containing the processor-dependent starting clock counts of the system.
455 : ! integer(IKD) :: clock = 0_IKD !< The scalar `integer` of kind double precision \IKD containing the total processor clock counts since `start`.
456 : ! integer(IKD) :: delta = 0_IKD !< The scalar `integer` of kind double precision \IKD containing total processor clock count since the last clock count measurement.
457 : ! integer(IKD) :: max = 0_IKD !< The scalar `integer` of kind double precision \IKD containing maximum value that the processor count may take, or `0` if there is no clock.
458 : ! integer(IKD) :: rate = 0_IKD !< The scalar `integer` of kind double precision \IKD containing number of clock counts per second, or `0` if there is no clock.
459 : !end type
460 :
461 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
462 :
463 : !> \brief
464 : !> This is the `timerSYS_type` class, containing attributes and static methods for
465 : !> setting up a timer based on the system-clock, using the Fortran intrinsic `system_clock()`.
466 : !>
467 : !> \details
468 : !> See the documentation of [constructTimerSYS](@ref pm_timer::constructTimerSYS)
469 : !> for the non-default constructor interface of this type.<br>
470 : !> See also the documentation details of [pm_timer](@ref pm_timer).
471 : !>
472 : !> \interface{timerSYS_type}
473 : !> \code{.F90}
474 : !>
475 : !> use pm_timer, only: timerSYS_type
476 : !> type(timerSYS_type) :: timer
477 : !>
478 : !> timer = timerSYS_type()
479 : !>
480 : !> \endcode
481 : !>
482 : !> \see
483 : !> [timer_type](@ref pm_timer::timer_type)<br>
484 : !> [timerDAT_type](@ref pm_timer::timerDAT_type)<br>
485 : !> [timerMPI_type](@ref pm_timer::timerMPI_type)<br>
486 : !> [timerOMP_type](@ref pm_timer::timerOMP_type)<br>
487 : !> [timerSYS_type](@ref pm_timer::timerSYS_type)<br>
488 : !>
489 : !> \example{timerSYS_type}
490 : !> \include{lineno} example/pm_timer/timerSYS_type/main.F90
491 : !> \compilef{timerSYS_type}
492 : !> \output{timerSYS_type}
493 : !> \include{lineno} example/pm_timer/timerSYS_type/main.out.F90
494 : !>
495 : !> \test
496 : !> [test_pm_timer](@ref test_pm_timer)
497 : !>
498 : !> \finmain{timerSYS_type}
499 : !>
500 : !> \author
501 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
502 : type, extends(timer_type) :: timerSYS_type
503 : !type(SysClockCount_type) :: Count !< An object of type [SysClockCount_type](@ref pm_timer::SysClockCount_type) containing information about the system clock.
504 : contains
505 : procedure, nopass :: time => getTimeSYS !< \public See [getTime_proc](@ref pm_timer::getTime_proc).
506 : procedure, nopass :: wait => setIdleSYS !< \public See [setIdle_proc](@ref pm_timer::setIdle_proc).
507 : end type
508 :
509 : interface timerSYS_type
510 : module procedure :: constructTimerSYS
511 : end interface
512 :
513 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
514 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
515 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
516 :
517 : contains
518 :
519 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
520 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
521 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
522 :
523 : !> \brief
524 : !> This is the constructor of the [timer_type](@ref pm_timer::timer_type) class.<br>
525 : !>
526 : !> \details
527 : !> Upon return, the constructor initializes all components of the timer object.<br>
528 : !> See also the documentation details of [pm_timer](@ref pm_timer).<br>
529 : !>
530 : !> \return
531 : !> `timer` : The output scalar object of class [timer_type](@ref pm_timer::timer_type).<br>
532 : !> Specifically, the output of type,<br>
533 : !> <ol>
534 : !> <li> [timerMPI_type](@ref pm_timer::timerMPI_type) if the ParaMonte library is built with preprocessor macro `MPI_ENABLED` defined.
535 : !> <li> [timerOMP_type](@ref pm_timer::timerOMP_type) if the ParaMonte library is built with preprocessor macro `OMP_ENABLED` defined.
536 : !> <li> [timerSYS_type](@ref pm_timer::timerSYS_type) if the ParaMonte library is built is serial mode.
537 : !> </ol>
538 : !>
539 : !> \interface{constructTimer}
540 : !> \code{.F90}
541 : !>
542 : !> use pm_timer, only: timer_type
543 : !> class(timer_type), allocatable :: timer
544 : !>
545 : !> timer = timer_type()
546 : !>
547 : !> \endcode
548 : !>
549 : !> \remark
550 : !> See the documentation of [timer_type](@ref pm_timer::timer_type) for example usage.
551 : !>
552 : !> \finmain{constructTimer}
553 : !>
554 : !> \author
555 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
556 3427 : function constructTimer() result(timer)
557 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
558 : !DEC$ ATTRIBUTES DLLEXPORT :: constructTimer
559 : #endif
560 : #if MPI_ENABLED
561 : type(timerMPI_type) :: timer
562 : timer = timerMPI_type()
563 : #elif OMP_ENABLED
564 : type(timerOMP_type) :: timer
565 : timer = timerOMP_type()
566 : #else
567 : type(timerSYS_type) :: timer
568 3427 : timer = timerSYS_type()
569 : #endif
570 3427 : end function
571 :
572 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
573 :
574 : !> \brief
575 : !> This is the constructor of the [timerCPU_type](@ref pm_timer::timerCPU_type) class.<br>
576 : !>
577 : !> \details
578 : !> Upon return, the constructor initializes all components of the timer object.<br>
579 : !> See also the documentation details of [pm_timer](@ref pm_timer).
580 : !>
581 : !> \return
582 : !> `timer` : The output scalar object of class [timerCPU_type](@ref pm_timer::timerCPU_type).
583 : !>
584 : !> \interface{constructTimerCPU}
585 : !> \code{.F90}
586 : !>
587 : !> use pm_timer, only: timerCPU_type
588 : !> type(timerCPU_type) :: timer
589 : !>
590 : !> timer = timerCPU_type()
591 : !>
592 : !> \endcode
593 : !>
594 : !> \remark
595 : !> See the documentation of [timerCPU_type](@ref pm_timer::timerCPU_type) for example usage.
596 : !>
597 : !> \finmain{constructTimerCPU}
598 : !>
599 : !> \author
600 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
601 6 : function constructTimerCPU() result(timer)
602 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
603 : !DEC$ ATTRIBUTES DLLEXPORT :: constructTimerCPU
604 : #endif
605 : type(timerCPU_type) :: timer
606 6 : timer%resol = getResTimerCPU()
607 6 : timer%start = timer%time()
608 6 : timer%delta = 0._RKD
609 6 : timer%clock = 0._RKD
610 6 : end function
611 :
612 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
613 :
614 : !> \brief
615 : !> This is the constructor of the [timerDAT_type](@ref pm_timer::timerDAT_type) class.<br>
616 : !>
617 : !> \details
618 : !> Upon return, the constructor initializes all components of the timer object.<br>
619 : !> See also the documentation details of [pm_timer](@ref pm_timer).
620 : !>
621 : !> \interface{constructTimerDAT}
622 : !> \code{.F90}
623 : !>
624 : !> use pm_timer, only: timerDAT_type
625 : !> type(timerDAT_type) :: timer
626 : !>
627 : !> timer = timerDAT_type()
628 : !>
629 : !> \endcode
630 : !>
631 : !> \return
632 : !> `timer` : The output scalar object of class [timerDAT_type](@ref pm_timer::timerDAT_type).
633 : !>
634 : !> \remark
635 : !> See the documentation of [timerDAT_type](@ref pm_timer::timerDAT_type) for example usage.
636 : !>
637 : !> \finmain{constructTimerDAT}
638 : !>
639 : !> \author
640 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
641 6 : function constructTimerDAT() result(timer)
642 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
643 : !DEC$ ATTRIBUTES DLLEXPORT :: constructTimerDAT
644 : #endif
645 : type(timerDAT_type) :: timer
646 6 : timer%resol = 0.001_RKD
647 6 : timer%start = timer%time()
648 6 : timer%delta = 0._RKD
649 6 : timer%clock = 0._RKD
650 6 : end function
651 :
652 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
653 :
654 : !> \brief
655 : !> This is the constructor of the [timerMPI_type](@ref pm_timer::timerMPI_type) class.<br>
656 : !>
657 : !> \details
658 : !> Upon return, the constructor initializes all components of the timer object.<br>
659 : !> See also the documentation details of [pm_timer](@ref pm_timer).
660 : !>
661 : !> \interface{constructTimerMPI}
662 : !> \code{.F90}
663 : !>
664 : !> use pm_timer, only: timerMPI_type
665 : !> type(timerMPI_type) :: timer
666 : !>
667 : !> timer = timerMPI_type()
668 : !>
669 : !> \endcode
670 : !>
671 : !> \return
672 : !> `timer` : The output scalar object of class [timerMPI_type](@ref pm_timer::timerMPI_type).
673 : !>
674 : !> \warning
675 : !> This constructor requires the MPI library to have been initialized and not finalized.<br>
676 : !> \vericon
677 : !>
678 : !> \warning
679 : !> This constructor is available only if the library is built with the preprocessor macro `MPI_ENABLED=1`.<br>
680 : !>
681 : !> \remark
682 : !> See the documentation of [timerMPI_type](@ref pm_timer::timerMPI_type) for example usage.
683 : !>
684 : !> \finmain{constructTimerMPI}
685 : !>
686 : !> \author
687 : !> \AmirShahmoradi, March 22, 2012, 02:51 AM, National Institute for Fusion Studies, The University of Texas at Austin
688 : #if MPI_ENABLED
689 : function constructTimerMPI() result(timer)
690 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
691 : !DEC$ ATTRIBUTES DLLEXPORT :: constructTimerMPI
692 : #endif
693 : use mpi !mpi_f08, only: mpi_initialized, mpi_init, mpi_wtick
694 : type(timerMPI_type) :: timer
695 : logical :: initialized
696 : integer :: ierrMPI
697 : !#if CHECK_ENABLED
698 : ! block
699 : ! use pm_err, only: setAsserted
700 : ! use pm_val2str, only: getStr
701 : ! logical :: initialized, finalized
702 : ! integer :: ierrMPI
703 : ! call mpi_initialized(initialized, ierrMPI)
704 : ! call setAsserted(ierrMPI /= 0 .and. initialized, MODULE_NAME//SK_"@constructTimerMPI(): The MPI library must be initialized before attempting to call the MPI timer.")
705 : ! call mpi_finalized(finalized, ierrMPI)
706 : ! call setAsserted(ierrMPI /= 0 .and. finalized, MODULE_NAME//SK_"@constructTimerMPI(): The MPI library must not be finalized prior to calling the MPI timer.")
707 : ! end block
708 : !#endif
709 : call mpi_initialized(initialized, ierrMPI)
710 : if (.not. initialized .and. ierrMPI == 0) call mpi_init(ierrMPI)
711 : if (ierrMPI /= 0) error stop MODULE_NAME//SK_"@constructTimerMPI(): Failed to initialize the MPI library."
712 : timer%resol = real(mpi_wtick(), RKD)
713 : timer%start = timer%time()
714 : timer%delta = 0._RKD
715 : timer%clock = 0._RKD
716 : end function
717 : #endif
718 :
719 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
720 :
721 : !> \brief
722 : !> This is the constructor of the [timerOMP_type](@ref pm_timer::timerOMP_type) class.<br>
723 : !>
724 : !> \details
725 : !> Upon return, the constructor initializes all components of the timer object.<br>
726 : !> See also the documentation details of [pm_timer](@ref pm_timer).
727 : !>
728 : !> \interface{constructTimerOMP}
729 : !> \code{.F90}
730 : !>
731 : !> use pm_timer, only: timerOMP_type
732 : !> type(timerOMP_type) :: timer
733 : !>
734 : !> timer = timerOMP_type()
735 : !>
736 : !> \endcode
737 : !>
738 : !> \return
739 : !> `timer` : The output scalar object of class [timerOMP_type](@ref pm_timer::timerOMP_type).
740 : !>
741 : !> \warning
742 : !> This constructor is available only if the library is built with the preprocessor macro `OMP_ENABLED=1`.<br>
743 : !>
744 : !> \remark
745 : !> See the documentation of [timerOMP_type](@ref pm_timer::timerOMP_type) for example usage.
746 : !>
747 : !> \finmain{constructTimerOMP}
748 : !>
749 : !> \author
750 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
751 : #if OMP_ENABLED
752 : function constructTimerOMP() result(timer)
753 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
754 : !DEC$ ATTRIBUTES DLLEXPORT :: constructTimerOMP
755 : #endif
756 : use omp_lib
757 : type(timerOMP_type) :: timer
758 : timer%resol = real(omp_get_wtick(), RKD)
759 : timer%start = timer%time()
760 : timer%delta = 0._RKD
761 : timer%clock = 0._RKD
762 : end function
763 : #endif
764 :
765 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
766 :
767 : !> \brief
768 : !> This is the constructor of the [timerSYS_type](@ref pm_timer::timerSYS_type) class.<br>
769 : !>
770 : !> \details
771 : !> Upon return, the constructor initializes all components of the timer object.<br>
772 : !> See also the documentation details of [pm_timer](@ref pm_timer).
773 : !>
774 : !> \interface{constructTimerSYS}
775 : !> \code{.F90}
776 : !>
777 : !> use pm_timer, only: timerSYS_type
778 : !> type(timerSYS_type) :: timer
779 : !>
780 : !> timer = timerSYS_type()
781 : !>
782 : !> \endcode
783 : !>
784 : !> \return
785 : !> `timer` : The output scalar object of class [timerSYS_type](@ref pm_timer::timerSYS_type).
786 : !>
787 : !> \remark
788 : !> See the documentation of [timerSYS_type](@ref pm_timer::timerSYS_type) for example usage.
789 : !>
790 : !> \finmain{constructTimerSYS}
791 : !>
792 : !> \author
793 : !> \AmirShahmoradi, March 22, 2012, 03:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
794 3432 : function constructTimerSYS() result(timer)
795 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
796 : !DEC$ ATTRIBUTES DLLEXPORT :: constructTimerSYS
797 : #endif
798 : type(timerSYS_type) :: timer
799 3432 : timer%resol = getResTimerSYS()
800 3432 : timer%start = timer%time()
801 3432 : timer%delta = 0._RKD
802 3432 : timer%clock = 0._RKD
803 3432 : end function
804 :
805 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
806 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
807 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
808 :
809 : !> \brief
810 : !> Generate and return the time in units of seconds since the
811 : !> specified time `since` or since an arbitrary processor-dependent time.
812 : !>
813 : !> \details
814 : !> This function relies on either of the following lower-level module functions depending on the ParaMonte library build:<br>
815 : !> <ol>
816 : !> <li> [getTimeMPI](@ref pm_timer::getTimeMPI) for MPI-enabled parallel builds of the ParaMonte library.<br>
817 : !> <li> [getTimeOMP](@ref pm_timer::getTimeOMP) for OpenMP-enabled parallel builds of the ParaMonte library.<br>
818 : !> <li> [getTimeSYS](@ref pm_timer::getTimeSYS) for all other parallel or serial builds of the ParaMonte library.<br>
819 : !> </ol>
820 : !>
821 : !> \param[in] since : The input scalar of type `real` of kind double precision \RKD representing the time origin (epoch) in seconds.<br>
822 : !> (**optional**. If not present, its value is processor-dependent but constant at runtime.)
823 : !>
824 : !> \return
825 : !> `timeInSec` : The output scalar of type `real` of kind double precision \RKD, containing the time in
826 : !> seconds since a specified time `since` or since an arbitrary processor-dependent time.
827 : !>
828 : !> \interface{getTime}
829 : !> \code{.F90}
830 : !>
831 : !> use pm_timer, only: getTime
832 : !> real(RKD) :: timeInSec
833 : !> logical(LK) :: failed
834 : !>
835 : !> timeInSec = getTime(since = since)
836 : !>
837 : !> \endcode
838 : !>
839 : !> \impure
840 : !>
841 : !> \see
842 : !> [timer_type](@ref pm_timer::timer_type)<br>
843 : !>
844 : !> \example{getTime}
845 : !> \include{lineno} example/pm_timer/getTime/main.F90
846 : !> \compilef{getTime}
847 : !> \output{getTime}
848 : !> \include{lineno} example/pm_timer/getTime/main.out.F90
849 : !>
850 : !> \test
851 : !> [test_pm_timer](@ref test_pm_timer)
852 : !>
853 : !> \finmain{getTime}
854 : !>
855 : !> \author
856 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
857 4 : function getTime(since) result(timeInSec)
858 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
859 : !DEC$ ATTRIBUTES DLLEXPORT :: getTime
860 : #endif
861 : real(RKD), intent(in), optional :: since
862 : real(RKD) :: timeInSec
863 : #if MPI_ENABLED
864 : timeInSec = getTimeMPI(since)
865 : #elif OMP_ENABLED
866 : timeInSec = getTimeOMP(since)
867 : #else
868 4 : timeInSec = getTimeSYS(since)
869 : #endif
870 4 : end function
871 :
872 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
873 :
874 : !> \brief
875 : !> Generate and return the CPU time in units of seconds since the
876 : !> specified time `since` or since an arbitrary processor-dependent time.
877 : !>
878 : !> \details
879 : !> See also the documentation details of [pm_timer](@ref pm_timer).
880 : !>
881 : !> \param[in] since : The input scalar of type `real` of kind double precision \RKD representing the time origin (epoch) in seconds.<br>
882 : !> (**optional**. If not present, its value is processor-dependent but constant at runtime.)
883 : !>
884 : !> \return
885 : !> `timeInSec` : The output scalar of type `real` of kind double precision \RKD, containing the time in
886 : !> seconds since a specified time `since` or since an arbitrary processor-dependent time.
887 : !>
888 : !> \interface{getTimeCPU}
889 : !> \code{.F90}
890 : !>
891 : !> use pm_timer, only: getTimeCPU
892 : !> real(RKD) :: timeInSec
893 : !> logical(LK) :: failed
894 : !>
895 : !> timeInSec = getTimeCPU(since = since)
896 : !>
897 : !> \endcode
898 : !>
899 : !> \impure
900 : !>
901 : !> \see
902 : !> [timerCPU_type](@ref pm_timer::timerCPU_type)<br>
903 : !>
904 : !> \example{getTimeCPU}
905 : !> \include{lineno} example/pm_timer/getTimeCPU/main.F90
906 : !> \compilef{getTimeCPU}
907 : !> \output{getTimeCPU}
908 : !> \include{lineno} example/pm_timer/getTimeCPU/main.out.F90
909 : !>
910 : !> \test
911 : !> [test_pm_timer](@ref test_pm_timer)
912 : !>
913 : !> \finmain{getTimeCPU}
914 : !>
915 : !> \author
916 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
917 490772 : function getTimeCPU(since) result(timeInSec)
918 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
919 : !DEC$ ATTRIBUTES DLLEXPORT :: getTimeCPU
920 : #endif
921 : real(RKD), intent(in), optional :: since
922 : real(RKD) :: timeInSec
923 490772 : call cpu_time(timeInSec)
924 490772 : CHECK_ASSERTION(__LINE__, timeInSec >= 0._RKD, SK_"@getTimeCPU(): The CPU does not have a clock.") ! fpp
925 490772 : if (present(since)) timeInSec = timeInSec - since
926 490772 : end function
927 :
928 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
929 :
930 : !> \brief
931 : !> Generate and return the calendrical clock time in units of seconds since the
932 : !> specified time `since` or since the Gregorian calendrical time origin.
933 : !>
934 : !> \details
935 : !> See also the documentation details of [pm_timer](@ref pm_timer).
936 : !>
937 : !> \param[in] since : The input scalar of type `real` of kind double precision \RKD representing the time origin (epoch) in seconds.<br>
938 : !> (**optional**. If not present, its value is the Gregorian calendrical time origin.)
939 : !>
940 : !> \return
941 : !> `timeInSec` : The output scalar of type `real` of kind double precision \RKD,
942 : !> containing the time in seconds since the specified time origin.
943 : !>
944 : !> \interface{getTimeDAT}
945 : !> \code{.F90}
946 : !>
947 : !> use pm_timer, only: getTimeDAT
948 : !> real(RKD) :: timeInSec
949 : !>
950 : !> timeInSec = getTimeDAT(since = since)
951 : !>
952 : !> \endcode
953 : !>
954 : !> \impure
955 : !>
956 : !> \see
957 : !> [timerDAT_type](@ref pm_timer::timerDAT_type)<br>
958 : !>
959 : !> \example{getTimeDAT}
960 : !> \include{lineno} example/pm_timer/getTimeDAT/main.F90
961 : !> \compilef{getTimeDAT}
962 : !> \output{getTimeDAT}
963 : !> \include{lineno} example/pm_timer/getTimeDAT/main.out.F90
964 : !>
965 : !> \test
966 : !> [test_pm_timer](@ref test_pm_timer)
967 : !>
968 : !> \finmain{getTimeDAT}
969 : !>
970 : !> \author
971 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
972 118850 : function getTimeDAT(since) result(timeInSec)
973 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
974 : !DEC$ ATTRIBUTES DLLEXPORT :: getTimeDAT
975 : #endif
976 : use pm_kind, only: IK, RKD
977 : use pm_dateTime, only: getJulianDay
978 : use pm_dateTime, only: SECONDS_PER_DAY
979 : real(RKD), intent(in), optional :: since
980 : real(RKD) :: timeInSec
981 : integer(IK) :: values(8)
982 118850 : call date_and_time(values = values)
983 118850 : timeInSec = getJulianDay(values) * SECONDS_PER_DAY
984 118850 : if (present(since)) timeInSec = timeInSec - since
985 118850 : end function
986 :
987 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
988 :
989 : !> \brief
990 : !> Generate and return the MPI clock time in units of seconds since the
991 : !> specified time `since` or since an arbitrary time origin set by the MPI library.
992 : !>
993 : !> \details
994 : !> See also the documentation details of [pm_timer](@ref pm_timer).
995 : !>
996 : !> \param[in] since : The input scalar of type `real` of kind double precision \RKD representing the time origin (epoch) in seconds.<br>
997 : !> (**optional**. If not present, its value is an arbitrary time origin set by the MPI library.)
998 : !>
999 : !> \return
1000 : !> `timeInSec` : The output scalar of type `real` of kind double precision \RKD,
1001 : !> containing the time in seconds since the specified time origin.
1002 : !>
1003 : !> \interface{getTimeMPI}
1004 : !> \code{.F90}
1005 : !>
1006 : !> use pm_timer, only: getTimeMPI
1007 : !> real(RKD) :: timeInSec
1008 : !>
1009 : !> timeInSec = getTimeMPI(since = since)
1010 : !>
1011 : !> \endcode
1012 : !>
1013 : !> \warning
1014 : !> This procedure exists only if the library is built with the preprocessor macro `MPI_ENABLED=1`.
1015 : !>
1016 : !> \impure
1017 : !>
1018 : !> \see
1019 : !> [timerMPI_type](@ref pm_timer::timerMPI_type)<br>
1020 : !>
1021 : !> \example{getTimeMPI}
1022 : !> \include{lineno} example/pm_timer/getTimeMPI/main.F90
1023 : !> \compilef{getTimeMPI}
1024 : !> \output{getTimeMPI}
1025 : !> \include{lineno} example/pm_timer/getTimeMPI/main.out.F90
1026 : !>
1027 : !> \test
1028 : !> [test_pm_timer](@ref test_pm_timer)
1029 : !>
1030 : !> \finmain{getTimeMPI}
1031 : !>
1032 : !> \author
1033 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1034 : #if MPI_ENABLED
1035 : function getTimeMPI(since) result(timeInSec)
1036 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1037 : !DEC$ ATTRIBUTES DLLEXPORT :: getTimeMPI
1038 : #endif
1039 : use mpi !mpi_f08, only: mpi_wtime
1040 : use pm_kind, only: IK, RKD
1041 : real(RKD), intent(in), optional :: since
1042 : real(RKD) :: timeInSec
1043 : timeInSec = mpi_wtime()
1044 : if (present(since)) timeInSec = timeInSec - since
1045 : end function
1046 : #endif
1047 :
1048 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1049 :
1050 : !> \brief
1051 : !> Generate and return the OpenMP clock time in units of seconds since the
1052 : !> specified time `since` or since an arbitrary time origin set by the OpenMP library.
1053 : !>
1054 : !> \details
1055 : !> See also the documentation details of [pm_timer](@ref pm_timer).
1056 : !>
1057 : !> \param[in] since : The input scalar of type `real` of kind double precision \RKD representing the time origin (epoch) in seconds.<br>
1058 : !> (**optional**. If not present, its value is an arbitrary time origin set by the OpenMP library.)
1059 : !>
1060 : !> \return
1061 : !> `timeInSec` : The output scalar of type `real` of kind double precision \RKD,
1062 : !> containing the time in seconds since the specified time origin.
1063 : !>
1064 : !> \interface{getTimeOMP}
1065 : !> \code{.F90}
1066 : !>
1067 : !> use pm_timer, only: getTimeOMP
1068 : !> real(RKD) :: timeInSec
1069 : !>
1070 : !> timeInSec = getTimeOMP(since = since)
1071 : !>
1072 : !> \endcode
1073 : !>
1074 : !> \warning
1075 : !> This procedure exists only if the library is built with the preprocessor macro `OMP_ENABLED=1`.
1076 : !>
1077 : !> \impure
1078 : !>
1079 : !> \see
1080 : !> [timerOMP_type](@ref pm_timer::timerOMP_type)<br>
1081 : !>
1082 : !> \example{getTimeOMP}
1083 : !> \include{lineno} example/pm_timer/getTimeOMP/main.F90
1084 : !> \compilef{getTimeOMP}
1085 : !> \output{getTimeOMP}
1086 : !> \include{lineno} example/pm_timer/getTimeOMP/main.out.F90
1087 : !>
1088 : !> \test
1089 : !> [test_pm_timer](@ref test_pm_timer)
1090 : !>
1091 : !> \finmain{getTimeOMP}
1092 : !>
1093 : !> \author
1094 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1095 : #if OMP_ENABLED
1096 : function getTimeOMP(since) result(timeInSec)
1097 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1098 : !DEC$ ATTRIBUTES DLLEXPORT :: getTimeOMP
1099 : #endif
1100 : use pm_kind, only: IK, RKD
1101 : use omp_lib, only: omp_get_wtime
1102 : real(RKD), intent(in), optional :: since
1103 : real(RKD) :: timeInSec
1104 : timeInSec = omp_get_wtime()
1105 : if (present(since)) timeInSec = timeInSec - since
1106 : end function
1107 : #endif
1108 :
1109 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1110 :
1111 : !> \brief
1112 : !> Generate and return the system clock time in units of seconds since the
1113 : !> specified time `since` or since an arbitrary processor-dependent time origin.
1114 : !>
1115 : !> \details
1116 : !> See also the documentation details of [pm_timer](@ref pm_timer).
1117 : !>
1118 : !> \param[in] since : The input scalar of type `real` of kind double precision \RKD representing the time origin (epoch) in seconds.<br>
1119 : !> (**optional**. If not present, its value is processor-dependent but constant at runtime.)
1120 : !>
1121 : !> \return
1122 : !> `timeInSec` : The output scalar of type `real` of kind double precision \RKD,
1123 : !> containing the time in seconds since a processor-defined origin of time.
1124 : !>
1125 : !> \interface{getTimeSYS}
1126 : !> \code{.F90}
1127 : !>
1128 : !> use pm_timer, only: getTimeSYS
1129 : !> real(RKD) :: timeInSec
1130 : !>
1131 : !> timeInSec = getTimeSYS(since = since)
1132 : !>
1133 : !> \endcode
1134 : !>
1135 : !> \impure
1136 : !>
1137 : !> \see
1138 : !> [timerSYS_type](@ref pm_timer::timerSYS_type)<br>
1139 : !>
1140 : !> \example{getTimeSYS}
1141 : !> \include{lineno} example/pm_timer/getTimeSYS/main.F90
1142 : !> \compilef{getTimeSYS}
1143 : !> \output{getTimeSYS}
1144 : !> \include{lineno} example/pm_timer/getTimeSYS/main.out.F90
1145 : !>
1146 : !> \test
1147 : !> [test_pm_timer](@ref test_pm_timer)
1148 : !>
1149 : !> \todo
1150 : !> \pmed
1151 : !> The performance of this procedure could be improved by avoiding the unnecessary retrieval of
1152 : !> the count rate and computing its inverse to get the period (which is fixed in the entire runtime).<br>
1153 : !> This requires redefining this procedure as a type-bound procedure of a parent timer type.<br>
1154 : !> However, preliminary benchmarks indicate that the division (instead of multiplication) slows down
1155 : !> the computation by at most 25 percent. On modern architecture, the extra cost is likely on the
1156 : !> of 25 CPU cycles, approximately 25 nanoseconds. This is likely negligible in most practical
1157 : !> scenarios. Here is a benchmark script for the cost of division vs. multiplication,
1158 : !> \code{.F90}
1159 : !>
1160 : !> implicit none
1161 : !> double precision :: t0, t1, x(2), summ
1162 : !> integer :: v0(8), v1(8)
1163 : !> integer :: i
1164 : !>
1165 : !> call cpu_time(t0)
1166 : !> !call date_and_time(values = v0)
1167 : !> do i = 1, 10**7
1168 : !> call random_number(x)
1169 : !> summ = summ + x(2) * x(1)
1170 : !> end do
1171 : !> call cpu_time(t1)
1172 : !> print *, summ, "prod", t1 - t0
1173 : !> !call date_and_time(values = v1)
1174 : !> !print *, summ, "prod", v1(7) - v0(7) + (v1(8) - v0(8)) / 1.d3
1175 : !>
1176 : !> summ = 0.d0
1177 : !> call cpu_time(t0)
1178 : !> !call date_and_time(values = v0)
1179 : !> do i = 1, 10**7
1180 : !> call random_number(x)
1181 : !> summ = summ + x(2) / x(1)
1182 : !> end do
1183 : !> !call date_and_time(values = v1)
1184 : !> call cpu_time(t1)
1185 : !> print *, summ, t1 - t0
1186 : !> !print *, summ, "divi", v1(7) - v0(7) + (v1(8) - v0(8)) / 1.d3
1187 : !>
1188 : !> end
1189 : !>
1190 : !> \endcode
1191 : !> The above code yields highly similar time costs for multiplication vs. division operation in most benchmarks.
1192 : !>
1193 : !> \finmain{getTimeSYS}
1194 : !>
1195 : !> \author
1196 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1197 462309760 : function getTimeSYS(since) result(timeInSec)
1198 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1199 : !DEC$ ATTRIBUTES DLLEXPORT :: getTimeSYS
1200 : #endif
1201 : real(RKD), intent(in), optional :: since
1202 : real(RKD) :: timeInSec
1203 : real(RKD) :: count_rate
1204 : integer(IKD) :: count
1205 : !integer(IKD) :: cmax
1206 462309760 : call system_clock(count, count_rate)!, count_max = cmax)
1207 462309760 : CHECK_ASSERTION(__LINE__, count > 0_IKD .and. count_rate > 0._RKD, SK_"@getTimeSYS(): The system does not have a clock.") ! .and. cmax > 0_IKD
1208 462309760 : timeInSec = count / count_rate
1209 462309760 : if (present(since)) timeInSec = timeInSec - since
1210 462309760 : end function
1211 :
1212 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1213 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1214 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1215 :
1216 : !> \brief
1217 : !> Generate and return the time resolution by the build-time default timer of the ParaMonte library
1218 : !> as used in [timer_type](@ref pm_timer::timer_type), in units of seconds.
1219 : !>
1220 : !> \details
1221 : !> This function relies on either of the following lower-level module functions depending on the ParaMonte library build:<br>
1222 : !> <ol>
1223 : !> <li> [getResTimerMPI](@ref pm_timer::getResTimerMPI) for MPI-enabled parallel builds of the ParaMonte library.<br>
1224 : !> <li> [getResTimerOMP](@ref pm_timer::getResTimerOMP) for OpenMP-enabled parallel builds of the ParaMonte library.<br>
1225 : !> <li> [getResTimerSYS](@ref pm_timer::getResTimerSYS) for all other parallel or serial builds of the ParaMonte library.<br>
1226 : !> </ol>
1227 : !>
1228 : !> \return
1229 : !> `resolution` : The output scalar of type `real` of kind double precision \RKD,
1230 : !> containing the time resolution of the clock in units of seconds.
1231 : !>
1232 : !> \interface{getResTimer}
1233 : !> \code{.F90}
1234 : !>
1235 : !> use pm_timer, only: getResTimer
1236 : !> real(RKD) :: resolution
1237 : !>
1238 : !> resolution = getResTimer()
1239 : !>
1240 : !> \endcode
1241 : !>
1242 : !> \impure
1243 : !>
1244 : !> \see
1245 : !> [timer_type](@ref pm_timer::timer_type)<br>
1246 : !>
1247 : !> \example{getResTimer}
1248 : !> \include{lineno} example/pm_timer/getResTimer/main.F90
1249 : !> \compilef{getResTimer}
1250 : !> \output{getResTimer}
1251 : !> \include{lineno} example/pm_timer/getResTimer/main.out.F90
1252 : !>
1253 : !> \test
1254 : !> [test_pm_timer](@ref test_pm_timer)
1255 : !>
1256 : !> \finmain{getResTimer}
1257 : !>
1258 : !> \author
1259 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1260 1 : function getResTimer() result(resolution)
1261 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1262 : !DEC$ ATTRIBUTES DLLEXPORT :: getResTimer
1263 : #endif
1264 : real(RKD) :: resolution
1265 : #if MPI_ENABLED
1266 : resolution = getResTimerMPI()
1267 : #elif OMP_ENABLED
1268 : resolution = getResTimerOMP()
1269 : #else
1270 1 : resolution = getResTimerSYS()
1271 : #endif
1272 1 : end function
1273 :
1274 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1275 :
1276 : !> \brief
1277 : !> Generate and return the time resolution of the Fortran intrinsic `cpu_time()`,
1278 : !> used in [timerCPU_type](@ref pm_timer::timerCPU_type), in units of seconds.
1279 : !>
1280 : !> \details
1281 : !> See also the documentation details of [pm_timer](@ref pm_timer).
1282 : !>
1283 : !> \return
1284 : !> `resolution` : The output scalar of type `real` of kind double precision \RKD,
1285 : !> containing the time resolution of the CPU clock in units of seconds.
1286 : !>
1287 : !> \interface{getResTimerCPU}
1288 : !> \code{.F90}
1289 : !>
1290 : !> use pm_timer, only: getResTimerCPU
1291 : !> real(RKD) :: resolution
1292 : !>
1293 : !> resolution = getResTimerCPU()
1294 : !>
1295 : !> \endcode
1296 : !>
1297 : !> \impure
1298 : !>
1299 : !> \see
1300 : !> [getResTimerSYS](@ref pm_timer::getResTimerSYS)<br>
1301 : !>
1302 : !> \example{getResTimerCPU}
1303 : !> \include{lineno} example/pm_timer/getResTimerCPU/main.F90
1304 : !> \compilef{getResTimerCPU}
1305 : !> \output{getResTimerCPU}
1306 : !> \include{lineno} example/pm_timer/getResTimerCPU/main.out.F90
1307 : !>
1308 : !> \test
1309 : !> [test_pm_timer](@ref test_pm_timer)
1310 : !>
1311 : !> \finmain{getResTimerCPU}
1312 : !>
1313 : !> \author
1314 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1315 7 : function getResTimerCPU() result(resolution)
1316 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1317 : !DEC$ ATTRIBUTES DLLEXPORT :: getResTimerCPU
1318 : #endif
1319 : real(RKD) :: resolution
1320 : real(RKD) :: resolution_
1321 7 : call cpu_time(time = resolution_)
1322 7 : CHECK_ASSERTION(__LINE__, resolution_ >= 0._RKD, SK_"@getResTimerCPU(): The CPU does not have a clock.") ! fpp
1323 : do
1324 9 : call cpu_time(time = resolution)
1325 9 : if (resolution == resolution_) cycle
1326 7 : resolution = resolution - resolution_
1327 2 : exit
1328 : end do
1329 7 : end function
1330 :
1331 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1332 :
1333 : !> \brief
1334 : !> Generate and return the time resolution of the Fortran intrinsic `date_and_time()`,
1335 : !> used in [timerDAT_type](@ref pm_timer::timerDAT_type), in units of seconds.
1336 : !>
1337 : !> \details
1338 : !> See also the documentation details of [pm_timer](@ref pm_timer).
1339 : !>
1340 : !> \return
1341 : !> `resolution` : The output scalar of type `real` of kind double precision \RKD,
1342 : !> containing the time resolution of `date_and_time()` in units of seconds.
1343 : !>
1344 : !> \interface{getResTimerDAT}
1345 : !> \code{.F90}
1346 : !>
1347 : !> use pm_timer, only: getResTimerDAT
1348 : !> real(RKD) :: resolution
1349 : !>
1350 : !> resolution = getResTimerDAT()
1351 : !>
1352 : !> \endcode
1353 : !>
1354 : !> \pure
1355 : !>
1356 : !> \see
1357 : !> [getResTimerSYS](@ref pm_timer::getResTimerSYS)<br>
1358 : !>
1359 : !> \example{getResTimerDAT}
1360 : !> \include{lineno} example/pm_timer/getResTimerDAT/main.F90
1361 : !> \compilef{getResTimerDAT}
1362 : !> \output{getResTimerDAT}
1363 : !> \include{lineno} example/pm_timer/getResTimerDAT/main.out.F90
1364 : !>
1365 : !> \test
1366 : !> [test_pm_timer](@ref test_pm_timer)
1367 : !>
1368 : !> \finmain{getResTimerDAT}
1369 : !>
1370 : !> \author
1371 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1372 1 : pure function getResTimerDAT() result(resolution)
1373 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1374 : !DEC$ ATTRIBUTES DLLEXPORT :: getResTimerDAT
1375 : #endif
1376 : real(RKD) :: resolution
1377 : resolution = 0.001_RKD
1378 1 : end function
1379 :
1380 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1381 :
1382 : !> \brief
1383 : !> Generate and return the time resolution of the MPI intrinsic timer,
1384 : !> used in [timerMPI_type](@ref pm_timer::timerMPI_type), in units of seconds.
1385 : !>
1386 : !> \details
1387 : !> See also the documentation details of [pm_timer](@ref pm_timer).
1388 : !>
1389 : !> \return
1390 : !> `resolution` : The output scalar of type `real` of kind double precision \RKD,
1391 : !> containing the time resolution of `mpi_wtime()` in units of seconds.
1392 : !>
1393 : !> \interface{getResTimerMPI}
1394 : !> \code{.F90}
1395 : !>
1396 : !> use pm_timer, only: getResTimerMPI
1397 : !> real(RKD) :: resolution
1398 : !>
1399 : !> resolution = getResTimerMPI()
1400 : !>
1401 : !> \endcode
1402 : !>
1403 : !> \warning
1404 : !> This procedure exists only if the library is built with the preprocessor macro `MPI_ENABLED=1`.
1405 : !>
1406 : !> \see
1407 : !> [getResTimerSYS](@ref pm_timer::getResTimerSYS)<br>
1408 : !>
1409 : !> \example{getResTimerMPI}
1410 : !> \include{lineno} example/pm_timer/getResTimerMPI/main.F90
1411 : !> \compilef{getResTimerMPI}
1412 : !> \output{getResTimerMPI}
1413 : !> \include{lineno} example/pm_timer/getResTimerMPI/main.out.F90
1414 : !>
1415 : !> \test
1416 : !> [test_pm_timer](@ref test_pm_timer)
1417 : !>
1418 : !> \finmain{getResTimerMPI}
1419 : !>
1420 : !> \author
1421 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1422 : #if MPI_ENABLED
1423 : function getResTimerMPI() result(resolution)
1424 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1425 : !DEC$ ATTRIBUTES DLLEXPORT :: getResTimerMPI
1426 : #endif
1427 : use mpi !mpi_f08, only: mpi_wtick
1428 : real(RKD) :: resolution
1429 : resolution = real(mpi_wtick(), RKD)
1430 : end function
1431 : #endif
1432 :
1433 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1434 :
1435 : !> \brief
1436 : !> Generate and return the time resolution of the OpenMP intrinsic timer,
1437 : !> used in [timerOMP_type](@ref pm_timer::timerOMP_type), in units of seconds.
1438 : !>
1439 : !> \details
1440 : !> See also the documentation details of [pm_timer](@ref pm_timer).
1441 : !>
1442 : !> \return
1443 : !> `resolution` : The output scalar of type `real` of kind double precision \RKD,
1444 : !> containing the time resolution of `mpi_wtime()` in units of seconds.
1445 : !>
1446 : !> \interface{getResTimerOMP}
1447 : !> \code{.F90}
1448 : !>
1449 : !> use pm_timer, only: getResTimerOMP
1450 : !> real(RKD) :: resolution
1451 : !>
1452 : !> resolution = getResTimerOMP()
1453 : !>
1454 : !> \endcode
1455 : !>
1456 : !> \warning
1457 : !> This procedure exists only if the library is built with the preprocessor macro `OMP_ENABLED=1`.
1458 : !>
1459 : !> \see
1460 : !> [getResTimerSYS](@ref pm_timer::getResTimerSYS)<br>
1461 : !>
1462 : !> \example{getResTimerOMP}
1463 : !> \include{lineno} example/pm_timer/getResTimerOMP/main.F90
1464 : !> \compilef{getResTimerOMP}
1465 : !> \output{getResTimerOMP}
1466 : !> \include{lineno} example/pm_timer/getResTimerOMP/main.out.F90
1467 : !>
1468 : !> \test
1469 : !> [test_pm_timer](@ref test_pm_timer)
1470 : !>
1471 : !> \finmain{getResTimerOMP}
1472 : !>
1473 : !> \author
1474 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1475 : #if OMP_ENABLED
1476 : function getResTimerOMP() result(resolution)
1477 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1478 : !DEC$ ATTRIBUTES DLLEXPORT :: getResTimerOMP
1479 : #endif
1480 : use omp_lib, only: omp_get_wtick
1481 : real(RKD) :: resolution
1482 : resolution = real(omp_get_wtick(), RKD)
1483 : end function
1484 : #endif
1485 :
1486 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1487 :
1488 : !> \brief
1489 : !> Generate and return the time resolution of the Fortran intrinsic `system_clock()`,
1490 : !> used in [timerSYS_type](@ref pm_timer::timerSYS_type), in units of seconds.
1491 : !>
1492 : !> \details
1493 : !> See also the documentation details of [pm_timer](@ref pm_timer).
1494 : !>
1495 : !> \return
1496 : !> `resolution` : The output scalar of type `real` of kind double precision \RKD,
1497 : !> containing the time resolution of the system clock in units of seconds.
1498 : !>
1499 : !> \interface{getResTimerSYS}
1500 : !> \code{.F90}
1501 : !>
1502 : !> use pm_timer, only: getResTimerSYS
1503 : !> real(RKD) :: resolution
1504 : !>
1505 : !> resolution = getResTimerSYS()
1506 : !>
1507 : !> \endcode
1508 : !>
1509 : !> \impure
1510 : !>
1511 : !> \see
1512 : !> [timerSYS_type](@ref pm_timer::timerSYS_type)<br>
1513 : !>
1514 : !> \example{getResTimerSYS}
1515 : !> \include{lineno} example/pm_timer/getResTimerSYS/main.F90
1516 : !> \compilef{getResTimerSYS}
1517 : !> \output{getResTimerSYS}
1518 : !> \include{lineno} example/pm_timer/getResTimerSYS/main.out.F90
1519 : !>
1520 : !> \test
1521 : !> [test_pm_timer](@ref test_pm_timer)
1522 : !>
1523 : !> \finmain{getResTimerSYS}
1524 : !>
1525 : !> \author
1526 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1527 3434 : function getResTimerSYS() result(resolution)
1528 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1529 : !DEC$ ATTRIBUTES DLLEXPORT :: getResTimerSYS
1530 : #endif
1531 : real(RKD) :: resolution
1532 3434 : call system_clock(count_rate = resolution)
1533 3434 : CHECK_ASSERTION(__LINE__, resolution > 0._RKD, SK_"@getResTimerSYS(): The system does not have a clock.") ! fpp
1534 3434 : resolution = 1._RKD / resolution
1535 3434 : end function
1536 :
1537 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1538 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1539 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1540 :
1541 : !> \brief
1542 : !> Set the processor in sleep mode for the input requested number of seconds via the default timer of the ParaMonte library.
1543 : !>
1544 : !> \details
1545 : !> This function relies on either of the following lower-level module functions depending on the ParaMonte library build:<br>
1546 : !> <ol>
1547 : !> <li> [setIdleMPI](@ref pm_timer::setIdleMPI) for MPI-enabled parallel builds of the ParaMonte library.<br>
1548 : !> <li> [setIdleOMP](@ref pm_timer::setIdleOMP) for OpenMP-enabled parallel builds of the ParaMonte library.<br>
1549 : !> <li> [setIdleSYS](@ref pm_timer::setIdleSYS) for all other parallel or serial builds of the ParaMonte library.<br>
1550 : !> </ol>
1551 : !>
1552 : !> \param[in] seconds : The input scalar of type `real` of kind double precision \RKD
1553 : !> representing the number of seconds to put the system into sleep.
1554 : !>
1555 : !> \interface{setIdle}
1556 : !> \code{.F90}
1557 : !>
1558 : !> use pm_timer, only: setIdle
1559 : !> real(RKD) :: seconds
1560 : !>
1561 : !> call setIdle(seconds)
1562 : !>
1563 : !> \endcode
1564 : !>
1565 : !> \impure
1566 : !>
1567 : !> \see
1568 : !> [setIdleCPU](@ref pm_timer::setIdleCPU)<br>
1569 : !> [setIdleDAT](@ref pm_timer::setIdleDAT)<br>
1570 : !> [setIdleMPI](@ref pm_timer::setIdleMPI)<br>
1571 : !> [setIdleOMP](@ref pm_timer::setIdleOMP)<br>
1572 : !> [setIdleSYS](@ref pm_timer::setIdleSYS)<br>
1573 : !> [getTimeCPU](@ref pm_timer::getTimeCPU)<br>
1574 : !> [getTimeDAT](@ref pm_timer::getTimeDAT)<br>
1575 : !> [getTimeMPI](@ref pm_timer::getTimeMPI)<br>
1576 : !> [getTimeOMP](@ref pm_timer::getTimeOMP)<br>
1577 : !> [getTimeSYS](@ref pm_timer::getTimeSYS)<br>
1578 : !>
1579 : !> \example{setIdle}
1580 : !> \include{lineno} example/pm_timer/setIdle/main.F90
1581 : !> \compilef{setIdle}
1582 : !> \output{setIdle}
1583 : !> \include{lineno} example/pm_timer/setIdle/main.out.F90
1584 : !>
1585 : !> \test
1586 : !> [test_pm_timer](@ref test_pm_timer)
1587 : !>
1588 : !> \finmain{setIdle}
1589 : !>
1590 : !> \author
1591 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1592 2 : subroutine setIdle(seconds)
1593 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1594 : !DEC$ ATTRIBUTES DLLEXPORT :: setIdle
1595 : #endif
1596 : use pm_kind, only: RKD
1597 : real(RKD) , intent(in) :: seconds
1598 : #if MPI_ENABLED
1599 : call setIdleMPI(seconds)
1600 : #elif OMP_ENABLED
1601 : call setIdleOMP(seconds)
1602 : #else
1603 2 : call setIdleSYS(seconds)
1604 : #endif
1605 2 : end subroutine
1606 :
1607 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1608 :
1609 : !> \brief
1610 : !> Set the processor in sleep mode for the input requested number of seconds via the Fortran intrinsic `cpu_time()`.
1611 : !>
1612 : !> \param[in] seconds : The input scalar of type `real` of kind double precision \RKD
1613 : !> representing the number of seconds to put the system into sleep.
1614 : !>
1615 : !> \interface{setIdleCPU}
1616 : !> \code{.F90}
1617 : !>
1618 : !> use pm_timer, only: setIdleCPU
1619 : !> real(RKD) :: seconds
1620 : !>
1621 : !> call setIdleCPU(seconds)
1622 : !>
1623 : !> \endcode
1624 : !>
1625 : !> \warning
1626 : !> If the system does not possess a clock, the results are unpredictable and the procedure may enter an indefinite loop.<br>
1627 : !> However, this condition rarely occurs and can be readily check only once on each new platform before a production run.<br>
1628 : !> The lack of a system clock will lead to runtime error only if the library is built with the preprocessor macro `CHECK_ENABLED=1`.
1629 : !>
1630 : !> \impure
1631 : !>
1632 : !> \see
1633 : !> [setIdleCPU](@ref pm_timer::setIdleCPU)<br>
1634 : !> [setIdleDAT](@ref pm_timer::setIdleDAT)<br>
1635 : !> [setIdleMPI](@ref pm_timer::setIdleMPI)<br>
1636 : !> [setIdleOMP](@ref pm_timer::setIdleOMP)<br>
1637 : !> [setIdleSYS](@ref pm_timer::setIdleSYS)<br>
1638 : !> [getTimeCPU](@ref pm_timer::getTimeCPU)<br>
1639 : !> [getTimeDAT](@ref pm_timer::getTimeDAT)<br>
1640 : !> [getTimeMPI](@ref pm_timer::getTimeMPI)<br>
1641 : !> [getTimeOMP](@ref pm_timer::getTimeOMP)<br>
1642 : !> [getTimeSYS](@ref pm_timer::getTimeSYS)<br>
1643 : !>
1644 : !> \example{setIdleCPU}
1645 : !> \include{lineno} example/pm_timer/setIdleCPU/main.F90
1646 : !> \compilef{setIdleCPU}
1647 : !> \output{setIdleCPU}
1648 : !> \include{lineno} example/pm_timer/setIdleCPU/main.out.F90
1649 : !>
1650 : !> \test
1651 : !> [test_pm_timer](@ref test_pm_timer)
1652 : !>
1653 : !> \finmain{setIdleCPU}
1654 : !>
1655 : !> \author
1656 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1657 4 : subroutine setIdleCPU(seconds)
1658 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1659 : !DEC$ ATTRIBUTES DLLEXPORT :: setIdleCPU
1660 : #endif
1661 : use pm_kind, only: RKD
1662 : real(RKD) , intent(in) :: seconds
1663 : real(RKD) :: since
1664 4 : since = getTimeCPU()
1665 : do
1666 439264 : if (getTimeCPU() - since > seconds) exit
1667 : end do
1668 4 : end subroutine
1669 :
1670 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1671 :
1672 : !> \brief
1673 : !> Set the processor in sleep mode for the input requested number of seconds via the Fortran intrinsic `date_and_time()`.
1674 : !>
1675 : !> \param[in] seconds : The input scalar of type `real` of kind double precision \RKD
1676 : !> representing the number of seconds to put the system into sleep.
1677 : !>
1678 : !> \interface{setIdleDAT}
1679 : !> \code{.F90}
1680 : !>
1681 : !> use pm_timer, only: setIdleDAT
1682 : !> real(RKD) :: seconds
1683 : !>
1684 : !> call setIdleDAT(seconds)
1685 : !>
1686 : !> \endcode
1687 : !>
1688 : !> \warning
1689 : !> If the system does not possess a clock, the results are unpredictable and the procedure may enter an indefinite loop.<br>
1690 : !> However, this condition rarely occurs and can be readily check only once on each new platform before a production run.<br>
1691 : !> The lack of a system clock will lead to runtime error only if the library is built with the preprocessor macro `CHECK_ENABLED=1`.
1692 : !>
1693 : !> \impure
1694 : !>
1695 : !> \see
1696 : !> [setIdleCPU](@ref pm_timer::setIdleCPU)<br>
1697 : !> [setIdleDAT](@ref pm_timer::setIdleDAT)<br>
1698 : !> [setIdleMPI](@ref pm_timer::setIdleMPI)<br>
1699 : !> [setIdleOMP](@ref pm_timer::setIdleOMP)<br>
1700 : !> [setIdleSYS](@ref pm_timer::setIdleSYS)<br>
1701 : !> [getTimeCPU](@ref pm_timer::getTimeCPU)<br>
1702 : !> [getTimeDAT](@ref pm_timer::getTimeDAT)<br>
1703 : !> [getTimeMPI](@ref pm_timer::getTimeMPI)<br>
1704 : !> [getTimeOMP](@ref pm_timer::getTimeOMP)<br>
1705 : !> [getTimeSYS](@ref pm_timer::getTimeSYS)<br>
1706 : !>
1707 : !> \example{setIdleDAT}
1708 : !> \include{lineno} example/pm_timer/setIdleDAT/main.F90
1709 : !> \compilef{setIdleDAT}
1710 : !> \output{setIdleDAT}
1711 : !> \include{lineno} example/pm_timer/setIdleDAT/main.out.F90
1712 : !>
1713 : !> \test
1714 : !> [test_pm_timer](@ref test_pm_timer)
1715 : !>
1716 : !> \finmain{setIdleDAT}
1717 : !>
1718 : !> \author
1719 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1720 4 : subroutine setIdleDAT(seconds)
1721 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1722 : !DEC$ ATTRIBUTES DLLEXPORT :: setIdleDAT
1723 : #endif
1724 : use pm_kind, only: RKD
1725 : real(RKD) , intent(in) :: seconds
1726 : real(RKD) :: since
1727 4 : since = getTimeDAT()
1728 : do
1729 118230 : if (getTimeDAT() - since > seconds) exit
1730 : end do
1731 4 : end subroutine
1732 :
1733 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1734 :
1735 : !> \brief
1736 : !> Set the processor in sleep mode for the input requested number of seconds via the MPI intrinsic timer `MPI_Wtime()`.
1737 : !>
1738 : !> \param[in] seconds : The input scalar of type `real` of kind double precision \RKD
1739 : !> representing the number of seconds to put the system into sleep.
1740 : !>
1741 : !> \interface{setIdleMPI}
1742 : !> \code{.F90}
1743 : !>
1744 : !> use pm_timer, only: setIdleMPI
1745 : !> real(RKD) :: seconds
1746 : !>
1747 : !> call setIdleMPI(seconds)
1748 : !>
1749 : !> \endcode
1750 : !>
1751 : !> \warning
1752 : !> This procedure exists only if the library is built with the preprocessor macro `MPI_ENABLED=1`.<br>
1753 : !>
1754 : !> \warning
1755 : !> Calling this procedure requires the MPI library to have been initialized and not finalized.<br>
1756 : !> \vericon
1757 : !>
1758 : !> \impure
1759 : !>
1760 : !> \see
1761 : !> [setIdleCPU](@ref pm_timer::setIdleCPU)<br>
1762 : !> [setIdleDAT](@ref pm_timer::setIdleDAT)<br>
1763 : !> [setIdleMPI](@ref pm_timer::setIdleMPI)<br>
1764 : !> [setIdleOMP](@ref pm_timer::setIdleOMP)<br>
1765 : !> [setIdleSYS](@ref pm_timer::setIdleSYS)<br>
1766 : !> [getTimeCPU](@ref pm_timer::getTimeCPU)<br>
1767 : !> [getTimeDAT](@ref pm_timer::getTimeDAT)<br>
1768 : !> [getTimeMPI](@ref pm_timer::getTimeMPI)<br>
1769 : !> [getTimeOMP](@ref pm_timer::getTimeOMP)<br>
1770 : !> [getTimeSYS](@ref pm_timer::getTimeSYS)<br>
1771 : !>
1772 : !> \example{setIdleMPI}
1773 : !> \include{lineno} example/pm_timer/setIdleMPI/main.F90
1774 : !> \compilef{setIdleMPI}
1775 : !> \output{setIdleMPI}
1776 : !> \include{lineno} example/pm_timer/setIdleMPI/main.out.F90
1777 : !>
1778 : !> \test
1779 : !> [test_pm_timer](@ref test_pm_timer)
1780 : !>
1781 : !> \finmain{setIdleMPI}
1782 : !>
1783 : !> \author
1784 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1785 : #if MPI_ENABLED
1786 : subroutine setIdleMPI(seconds)
1787 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1788 : !DEC$ ATTRIBUTES DLLEXPORT :: setIdleMPI
1789 : #endif
1790 : use pm_kind, only: RKD
1791 : real(RKD) , intent(in) :: seconds
1792 : real(RKD) :: since
1793 : since = getTimeMPI()
1794 : do
1795 : if (getTimeMPI() - since > seconds) exit
1796 : end do
1797 : end subroutine
1798 : #endif
1799 :
1800 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1801 :
1802 : !> \brief
1803 : !> Set the processor in sleep mode for the input requested number of seconds via the OpenMP intrinsic timer `omp_get_wtime()`.
1804 : !>
1805 : !> \param[in] seconds : The input scalar of type `real` of kind double precision \RKD
1806 : !> representing the number of seconds to put the system into sleep.
1807 : !>
1808 : !> \interface{setIdleOMP}
1809 : !> \code{.F90}
1810 : !>
1811 : !> use pm_timer, only: setIdleOMP
1812 : !> real(RKD) :: seconds
1813 : !>
1814 : !> call setIdleOMP(seconds)
1815 : !>
1816 : !> \endcode
1817 : !>
1818 : !> \warning
1819 : !> This procedure exists only if the library is built with the preprocessor macro `OMP_ENABLED=1`.<br>
1820 : !>
1821 : !> \impure
1822 : !>
1823 : !> \see
1824 : !> [setIdleCPU](@ref pm_timer::setIdleCPU)<br>
1825 : !> [setIdleDAT](@ref pm_timer::setIdleDAT)<br>
1826 : !> [setIdleMPI](@ref pm_timer::setIdleMPI)<br>
1827 : !> [setIdleOMP](@ref pm_timer::setIdleOMP)<br>
1828 : !> [setIdleSYS](@ref pm_timer::setIdleSYS)<br>
1829 : !> [getTimeCPU](@ref pm_timer::getTimeCPU)<br>
1830 : !> [getTimeDAT](@ref pm_timer::getTimeDAT)<br>
1831 : !> [getTimeMPI](@ref pm_timer::getTimeMPI)<br>
1832 : !> [getTimeOMP](@ref pm_timer::getTimeOMP)<br>
1833 : !> [getTimeSYS](@ref pm_timer::getTimeSYS)<br>
1834 : !>
1835 : !> \example{setIdleOMP}
1836 : !> \include{lineno} example/pm_timer/setIdleOMP/main.F90
1837 : !> \compilef{setIdleOMP}
1838 : !> \output{setIdleOMP}
1839 : !> \include{lineno} example/pm_timer/setIdleOMP/main.out.F90
1840 : !>
1841 : !> \test
1842 : !> [test_pm_timer](@ref test_pm_timer)
1843 : !>
1844 : !> \finmain{setIdleOMP}
1845 : !>
1846 : !> \author
1847 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1848 : #if OMP_ENABLED
1849 : subroutine setIdleOMP(seconds)
1850 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1851 : !DEC$ ATTRIBUTES DLLEXPORT :: setIdleOMP
1852 : #endif
1853 : use pm_kind, only: RKD
1854 : real(RKD) , intent(in) :: seconds
1855 : real(RKD) :: since
1856 : since = getTimeOMP()
1857 : do
1858 : if (getTimeOMP() - since > seconds) exit
1859 : end do
1860 : end subroutine
1861 : #endif
1862 :
1863 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1864 :
1865 : !> \brief
1866 : !> Set the processor in sleep mode for the input requested number of seconds via the Fortran intrinsic `system_clock()`.
1867 : !>
1868 : !> \param[in] seconds : The input scalar of type `real` of kind double precision \RKD
1869 : !> representing the number of seconds to put the system into sleep.
1870 : !>
1871 : !> \interface{setIdleSYS}
1872 : !> \code{.F90}
1873 : !>
1874 : !> use pm_timer, only: setIdleSYS
1875 : !> real(RKD) :: seconds
1876 : !>
1877 : !> call setIdleSYS(seconds)
1878 : !>
1879 : !> \endcode
1880 : !>
1881 : !> \warning
1882 : !> If the system does not possess a clock, the results are unpredictable and the procedure may enter an indefinite loop.<br>
1883 : !> However, this condition rarely occurs and can be readily check only once on each new platform before a production run.<br>
1884 : !> The lack of a system clock will lead to runtime error only if the library is built with the preprocessor macro `CHECK_ENABLED=1`.
1885 : !>
1886 : !> \impure
1887 : !>
1888 : !> \see
1889 : !> [setIdleCPU](@ref pm_timer::setIdleCPU)<br>
1890 : !> [setIdleDAT](@ref pm_timer::setIdleDAT)<br>
1891 : !> [setIdleMPI](@ref pm_timer::setIdleMPI)<br>
1892 : !> [setIdleOMP](@ref pm_timer::setIdleOMP)<br>
1893 : !> [setIdleSYS](@ref pm_timer::setIdleSYS)<br>
1894 : !> [getTimeCPU](@ref pm_timer::getTimeCPU)<br>
1895 : !> [getTimeDAT](@ref pm_timer::getTimeDAT)<br>
1896 : !> [getTimeMPI](@ref pm_timer::getTimeMPI)<br>
1897 : !> [getTimeOMP](@ref pm_timer::getTimeOMP)<br>
1898 : !> [getTimeSYS](@ref pm_timer::getTimeSYS)<br>
1899 : !>
1900 : !> \example{setIdleSYS}
1901 : !> \include{lineno} example/pm_timer/setIdleSYS/main.F90
1902 : !> \compilef{setIdleSYS}
1903 : !> \output{setIdleSYS}
1904 : !> \include{lineno} example/pm_timer/setIdleSYS/main.out.F90
1905 : !>
1906 : !> \test
1907 : !> [test_pm_timer](@ref test_pm_timer)
1908 : !>
1909 : !> \finmain{setIdleSYS}
1910 : !>
1911 : !> \author
1912 : !> \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
1913 2419128 : subroutine setIdleSYS(seconds)
1914 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1915 : !DEC$ ATTRIBUTES DLLEXPORT :: setIdleSYS
1916 : #endif
1917 : use pm_kind, only: RKD
1918 : real(RKD) , intent(in) :: seconds
1919 : real(RKD) :: since
1920 2419128 : since = getTimeSYS()
1921 : do
1922 444217156 : if (getTimeSYS() - since > seconds) exit
1923 : end do
1924 2419128 : end subroutine
1925 :
1926 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1927 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1928 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1929 :
1930 : end module pm_timer ! LCOV_EXCL_LINE
|