]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliCodeTimer.cxx
- Change in ProcessFloatAll method for calculating statistical quantities for
[u/mrichter/AliRoot.git] / STEER / AliCodeTimer.cxx
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
7 * Permission to use, copy, modify and distribute this software and its   *
8 * documentation strictly for non-commercial purposes is hereby granted   *
9 * without fee, provided that the above copyright notice appears in all   *
10 * copies and that both the copyright notice and this permission notice   *
11 * appear in the supporting documentation. The authors make no claims     *
12 * about the suitability of this software for any purpose. It is          *
13 * provided "as is" without express or implied warranty.                  *
14 **************************************************************************/
15 //  $Id$
16
17 //_________________________________________________________________________
18 // Class to get organized with the way we're timing our methods...
19 //
20 // Typical usage is based on macros (like for AliLog related ones AliDebug...)
21 //
22 // The idea is to instrument the code with a few macro calls, and then,
23 // at the end of the execution, get a printout of *all* the timers, by using
24 // AliCodeTimer::Instance()->Print()
25 // instead of getting scattered outputs all over the place.
26 //
27 // To time a given method, use :
28 //
29 // void ClassA::MethodA(....)
30 // {
31 //    AliCodeTimerAuto("")
32 // }
33 //
34 // To get several timers within a same method, use : 
35 //
36 // void ClassA::MethodB(...)
37 // {
38 //   AliCodeTimerStart("doing something")
39 //   ....
40 //   AliCodeTimerStop("doing something")
41 //
42 //   AliCodeTimerStart("doing something else")
43 //   ....
44 //   AliCodeTimerStop("doing something else")
45 // }
46
47 #include "AliCodeTimer.h"
48
49 #include <TMap.h>
50 #include <TObjString.h>
51 #include <TStopwatch.h>
52 #include <Riostream.h>
53
54 /// \cond CLASSIMP
55 ClassImp(AliCodeTimer)
56 ClassImp(AliCodeTimer::AliPair)
57 /// \endcond
58
59 AliCodeTimer* AliCodeTimer::fgInstance(0x0);
60
61 //_____________________________________________________________________________
62 void
63 AliCodeTimer::AliPair::Print(Option_t* opt) const
64 {
65   // Print timer information
66   cout << opt << Form("%s R:%.4fs C:%.4fs (%d slices)",
67                       Name().Data(),Timer()->RealTime(),
68                       Timer()->CpuTime(),Timer()->Counter()-1) << endl;
69 }
70
71
72 //_____________________________________________________________________________
73 AliCodeTimer::AliCodeTimer() : TObject(), fTimers(new TMap)
74 {
75   /// Ctor
76   fTimers->SetOwner(kTRUE);
77 }
78
79 //_____________________________________________________________________________
80 AliCodeTimer::~AliCodeTimer()
81 {
82   /// Dtor
83   Reset();
84   delete fTimers;
85 }
86
87 //_____________________________________________________________________________
88 AliCodeTimer*
89 AliCodeTimer::Instance()
90 {
91   // single instance of this class
92   if (!fgInstance) fgInstance = new AliCodeTimer;
93   return fgInstance;
94 }
95
96 //_____________________________________________________________________________
97 void AliCodeTimer::Continue(const char* classname, const char* methodname, 
98                             const char* message)
99 {
100   /// Resume a previously stop timer
101   TStopwatch* t = Stopwatch(classname,methodname,message);
102   if (t)
103   {
104     t->Continue();
105   }
106   else
107   {
108     AliError(Form("No timer for %s/%s/%s",classname,methodname,message));
109   }
110 }
111
112 //_____________________________________________________________________________
113 Double_t AliCodeTimer::CpuTime(const char* classname, 
114                                const char* methodname,
115                                const char* message) const
116 {
117   /// Return cpu time for a given timer
118   TStopwatch* t = Stopwatch(classname,methodname,message);
119   if (t)
120   {
121     return t->CpuTime();
122   }
123   else
124   {
125     return 0;
126   }
127 }
128
129 //_____________________________________________________________________________
130 TMap*
131 AliCodeTimer::MethodMap(const char* classname) const
132 {
133   /// Return the map for a given "classname"
134   return static_cast<TMap*>(fTimers->GetValue(classname));
135 }
136
137 //_____________________________________________________________________________
138 TObjArray*
139 AliCodeTimer::MessageArray(const char* classname, const char* methodname) const
140 {
141   /// Return the array for a given AliPair (classname,methodname)
142   TMap* m = MethodMap(classname);
143   if ( m ) 
144   {
145     return static_cast<TObjArray*>(m->GetValue(methodname));
146   }
147   return 0;
148 }
149
150 //_____________________________________________________________________________
151 void AliCodeTimer::PrintMethod(const char* classname, const char* methodname) const
152 {
153   /// Print all the timers for a given method
154   TObjArray* messages = MessageArray(classname,methodname);
155   messages->Sort();
156   
157   cout << "   " << methodname << " ";
158   
159   if ( messages->GetLast() == 0 ) 
160   {
161     AliPair* p = static_cast<AliPair*>(messages->First());
162     p->Print();
163   }
164   else
165   {
166     cout << endl;
167     
168     TIter next(messages);
169     AliPair* p;
170   
171     while ( ( p = static_cast<AliPair*>(next()) ) ) 
172     {
173       p->Print("        ");
174     }   
175   }
176 }
177
178 //_____________________________________________________________________________
179 void AliCodeTimer::PrintClass(const char* classname) const
180 {
181   /// Print all the timers for a given class
182   TMap* methods = MethodMap(classname);
183   TIter next(methods);
184   TObjString* methodname;
185   TObjArray methodNameArray;
186   
187   while ( ( methodname = static_cast<TObjString*>(next()) ) ) 
188   {
189     methodNameArray.Add(methodname);
190   }
191   
192   cout << classname << endl;
193   
194   methodNameArray.Sort();
195   
196   TIter mnext(&methodNameArray);
197   
198   while ( ( methodname = static_cast<TObjString*>(mnext()) ) ) 
199   {
200     PrintMethod(classname,methodname->String().Data());
201   }
202 }
203   
204 //_____________________________________________________________________________
205 void AliCodeTimer::Print(Option_t* /*opt*/) const
206 {
207   /// Print all the timers we hold
208   TIter next(fTimers);
209   TObjString* classname;
210   TObjArray classNameArray;
211   
212   while ( ( classname = static_cast<TObjString*>(next()) ) )
213   {
214     classNameArray.Add(classname);
215   }
216   
217   classNameArray.Sort();
218   
219   TIter cnext(&classNameArray);
220   while ( ( classname = static_cast<TObjString*>(cnext()) ) )
221   {
222     PrintClass(classname->String().Data());
223   }
224 }
225
226 //_____________________________________________________________________________
227 Double_t 
228 AliCodeTimer::RealTime(const char* classname, const char* methodname,
229                        const char* message) const
230 {
231   /// Return real time of a given time
232   TStopwatch* t = Stopwatch(classname,methodname,message);
233   if (t)
234   {
235     return t->RealTime();
236   }
237   else
238   {
239     return 0;
240   }
241 }
242
243 //_____________________________________________________________________________
244 void
245 AliCodeTimer::Reset()
246 {
247   /// Reset
248   TIter next(fTimers);
249   TObjString* classname;
250   
251   while ( ( classname = static_cast<TObjString*>(next()) ) ) 
252   {
253     TMap* m = static_cast<TMap*>(fTimers->GetValue(classname->String().Data()));
254     m->DeleteAll();
255   }
256   
257   fTimers->DeleteAll();
258 }
259
260 //_____________________________________________________________________________
261 void 
262 AliCodeTimer::Start(const char* classname, const char* methodname,
263                     const char* message)
264 {
265   /// Start a given time
266   TStopwatch* t = Stopwatch(classname,methodname,message);
267   if (!t)
268   {
269     TMap* m = MethodMap(classname);
270     if (!m)
271     {
272       m = new TMap;
273       m->SetOwner(kTRUE);
274       fTimers->Add(new TObjString(classname),m);
275     }      
276     TObjArray* messages = MessageArray(classname,methodname);
277     if (!messages)
278     {
279       messages = new TObjArray;
280       messages->SetOwner(kTRUE);
281       m->Add(new TObjString(methodname),messages);
282     }
283     t = new TStopwatch;
284     t->Start(kTRUE);
285     t->Stop();
286     messages->Add(new AliPair(new TObjString(message),t));
287   }
288   t->Start(kFALSE);
289 }
290
291 //_____________________________________________________________________________
292 void 
293 AliCodeTimer::Stop(const char* classname, const char* methodname,
294                    const char* message)
295 {
296   /// Stop a given timer
297   TStopwatch* t = Stopwatch(classname,methodname,message);
298   if (!t)
299   {
300     AliError(Form("No timer for %s/%s/%s",classname,methodname,message));
301   }
302   else
303   {
304     t->Stop();
305   }
306 }
307
308 //_____________________________________________________________________________
309 TStopwatch* 
310 AliCodeTimer::Stopwatch(const char* classname, const char* methodname,
311                         const char* message) const
312 {
313   /// Return the internal TStopwatch for a given timer
314   TObjArray* a = MessageArray(classname,methodname);
315   if ( a ) 
316   {
317     if (message)
318     {
319       TIter next(a);
320       AliPair* p;
321       while ( ( p = static_cast<AliPair*>(next()) ) ) 
322       {
323         TString s = p->Name();
324         if ( s == TString(message) ) 
325         {
326           return p->Timer();
327         }
328       }
329     }
330     else
331     {
332       return static_cast<TStopwatch*>(a->First());
333     }
334   }
335   return 0x0;
336 }