]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - HLT/BASE/AliHLTTTreeProcessor.cxx
configurable histogram titles added (Hege)
[u/mrichter/AliRoot.git] / HLT / BASE / AliHLTTTreeProcessor.cxx
... / ...
CommitLineData
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 "AliHLTErrorGuard.h"
30#include "TDirectory.h"
31#include "TDatime.h"
32#include "TString.h"
33#include "TTree.h"
34#include "TH1.h"
35#include "TStopwatch.h"
36#include "TUUID.h"
37#include "TSystem.h"
38#include "TRandom3.h"
39
40/** ROOT macro for the implementation of ROOT specific class methods */
41ClassImp(AliHLTTTreeProcessor)
42
43AliHLTTTreeProcessor::AliHLTTTreeProcessor()
44 : AliHLTProcessor(),
45 fDefinitions(),
46 fTree(0),
47 fMaxEntries(kMaxEntries),
48 fPublishInterval(kInterval),
49 fLastTime(0),
50 fpEventTimer(NULL),
51 fpCycleTimer(NULL),
52 fMaxMemory(700000),
53 fMaxEventTime(0),
54 fNofEventsForce(0),
55 fForcedEventsCount(0),
56 fSkippedEventsCount(0),
57 fNewEventsCount(0),
58 fUniqueId(0),
59 fIgnoreCycleTime(10),
60 fCycleTimeFactor(1.0)
61{
62 // see header file for class documentation
63 // or
64 // refer to README to build package
65 // or
66 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
67}
68
69const AliHLTUInt32_t AliHLTTTreeProcessor::fgkTimeScale=1000000; // ticks per second
70
71AliHLTTTreeProcessor::~AliHLTTTreeProcessor()
72{
73 // see header file for class documentation
74}
75
76AliHLTComponentDataType AliHLTTTreeProcessor::GetOutputDataType()
77{
78 // get the component output data type
79 return kAliHLTDataTypeHistogram;
80}
81
82void AliHLTTTreeProcessor::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
83{
84 // get the output size estimator
85 //
86 if (!fDefinitions.size()) {
87 HLTError("Can not calculate output data size, no histogram definitions were provided");
88 return;
89 }
90
91 constBase = 0;
92 for (list_const_iterator i = fDefinitions.begin(); i != fDefinitions.end(); ++i)
93 constBase += i->GetSize();
94
95 inputMultiplier = 1.;
96}
97
98int AliHLTTTreeProcessor::DoInit(int argc, const char** argv)
99{
100 // init component
101 // ask child to create the tree.
102 int iResult = 0;
103
104 // component configuration
105 //Stage 1: default initialization.
106 //"Default" (for derived component) histograms.
107 FillHistogramDefinitions();
108 //Default values.
109 fMaxEntries = kMaxEntries;
110 fPublishInterval = kInterval;
111 fLastTime = 0;
112 //Stage 2: OCDB.
113 TString cdbPath("HLT/ConfigHLT/");
114 cdbPath += GetComponentID();
115 //
116 iResult = ConfigureFromCDBTObjString(cdbPath);
117 //
118 if (iResult < 0)
119 return iResult;
120 //Stage 3: command line arguments.
121 if (argc && (iResult = ConfigureFromArgumentString(argc, argv)) < 0)
122 return iResult;
123
124 // calculating a unique id from the hostname and process id
125 // used for identifying output of multiple components
126 TUUID guid = GenerateGUID();
127 union
128 {
129 UChar_t buf[16];
130 UInt_t bufAsInt[4];
131 };
132 guid.GetUUID(buf);
133 fUniqueId = bufAsInt[0];
134
135 if (!fTree) {
136 // originally foreseen to pass the arguments to the function, however
137 // this is not appropriate. Argument scan via overloaded function
138 // ScanConfigurationArgument
139 std::auto_ptr<TTree> ptr(CreateTree(0, NULL));
140 if (ptr.get()) {
141 ptr->SetDirectory(0);
142 ptr->SetCircular(fMaxEntries);
143 fTree = ptr.release();
144 } else //No way to process error correctly - error is unknown here.
145 return -EINVAL;
146 } else {
147 HLTError("fTree pointer must be null before DoInit call");
148 return -EINVAL;
149 }
150
151 if (iResult>=0 && fMaxEventTime>0) {
152 fpEventTimer=new TStopwatch;
153 if (fpEventTimer) {
154 fpEventTimer->Reset();
155 }
156 fpCycleTimer=new TStopwatch;
157 if (fpCycleTimer) {
158 fpCycleTimer->Reset();
159 }
160 }
161 fSkippedEventsCount=0;
162
163 return iResult;
164}
165
166int AliHLTTTreeProcessor::DoDeinit()
167{
168 // cleanup component
169 delete fTree;
170 fTree = 0;
171 fDefinitions.clear();
172
173 if (fpEventTimer) delete fpEventTimer;
174 fpEventTimer=NULL;
175 if (fpCycleTimer) delete fpCycleTimer;
176 fpCycleTimer=NULL;
177
178 return 0;
179}
180
181int AliHLTTTreeProcessor::DoEvent(const AliHLTComponentEventData& evtData, AliHLTComponentTriggerData& trigData)
182{
183 //Process event and publish histograms.
184 AliHLTUInt32_t eventType=0;
185 if (!IsDataEvent(&eventType) && eventType!=gkAliEventTypeEndOfRun) return 0;
186
187 //I'm pretty sure, that if fTree == 0 (DoInit failed) DoEvent is not called.
188 //But interface itself does not force you to call DoInit before DoEvent, so,
189 //I make this check explicit.
190 if (!fTree) {
191 HLTError("fTree is a null pointer, try to call AliHLTTTreeProcessor::DoInit first.");
192 return -EINVAL;//-ENULLTREE? :)
193 }
194
195 AliHLTUInt32_t averageEventTime=0;
196 AliHLTUInt32_t averageCycleTime=0;
197
198 int fillingtime=0;
199 int publishtime=0;
200 bool bDoFilling=true;
201 bool bDoPublishing=false;
202 const int cycleResetInterval=1000;
203 if (fpEventTimer && fpCycleTimer) {
204 averageEventTime=AliHLTUInt32_t(fpEventTimer->RealTime()*fgkTimeScale)/(GetEventCount()+1);
205 fillingtime=int(fpEventTimer->RealTime()*fgkTimeScale);
206 publishtime=fillingtime;
207 fpEventTimer->Start(kFALSE);
208 fpCycleTimer->Stop();
209 averageCycleTime=AliHLTUInt32_t(fpCycleTimer->RealTime()*fgkTimeScale)/((GetEventCount()%cycleResetInterval)+1);
210 // adapt processing to 3/4 of the max time
211 bDoFilling=4*averageEventTime<3*fMaxEventTime ||
212 (averageEventTime<fCycleTimeFactor*averageCycleTime && fpCycleTimer->RealTime()>fIgnoreCycleTime);
213 if (fNofEventsForce>0 && fForcedEventsCount<fNofEventsForce) {
214 fForcedEventsCount++;
215 bDoFilling=true;
216 }
217 }
218
219 // FIXME: there is still an unclear increase in memory consumption, even if the number of entries
220 // in the tree is restricted. Valgrind studies did not show an obvious memory leak. This is likely
221 // to be caused by something deep in the Root TTree functionality and needs to be studied in detail.
222 ProcInfo_t ProcInfo;
223 gSystem->GetProcInfo(&ProcInfo);
224 if (ProcInfo.fMemResident>fMaxMemory) bDoFilling=false;
225
226 // process input data blocks and fill the tree
227 int iResult = 0;
228 if (eventType!=gkAliEventTypeEndOfRun) {
229 if (bDoFilling) {iResult=FillTree(fTree, evtData, trigData); fNewEventsCount++;}
230 else fSkippedEventsCount++;
231 }
232 if (fpEventTimer) {
233 fpEventTimer->Stop();
234 fillingtime=int(fpEventTimer->RealTime()*fgkTimeScale)-fillingtime;
235 if (fillingtime<0) fillingtime=0;
236 fpEventTimer->Start(kFALSE);
237 }
238
239 if (iResult < 0) {
240 ALIHLTERRORGUARD(5, "FillTree failed with %d, first event %d", iResult, GetEventCount());
241 return iResult;
242 }
243
244 const TDatime time;
245
246 if (( time.Get() - fLastTime > fPublishInterval && fNewEventsCount>0) ||
247 eventType==gkAliEventTypeEndOfRun) {
248 if ((bDoPublishing=fLastTime>0)) { // publish earliest after the first interval but set the timer
249
250 for (list_const_iterator i = fDefinitions.begin(); i != fDefinitions.end(); ++i) {
251 if (TH1* h = CreateHistogram(*i)) {
252 //I do not care about errors here - since I'm not able
253 //to rollback changes.
254 // TODO: in case of -ENOSPC et the size of the last object by calling
255 // GetLastObjectSize() and accumulate the necessary output buffer size
256 PushBack(h, GetOriginDataType(), GetDataSpec());
257 delete h;
258 }
259 }
260 unsigned eventcount=GetEventCount()+1;
261 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);
262 fNewEventsCount=0;
263 HLTBenchmark("current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual);
264 }
265
266 fLastTime=time.Get();
267 if (fLastTime==0) {
268 // choose a random offset at beginning to equalize traffic for multiple instances
269 // of the component
270 gRandom->SetSeed(fUniqueId);
271 fLastTime-=gRandom->Integer(fPublishInterval);
272 }
273 }
274
275 if (fpEventTimer) {
276 fpEventTimer->Stop();
277 publishtime=int(fpEventTimer->RealTime()*fgkTimeScale)-publishtime;
278 if (publishtime>fillingtime) publishtime-=fillingtime;
279 else publishtime=0;
280
281 averageEventTime=AliHLTUInt32_t(fpEventTimer->RealTime()*fgkTimeScale)/(GetEventCount()+1);
282
283 // info output once every 5 seconds
284 static UInt_t lastTime=0;
285 if (time.Get()-lastTime>5 ||
286 eventType==gkAliEventTypeEndOfRun ||
287 bDoPublishing) {
288 lastTime=time.Get();
289 unsigned eventcount=GetEventCount()+1;
290 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);
291 }
292 }
293 if (fpCycleTimer) {
294 bool bReset=(GetEventCount()%cycleResetInterval)==0;
295 fpCycleTimer->Start(bReset);
296 }
297
298 return iResult;
299}
300
301int AliHLTTTreeProcessor::ScanConfigurationArgument(int argc, const char** argv)
302{
303 // scan one argument and its parameters from the list
304 // return number of processed entries.
305 // possible arguments:
306 // -maxentries number
307 // -interval number
308 // -histogram name -size number -expression expression [-title expression ] -cut expression ][-opt option]
309 // As soon as "-histogram" found, -size and -expression and -outtype are required,
310 // cut and option can be omitted.
311 if (argc <= 0)
312 return 0;
313
314 std::list<AliHLTHistogramDefinition> newDefs;
315 AliHLTHistogramDefinition def;
316
317 int i = 0;
318 int maxEntries = fMaxEntries;
319
320 while (i < argc) {
321 const TString argument(argv[i]);
322
323 if (argument.CompareTo("-maxentries") == 0) { //1. Max entries argument for TTree.
324 if (i + 1 == argc) {
325 HLTError("Numeric value for '-maxentries' is expected");
326 return -EPROTO;
327 }
328 //Next must be a number.
329 //TString returns 0 (number) even if string contains non-numeric symbols.
330 maxEntries = TString(argv[i + 1]).Atoi();
331 if (maxEntries <= 0) {
332 HLTError("Bad value for '-maxentries': %d", maxEntries);
333 return -EPROTO;
334 }
335
336 i += 2;
337 } else if (argument.CompareTo("-interval") == 0) { //2. Interval argument for publishing.
338 if (i + 1 == argc) {
339 HLTError("Numeric value for '-interval' is expected");
340 return -EPROTO;
341 }
342
343 const Int_t interval = TString(argv[i + 1]).Atoi();
344 if (interval < 0) {
345 HLTError("Bad value for '-interval' argument: %d", interval);
346 return -EPROTO;
347 }
348
349 fPublishInterval = interval;
350
351 i += 2;
352 } else if (argument.CompareTo("-maxeventtime") == 0) { // max average processing time in us
353 if (i + 1 == argc) {
354 HLTError("Numeric value for '-maxeventtime' is expected");
355 return -EPROTO;
356 }
357
358 const Int_t time = TString(argv[i + 1]).Atoi();
359 if (time < 0) {
360 HLTError("Bad value for '-maxeventtime' argument: %d", time);
361 return -EPROTO;
362 }
363
364 fMaxEventTime = time;
365
366 i += 2;
367 } else if (argument.CompareTo("-forced-events") == 0) { // number of forced events
368 if (i + 1 == argc) {
369 HLTError("Numeric value for '-forced-events' is expected");
370 return -EPROTO;
371 }
372
373 const Int_t count = TString(argv[i + 1]).Atoi();
374 if (count < 0) {
375 HLTError("Bad value for '-forced-events' argument: %d", count);
376 return -EPROTO;
377 }
378
379 fNofEventsForce = count;
380 fForcedEventsCount=0;
381
382 i += 2;
383 } else if (argument.CompareTo("-ignore-cycletime") == 0) { // ignore cycle time for n sec
384 if (i + 1 == argc) {
385 HLTError("Numeric value for '-ignore-cycletime' is expected");
386 return -EPROTO;
387 }
388
389 const Int_t time = TString(argv[i + 1]).Atoi();
390 if (time < 0) {
391 HLTError("Bad value for '-ignore-cycletime' argument: %d", time);
392 return -EPROTO;
393 }
394
395 fIgnoreCycleTime = time;
396 i += 2;
397 } else if (argument.CompareTo("-maxmemory") == 0) { // maximum of memory in kByte to be used by the component
398 if (i + 1 == argc) {
399 HLTError("Numeric value for '-maxmemory' is expected");
400 return -EPROTO;
401 }
402
403 const Int_t mem = TString(argv[i + 1]).Atoi();
404 if (mem < 0) {
405 HLTError("Bad value for '-maxmemory' argument: %d", time);
406 return -EPROTO;
407 }
408
409 fMaxMemory = mem;
410 i += 2;
411 } else if (argument.CompareTo("-cycletime-factor") == 0) { // weight factor for cycle time
412 if (i + 1 == argc) {
413 HLTError("Numeric value for '-cycletime-factor' is expected");
414 return -EPROTO;
415 }
416
417 const Float_t factor = TString(argv[i + 1]).Atof();
418 if (factor < 0) {
419 HLTError("Bad value for '-cycletime-factor' argument: %f", factor);
420 return -EPROTO;
421 }
422
423 fCycleTimeFactor = factor;
424 i += 2;
425 } else if (argument.CompareTo("-histogram") == 0) { //3. Histogramm definition.
426 const int nParsed = ParseHistogramDefinition(argc, argv, i, def);
427 if (!nParsed)
428 return -EPROTO;
429
430 newDefs.push_back(def);
431
432 i += nParsed;
433 } else {
434 HLTError("Unknown argument %s", argument.Data());
435 return -EPROTO;
436 }
437 }
438
439 if (maxEntries != fMaxEntries) {
440 fMaxEntries = maxEntries;
441 if (fTree) {
442 fTree->Reset();
443 fTree->SetCircular(fMaxEntries);
444 }
445 }
446
447 if (newDefs.size())
448 fDefinitions.swap(newDefs);
449
450 return i;
451}
452
453TH1* AliHLTTTreeProcessor::CreateHistogram(const AliHLTHistogramDefinition& d)
454{
455
456 // create a histogram from the tree
457 if (!fTree) {
458 HLTError("fTree is a null pointer, try to call AliHLTTTreeProcessor::DoInit first.");
459 return 0;
460 }
461
462 TString histName(d.GetName());
463 if (!histName.Contains("(")) {
464 //Without number of bins, the histogram will be "fixed"
465 //and most of values can go to underflow/overflow bins,
466 //since kCanRebin will be false.
467 histName += TString::Format("(%d)", Int_t(kDefaultNBins));
468 }
469
470 const Long64_t rez = fTree->Project(histName.Data(), d.GetExpression().Data(), d.GetCut().Data(), d.GetDrawOption().Data());
471
472 if (rez == -1) {
473 HLTError("TTree::Project failed");
474 return 0;
475 }
476
477 //Now, cut off the binning part of a name
478 histName = histName(0, histName.Index("("));
479 TH1 * hist = dynamic_cast<TH1*>(gDirectory->Get(histName.Data()));
480 if (!hist) {
481 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()));
482 HLTError(msg.Data());
483 }else if (d.GetDrawOption().Length()) {
484 hist->SetOption(d.GetDrawOption().Data());
485 }
486
487 //Reformatting the histogram name
488 TString str2=d.GetCut().Data();
489 str2.ReplaceAll("Track_", "");
490 str2.ReplaceAll("&&", " ");
491 str2 = histName+" "+str2;
492 hist->SetTitle(str2);
493
494 if(d.GetTitle().Length()){
495
496 //removing underscore
497 size_t found;
498 string str=d.GetTitle().Data();
499 found=str.find_first_of("_");
500 if(!(d.GetExpression().CompareTo("Track_pt"))){
501 found=str.find_first_of("_",found+1);
502 }
503 str[found]=' ';
504 char axis[100];
505 sprintf(axis,"%s",str.c_str());
506
507 hist->SetXTitle(axis);
508 hist->GetXaxis()->CenterTitle();
509 }
510 return hist;
511}
512
513int AliHLTTTreeProcessor::ParseHistogramDefinition(int argc, const char** argv, int pos, AliHLTHistogramDefinition& dst)const
514{
515 //Histogram-definition:
516 // -histogram name -size number -expression expression [-title expression][-cut expression][-opt option]
517
518 //at pos we have '-histogram', at pos + 1 must be the name.
519 if (pos + 1 == argc) {
520 HLTError("Bad histogram definition, histogram name is expected");
521 return 0;
522 }
523
524 dst.SetName(argv[pos + 1]);
525 pos += 2;
526
527 //At pos must be '-size', and number at pos + 1.
528 if (pos == argc || TString(argv[pos]).CompareTo("-size")) {
529 HLTError("Bad histogram definition, '-size' is expected");
530 return 0;
531 }
532
533 if (pos + 1 == argc) {
534 HLTError("Bad histogram definition, size is expected");
535 return 0;
536 }
537
538 dst.SetSize(TString(argv[pos + 1]).Atoi());
539 if (dst.GetSize() <= 0) {
540 HLTError("Bad histogram definition, positive size is required");
541 return 0;
542 }
543
544 pos += 2;
545 //At pos must be '-expression', and expression at pos + 1.
546 if (pos == argc || TString(argv[pos]).CompareTo("-expression")) {
547 HLTError("Bad histogram definition, '-expression' is expected");
548 return 0;
549 }
550
551 if (pos + 1 == argc) {
552 HLTError("Bad histogram definition, expression is expected");
553 return 0;
554 }
555
556 dst.SetExpression(argv[pos + 1]);
557 pos += 2;
558
559 int processed = 6;
560 dst.SetTitle("");
561 dst.SetCut("");
562 dst.SetDrawOption("");
563
564 //remaining options can be the title, cut and Draw option.
565 //title must be first
566 if (pos + 1 >= argc){
567 return processed;
568 }
569 if (TString(argv[pos]).CompareTo("-title") == 0) {
570 dst.SetTitle(argv[pos + 1]);
571 pos += 2;
572 processed += 2;
573 }
574
575 //cut must be second.
576 if (pos + 1 >= argc)
577 return processed;
578
579 if (TString(argv[pos]).CompareTo("-cut") == 0) {
580 dst.SetCut(argv[pos + 1]);
581 pos += 2;
582 processed += 2;
583 }
584
585 if (pos + 1 >= argc)
586 return processed;
587
588 if (TString(argv[pos]).CompareTo("-opt") == 0) {
589 dst.SetDrawOption(argv[pos + 1]);
590 processed += 2;
591 }
592
593 return processed;
594}