]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/BASE/AliHLTTTreeProcessor.cxx
skipping publishing if there was no new event added during the last interval
[u/mrichter/AliRoot.git] / HLT / BASE / AliHLTTTreeProcessor.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Timur Pocheptsov <Timur.Pocheptsov@cern.ch>           *
8 //*                  Matthias Richter <Matthias.Richter@cern.ch>
9 //*                  for The ALICE HLT Project.                            *
10 //*                                                                        *
11 //* Permission to use, copy, modify and distribute this software and its   *
12 //* documentation strictly for non-commercial purposes is hereby granted   *
13 //* without fee, provided that the above copyright notice appears in all   *
14 //* copies and that both the copyright notice and this permission notice   *
15 //* appear in the supporting documentation. The authors make no claims     *
16 //* about the suitability of this software for any purpose. It is          *
17 //* provided "as is" without express or implied warranty.                  *
18 //**************************************************************************
19
20 /// @file   AliHLTTTreeProcessor.cxx
21 /// @author Timur Pocheptsov, Matthias Richter
22 /// @date   05.07.2010
23 /// @brief  Generic component for data collection in a TTree
24
25 #include <cerrno>
26 #include <memory>
27
28 #include "AliHLTTTreeProcessor.h"
29 #include "TDirectory.h"
30 #include "TDatime.h"
31 #include "TString.h"
32 #include "TTree.h"
33 #include "TH1.h"
34 #include "TStopwatch.h"
35
36 /** ROOT macro for the implementation of ROOT specific class methods */
37 ClassImp(AliHLTTTreeProcessor)
38
39 AliHLTTTreeProcessor::AliHLTTTreeProcessor()
40                         : AliHLTProcessor(), 
41                           fDefinitions(),
42                           fTree(0),
43                           fMaxEntries(kMaxEntries),
44                           fPublishInterval(kInterval),
45                           fLastTime(0),
46                           fpEventTimer(NULL),
47                           fpCycleTimer(NULL),
48                           fMaxEventTime(0),
49                           fNofEventsForce(0),
50                           fForcedEventsCount(0),
51                           fSkippedEventsCount(0),
52                           fNewEventsCount(0)
53 {
54   // see header file for class documentation
55   // or
56   // refer to README to build package
57   // or
58   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
59 }
60
61 const AliHLTUInt32_t AliHLTTTreeProcessor::fgkTimeScale=1000000; // ticks per second
62
63 AliHLTTTreeProcessor::~AliHLTTTreeProcessor()
64 {
65   // see header file for class documentation
66 }
67
68 AliHLTComponentDataType AliHLTTTreeProcessor::GetOutputDataType()
69 {
70   // get the component output data type
71   return kAliHLTDataTypeHistogram;
72 }
73
74 void AliHLTTTreeProcessor::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
75 {
76   // get the output size estimator
77   //
78   if (!fDefinitions.size()) {
79     HLTError("Can not calculate output data size, no histogram definitions were provided");
80     return;
81   }
82
83   constBase = 0;
84   for (list_const_iterator i = fDefinitions.begin(); i != fDefinitions.end(); ++i)
85     constBase += i->GetSize();
86
87   inputMultiplier = 1.;
88 }
89
90 int AliHLTTTreeProcessor::DoInit(int argc, const char** argv)
91 {
92   // init component
93   // ask child to create the tree.
94   int iResult = 0;
95
96   if (!fTree) {
97     std::auto_ptr<TTree> ptr(CreateTree(argc, argv));
98     if (ptr.get()) {
99       //Stage 1: default initialization.
100       ptr->SetDirectory(0);
101       //"Default" (for derived component) histograms.
102       FillHistogramDefinitions();
103       //Default values.
104       fMaxEntries = kMaxEntries;
105       fPublishInterval = kInterval;
106       fLastTime = 0;
107       //Stage 2: OCDB.
108       TString cdbPath("HLT/ConfigHLT/");
109       cdbPath += GetComponentID();
110       //
111       iResult = ConfigureFromCDBTObjString(cdbPath);
112       //
113       if (iResult < 0)
114         return iResult;
115       //Stage 3: command line arguments.
116       if (argc && (iResult = ConfigureFromArgumentString(argc, argv)) < 0)
117         return iResult;
118
119       ptr->SetCircular(fMaxEntries);
120       fTree = ptr.release();
121     } else //No way to process error correctly - error is unknown here.
122       return -EINVAL;
123   } else {
124     HLTError("fTree pointer must be null before DoInit call");
125     return -EINVAL;
126   }
127
128   if (iResult>=0 && fMaxEventTime>0) {
129     fpEventTimer=new TStopwatch;
130     if (fpEventTimer) {
131       fpEventTimer->Reset();
132     }
133     fpCycleTimer=new TStopwatch;
134     if (fpCycleTimer) {
135       fpCycleTimer->Reset();
136     }
137   }
138   fSkippedEventsCount=0;
139
140   return iResult;
141 }
142
143 int AliHLTTTreeProcessor::DoDeinit()
144 {
145   // cleanup component
146   delete fTree;
147   fTree = 0;
148   fDefinitions.clear();
149
150   if (fpEventTimer) delete fpEventTimer;
151   fpEventTimer=NULL;
152   if (fpCycleTimer) delete fpCycleTimer;
153   fpCycleTimer=NULL;
154
155   return 0;
156 }
157
158 int AliHLTTTreeProcessor::DoEvent(const AliHLTComponentEventData& evtData, AliHLTComponentTriggerData& trigData)
159 {
160   //Process event and publish histograms.
161   AliHLTUInt32_t eventType=0;
162   if (!IsDataEvent(&eventType) && eventType!=gkAliEventTypeEndOfRun) return 0;
163
164   //I'm pretty sure, that if fTree == 0 (DoInit failed) DoEvent is not called.
165   //But interface itself does not force you to call DoInit before DoEvent, so,
166   //I make this check explicit.
167   if (!fTree) {
168     HLTError("fTree is a null pointer, try to call AliHLTTTreeProcessor::DoInit first.");
169     return -EINVAL;//-ENULLTREE? :)
170   }
171
172   AliHLTUInt32_t averageEventTime=0;
173   AliHLTUInt32_t averageCycleTime=0;
174
175   int fillingtime=0;
176   int publishtime=0;
177   bool bDoFilling=true;
178   bool bDoPublishing=false;
179   const int cycleResetInterval=1000;
180   if (fpEventTimer && fpCycleTimer) {
181     averageEventTime=(fpEventTimer->RealTime()*fgkTimeScale)/(GetEventCount()+1);
182     fillingtime=fpEventTimer->RealTime()*fgkTimeScale;
183     publishtime=fillingtime;
184     fpEventTimer->Start(kFALSE);
185     fpCycleTimer->Stop();
186     averageCycleTime=(fpCycleTimer->RealTime()*fgkTimeScale)/((GetEventCount()%cycleResetInterval)+1);
187     // adapt processing to 3/4 of the max time
188     bDoFilling=4*averageEventTime<3*fMaxEventTime || averageEventTime<averageCycleTime;
189     if (fNofEventsForce>0 && fForcedEventsCount<fNofEventsForce) {
190       fForcedEventsCount++;
191       bDoFilling=true;
192     }
193   }
194
195   // process input data blocks and fill the tree
196   int iResult = 0;
197   if (eventType!=gkAliEventTypeEndOfRun) {
198     if (bDoFilling) {iResult=FillTree(fTree, evtData, trigData); fNewEventsCount++;}
199     else fSkippedEventsCount++;
200   }
201   if (fpEventTimer) {
202     fpEventTimer->Stop();
203     fillingtime=fpEventTimer->RealTime()*fgkTimeScale-fillingtime;
204     fpEventTimer->Start(kFALSE);
205   }
206
207   if (iResult < 0)
208     return iResult;
209
210   const TDatime time;
211
212   if (( time.Get() - fLastTime > fPublishInterval && fNewEventsCount>0) ||
213       eventType==gkAliEventTypeEndOfRun) {
214     bDoPublishing=true;
215     for (list_const_iterator i = fDefinitions.begin(); i != fDefinitions.end(); ++i) {
216       if (TH1* h = CreateHistogram(*i)) {
217         //I do not care about errors here - since I'm not able
218         //to rollback changes.
219         // TODO: in case of -ENOSPC et the size of the last object by calling
220         // GetLastObjectSize() and accumulate the necessary output buffer size
221         PushBack(h, GetOriginDataType(), GetDataSpec());
222       }
223     }
224       unsigned eventcount=GetEventCount()+1;
225       HLTBenchmark("publishing %d histograms, %d entries in tree, %d new events since last publishing, accumulated %d of %d events (%.1f%%)", fDefinitions.size(), fTree->GetEntriesFast(), fNewEventsCount, eventcount-fSkippedEventsCount, eventcount, eventcount>0?(100*float(eventcount-fSkippedEventsCount)/eventcount):0);
226     fNewEventsCount=0;
227
228     fLastTime = time.Get();
229   }
230
231   if (fpEventTimer) {
232     fpEventTimer->Stop();
233     publishtime=fpEventTimer->RealTime()*fgkTimeScale-publishtime;
234     if (publishtime>fillingtime) publishtime-=fillingtime;
235     else publishtime=0;
236
237     averageEventTime=(fpEventTimer->RealTime()*fgkTimeScale)/(GetEventCount()+1);
238
239     // info output once every 5 seconds
240     static UInt_t lastTime=0;
241     if (time.Get()-lastTime>5 ||
242         eventType==gkAliEventTypeEndOfRun ||
243         bDoPublishing) {
244       lastTime=time.Get();
245       unsigned eventcount=GetEventCount()+1;
246       HLTBenchmark("filling time %d us, publishing time %d, average total processing time %d us, cycle time %d us, accumulated %d of %d events (%.1f%%)", fillingtime, publishtime, averageEventTime, averageCycleTime, eventcount-fSkippedEventsCount, eventcount, eventcount>0?(100*float(eventcount-fSkippedEventsCount)/eventcount):0);
247     }
248   }
249   if (fpCycleTimer) {
250     bool bReset=(GetEventCount()%cycleResetInterval)==0;
251     fpCycleTimer->Start(bReset);
252   }
253
254   return iResult;
255 }
256
257 int AliHLTTTreeProcessor::ScanConfigurationArgument(int argc, const char** argv)
258 {
259   // scan one argument and its parameters from the list
260   // return number of processed entries.
261   // possible arguments: 
262   // -maxentries number
263   // -interval number
264   // -histogram name -size number -expression expression [-cut expression ][-opt option]
265   // As soon as "-histogram" found, -size and -expression and -outtype are required, 
266   // cut and option can be omitted.
267   if (argc <= 0)
268     return 0;
269
270   std::list<AliHLTHistogramDefinition> newDefs;
271   AliHLTHistogramDefinition def;
272
273   int i = 0;
274   int maxEntries = 0;
275
276   while (i < argc) {
277     const TString argument(argv[i]);
278
279     if (argument.CompareTo("-maxentries") == 0) { //1. Max entries argument for TTree.
280       if (i + 1 == argc) {
281         HLTError("Numeric value for '-maxentries' is expected");
282         return -EPROTO;
283       }
284       //Next must be a number.
285       //TString returns 0 (number) even if string contains non-numeric symbols.
286       maxEntries = TString(argv[i + 1]).Atoi();
287       if (maxEntries <= 0) {
288         HLTError("Bad value for '-maxentries': %d", maxEntries);
289         return -EPROTO;
290       }
291   
292       i += 2;
293     } else if (argument.CompareTo("-interval") == 0) { //2. Interval argument for publishing.
294       if (i + 1 == argc) {
295         HLTError("Numeric value for '-interval' is expected");
296         return -EPROTO;
297       }
298
299       const Int_t interval = TString(argv[i + 1]).Atoi();
300       if (interval <= 0) {
301         HLTError("Bad value for '-interval' argument: %d", interval);
302         return -EPROTO;
303       }
304
305       fPublishInterval = interval;
306
307       i += 2;
308     } else if (argument.CompareTo("-maxeventtime") == 0) { // max average processing time in us
309       if (i + 1 == argc) {
310         HLTError("Numeric value for '-maxeventtime' is expected");
311         return -EPROTO;
312       }
313
314       const Int_t time = TString(argv[i + 1]).Atoi();
315       if (time <= 0) {
316         HLTError("Bad value for '-maxeventtime' argument: %d", time);
317         return -EPROTO;
318       }
319
320       fMaxEventTime = time;
321
322       i += 2;
323     } else if (argument.CompareTo("-forced-events") == 0) { // number of forced events
324       if (i + 1 == argc) {
325         HLTError("Numeric value for '-forced-events' is expected");
326         return -EPROTO;
327       }
328
329       const Int_t count = TString(argv[i + 1]).Atoi();
330       if (count <= 0) {
331         HLTError("Bad value for '-forced-events' argument: %d", count);
332         return -EPROTO;
333       }
334
335       fNofEventsForce = count;
336       fForcedEventsCount=0;
337
338       i += 2;
339     } else if (argument.CompareTo("-histogram") == 0) { //3. Histogramm definition.
340       const int nParsed = ParseHistogramDefinition(argc, argv, i, def);
341       if (!nParsed)
342         return -EPROTO;
343
344       newDefs.push_back(def);
345
346       i += nParsed;   
347     } else {
348       HLTError("Unknown argument %s", argument.Data());
349       return -EPROTO;
350     }
351   }
352
353   if (maxEntries != fMaxEntries) {
354     fMaxEntries = maxEntries;
355     if (fTree) {
356       fTree->Reset();
357       fTree->SetCircular(fMaxEntries);
358     }
359   }
360
361   if (newDefs.size())
362     fDefinitions.swap(newDefs);
363
364   return i;
365 }
366
367 TH1* AliHLTTTreeProcessor::CreateHistogram(const AliHLTHistogramDefinition& d)
368 {
369
370   // create a histogram from the tree
371   if (!fTree) {
372     HLTError("fTree is a null pointer, try to call AliHLTTTreeProcessor::DoInit first.");
373     return 0;
374   }
375
376   TString histName(d.GetName());
377   if (!histName.Contains("(")) {
378     //Without number of bins, the histogram will be "fixed"
379     //and most of values can go to underflow/overflow bins,
380     //since kCanRebin will be false.
381     histName += TString::Format("(%d)", Int_t(kDefaultNBins));
382   }
383
384   const Long64_t rez = fTree->Project(histName.Data(), d.GetExpression().Data(), d.GetCut().Data(), d.GetDrawOption().Data());
385
386   if (rez == -1) {
387     HLTError("TTree::Project failed");
388     return 0;
389   }
390
391   //Now, cut off the binning part of a name
392   histName = histName(0, histName.Index("("));
393   TH1 * hist = dynamic_cast<TH1*>(gDirectory->Get(histName.Data()));
394   if (!hist) {
395     const TString msg(Form("Hist %s is a null pointer, selection was %s, strange name or hist's type\n", histName.Data(), d.GetExpression().Data()));
396     HLTError(msg.Data());
397   }else if (d.GetDrawOption().Length()) {
398     hist->SetOption(d.GetDrawOption().Data());
399   }
400
401   return hist;
402 }
403
404 int AliHLTTTreeProcessor::ParseHistogramDefinition(int argc, const char** argv, int pos, AliHLTHistogramDefinition& dst)const
405 {
406   //Histogram-definition:
407   //    -histogram name -size number -expression expression [-cut expression][-opt option]
408
409   //at pos we have '-histogram', at pos + 1 must be the name.
410   if (pos + 1 == argc) {
411     HLTError("Bad histogram definition, histogram name is expected");
412     return 0;
413   }
414
415   dst.SetName(argv[pos + 1]);
416   pos += 2;
417   
418   //At pos must be '-size', and number at pos + 1.
419   if (pos == argc || TString(argv[pos]).CompareTo("-size")) {
420     HLTError("Bad histogram definition, '-size' is expected");
421     return 0;
422   }
423
424   if (pos + 1 == argc) {
425     HLTError("Bad histogram definition, size is expected");
426     return 0;
427   }
428
429   dst.SetSize(TString(argv[pos + 1]).Atoi());
430   if (dst.GetSize() <= 0) {
431     HLTError("Bad histogram definition, positive size is required");
432     return 0;
433   }
434
435   pos += 2;
436   //At pos must be '-expression', and expression at pos + 1. 
437   if (pos == argc || TString(argv[pos]).CompareTo("-expression")) {
438     HLTError("Bad histogram definition, '-expression' is expected");
439     return 0;
440   }
441
442   if (pos + 1 == argc) {
443     HLTError("Bad histogram definition, expression is expected");
444     return 0;
445   }
446
447   dst.SetExpression(argv[pos + 1]);
448   pos += 2;
449
450   int processed = 6;
451   dst.SetCut("");
452   dst.SetDrawOption("");
453
454   //remaining options can be the cut and Draw option.
455   //cut must be first.
456   if (pos + 1 >= argc)
457     return processed;
458
459   if (TString(argv[pos]).CompareTo("-cut") == 0) {
460     dst.SetCut(argv[pos + 1]);
461     pos += 2;
462     processed += 2;
463   }
464
465   if (pos + 1 >= argc)
466     return processed;
467
468   if (TString(argv[pos]).CompareTo("-opt") == 0) {
469     dst.SetDrawOption(argv[pos + 1]);
470     processed += 2;
471   }
472
473   return processed;
474 }