source: timer.F90

Last change on this file was 3fbbfbb, checked in by Jan Meinke <j.meinke@…>, 14 years ago

Move to doxygen comments and smmp_p.

Doxygen comments in Fortran are !> ... !! ... !<. I'm planning move the API documentation from the
lyx file into the code. This should make it easier to get documentation for all the common block
variables as well.

Use import smmp_p to indicate the parallel version of the Python bindings.

  • Property mode set to 100644
File size: 4.7 KB
Line 
1! timer.f90
2!
3! (c)2007 Jan H. Meinke
4!
5#ifdef BGL
6module rts
7! interface definition of rts_get_timebase, Fortran2003, available in XLF Version 9
8 interface
9 function rts_get_timebase() bind(c)
10 use, intrinsic :: iso_c_binding
11 integer(c_long_long) :: rts_get_timebase
12 end function rts_get_timebase
13 end interface
14end module rts
15#endif
16
17!! This module contains a minimal set of routines to time sections of code.
18! It can keep track of several timers at once. Each timer is associated with
19! an integer ID. Individual timing events are not stored.
20!
21! Currently only timerIDs between 1 and 10 are valid.
22module timer
23#ifdef BGL
24 use rts
25 use, intrinsic :: iso_c_binding
26#endif
27 implicit none
28
29 !! Contains the minimum time ever measured, the maximum time ever measured,
30 ! the average time, the standard deviation, and the number of measurements.
31 type TimingStructure
32 integer :: counter
33#ifdef BGL
34 integer(c_long_long) :: t0, t1
35#else
36 real(kind=8) :: t0, t1
37#endif
38 real(kind=8) :: minTime, maxTime
39 real(kind=8) :: sumOfTimes, sumOfSquares
40 real(kind=8) :: average, sigma
41 end type TimingStructure
42
43 !! Contains the timing data for all timers.
44 type(TimingStructure), allocatable :: timingData(:)
45
46 integer :: initialTimerCount = 10
47#ifdef BGL
48 real(kind=8) :: clockspeed = 700.0d6
49#endif
50
51 contains
52 !! Initializes timer
53 ! This routine must be called before starting any timer.
54 !
55 ! @author Jan H. Meinke
56 subroutine init_timer()
57 integer :: i
58
59 ! Initialize timers
60 allocate (timingData(initialTimerCount))
61! forall(i = 1:initialTimerCount)
62 do i = 1, initialTimerCount
63 timingData(i) = TimingStructure(0, -1.0, -1.0, -1.0, -1.0, 0.0, &
64 0.0, -1.0, -1.0)
65! end forall
66 end do
67 end subroutine init_timer
68
69 !! Start a new measurement.
70 !
71 ! @param id ID of the timer
72 !
73 ! @author Jan H. Meinke
74 subroutine start_timer(id)
75 integer, intent(in) :: id
76 timingData(id)%t0 = getCurrentTime()
77 end subroutine start_timer
78
79 !! Take a measurement and stop the timer with ID id
80 !
81 ! @param id ID of the timer
82 !
83 ! @author Jan H. Meinke
84 subroutine stop_timer(id)
85 integer, intent(in) :: id
86 real(kind=8) :: dt
87
88 timingData(id)%t1 = getCurrentTime()
89#ifdef BGL
90 dt = (timingData(id)%t1 - timingData(id)%t0) / clockspeed
91#else
92 dt = (timingData(id)%t1 - timingData(id)%t0)
93#endif
94 timingData(id)%counter = timingData(id)%counter + 1
95 timingData(id)%maxTime = max(dt, timingData(id)%maxTime)
96 if (timingData(id)%minTime > 0) then
97 timingData(id)%minTime = min(dt, timingData(id)%minTime)
98 else
99 timingData(id)%minTime = dt
100 end if
101 timingData(id)%sumOfTimes = timingData(id)%sumOfTimes + dt
102 timingData(id)%sumOfSquares = timingData(id)%sumOfSquares + dt ** 2
103
104 end subroutine stop_timer
105
106 subroutine evaluate(id)
107 integer, intent(in) :: id
108 if (timingData(id)%counter > 0) then
109 timingData(id)%average = timingData(id)%sumOfTimes / timingData(id)%counter
110 end if
111 if (timingData(id)%counter > 1) then
112! write (*,*) timingData(id)%average, timingData(id)%sumOfSquares,&
113! timingData(id)%counter
114 timingData(id)%sigma = sqrt((abs(timingData(id)%average ** 2 &
115 - timingData(id)%sumOfSquares / timingData(id)%counter)) &
116 / (timingData(id)%counter - 1))
117 end if
118 end subroutine evaluate
119
120 !! Wrapper around the timing function used.
121 ! If we are using MPI, we can use MPI_WTime, otherwise we use cpu_time(time)
122#ifdef BGL
123 integer(c_long_long) function getCurrentTime()
124 integer(c_long_long) :: time
125 time = rts_get_timebase()
126 getCurrentTime = time
127 end function getCurrentTime
128#elif MPI
129 real(kind=8) function getCurrentTime()
130 real(kind=8) :: time
131 include 'mpif.h'
132 time = MPI_Wtime()
133 getCurrentTime = time
134 end function getCurrentTime
135#else
136 real(kind=8) function getCurrentTime()
137 real(kind=8) :: time
138 call cpu_time(time)
139 getCurrentTime = time
140 end function getCurrentTime
141#endif
142
143end module timer
Note: See TracBrowser for help on using the repository browser.