]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/BASE/util/AliHLTCompStatCollector.cxx
fixing warnings
[u/mrichter/AliRoot.git] / HLT / BASE / util / AliHLTCompStatCollector.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: Matthias Richter <Matthias.Richter@ift.uib.no>        *
8 //*                  for The ALICE HLT Project.                            *
9 //*                                                                        *
10 //* Permission to use, copy, modify and distribute this software and its   *
11 //* documentation strictly for non-commercial purposes is hereby granted   *
12 //* without fee, provided that the above copyright notice appears in all   *
13 //* copies and that both the copyright notice and this permission notice   *
14 //* appear in the supporting documentation. The authors make no claims     *
15 //* about the suitability of this software for any purpose. It is          *
16 //* provided "as is" without express or implied warranty.                  *
17 //**************************************************************************
18
19 /** @file   AliHLTCompStatCollector.cxx
20     @author Matthias Richter
21     @date   
22     @brief  Collector component for the component statistics information.
23 */
24
25 #include "AliHLTCompStatCollector.h"
26 #include "TFile.h"
27 #include "TStopwatch.h"
28 #include "TH1F.h"
29 #include "TH2F.h"
30 #include "TH2C.h"
31 #include "TTree.h"
32 #include "TFolder.h"
33 #include "TNamed.h"
34 #include "TString.h"
35 #include <cassert>
36
37 #define HLTSTAT_FOLDER_NAME               "HLTstat"
38 #define HLTSTAT_FOLDER_DESC               "ALICE HLT component statistics"
39 #define HLTSTAT_ENTRY_PARENT_FOLDER_NAME  "parents"
40 #define HLTSTAT_ENTRY_PARENT_FOLDER_DESC  "parent components"
41 #define HLTSTAT_ENTRY_PROPS_FOLDER_NAME   "props"
42 #define HLTSTAT_ENTRY_PROPS_FOLDER_DESC   "component properties"
43 #define HLTSTAT_ENTRY_PROPS_IDOBJ_NAME    "id"
44 #define HLTSTAT_ENTRY_PROPS_IDOBJ_DESC    "numerical id calculated from chain id"
45
46 /** ROOT macro for the implementation of ROOT specific class methods */
47 ClassImp(AliHLTCompStatCollector)
48
49 AliHLTCompStatCollector::AliHLTCompStatCollector()
50   :
51   AliHLTProcessor(),
52   fpTimer(NULL),
53   fpFolder(NULL),
54   fpStatTree(NULL),
55   fCycleTime(0),
56   fNofSets(0),
57   fArraySize(1000),
58   fPosition(0),
59   fpLevelArray(NULL),
60   fpSpecArray(NULL),
61   fpBlockNoArray(NULL),
62   fpIdArray(NULL),
63   fpTimeArray(NULL),
64   fpCTimeArray(NULL),
65   fpInputBlockCountArray(NULL),
66   fpTotalInputSizeArray(NULL),
67   fpNormalizedInputSizeArray(NULL),
68   fpOutputBlockCountArray(NULL),
69   fpTotalOutputSizeArray(NULL)
70   , fpInputOutputRatioArray(NULL)
71   , fpNormalizedInputOutputRatioArray(NULL)
72   , fpComponentCycleTimeArray(NULL)
73   , fpEventTypeArray(NULL)
74   , fpEventCountArray(NULL)
75   , fSizeEstimator(1000)
76   , fMode(kPublishObjects)
77   , fFileName()
78   , fFile(NULL)
79   , fLastTime(time(NULL))
80   , fPeriod(0)
81   , fEventModulo(0)
82 {
83   // see header file for class documentation
84   // or
85   // refer to README to build package
86   // or
87   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
88 }
89
90 AliHLTCompStatCollector::~AliHLTCompStatCollector()
91 {
92   // see header file for class documentation
93   ClearAll();
94 }
95
96 void AliHLTCompStatCollector::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
97 {
98   // see header file for class documentation
99   list.push_back(kAliHLTDataTypeComponentStatistics);
100 }
101
102 AliHLTComponentDataType AliHLTCompStatCollector::GetOutputDataType()
103 {
104   // see header file for class documentation
105   return kAliHLTMultipleDataType;
106 }
107
108 int AliHLTCompStatCollector::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
109 {
110   // see header file for class documentation
111   tgtList.clear();
112   tgtList.push_back(kAliHLTDataTypeHistogram);
113   tgtList.push_back(kAliHLTDataTypeTTree);
114   return tgtList.size();
115 }
116
117 void AliHLTCompStatCollector::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
118 {
119   // see header file for class documentation
120   constBase=fSizeEstimator;
121   inputMultiplier=100.0;
122 }
123
124 int AliHLTCompStatCollector::DoInit( int argc, const char** argv )
125 {
126   // see header file for class documentation
127   int iResult=0;
128   TString argument="";
129   int bMissingParam=0;
130   for (int i=0; i<argc && iResult>=0; i++) {
131     argument=argv[i];
132     if (argument.IsNull()) continue;
133
134     // -file
135     if (argument.CompareTo("-file")==0) {
136       if ((bMissingParam=(++i>=argc))) break;
137       fFileName=argv[i];
138       fMode|=kSaveObjects;
139
140     // -modulo
141     } else if (argument.CompareTo("-modulo")==0) {
142       if ((bMissingParam=(++i>=argc))) break;
143       TString param=argv[i];
144       if (param.IsDigit()) {
145         fEventModulo=param.Atoi();
146       } else {
147         HLTError("expecting number as parameter for option %s", argument.Data());
148         iResult=-EINVAL;
149       }
150
151     // -period
152     } else if (argument.CompareTo("-period")==0) {
153       if ((bMissingParam=(++i>=argc))) break;
154       TString param=argv[i];
155       if (param.IsDigit()) {
156         fPeriod=param.Atoi();
157       } else {
158         HLTError("expecting number as parameter for option %s", argument.Data());
159         iResult=-EINVAL;
160       }
161
162     // -publish
163     } else if (argument.CompareTo("-publish")==0) {
164       if ((bMissingParam=(++i>=argc))) break;
165       TString param=argv[i];
166       if (param.IsDigit()) {
167         if (param.Atoi()==1) fMode|=kPublishObjects;
168         else if (param.Atoi()==0) fMode&=~kPublishObjects;
169         else {
170           HLTError("expecting 0 or 1 as parameter for option %s", argument.Data());
171           iResult=-EINVAL;
172         }
173       } else {
174         HLTError("expecting number as parameter for option %s", argument.Data());
175         iResult=-EINVAL;
176       }
177
178     // -arraysize
179     } else if (argument.CompareTo("-arraysize")==0) {
180       if ((bMissingParam=(++i>=argc))) break;
181       TString param=argv[i];
182       if (param.IsDigit()) {
183         fArraySize=param.Atoi();
184       } else {
185         HLTError("expecting number as parameter for option %s", argument.Data());
186         iResult=-EINVAL;
187       }
188
189     } else {
190       HLTError("unknown argument %s", argument.Data());
191       iResult=-EINVAL;
192     }
193   }
194   if (bMissingParam) {
195     HLTError("missing parameter for argument %s", argument.Data());
196     iResult=-EINVAL;
197   }
198
199   if (iResult>=0) {
200     fpLevelArray=new UInt_t[fArraySize];
201     fpSpecArray=new UInt_t[fArraySize];
202     fpBlockNoArray=new UInt_t[fArraySize];
203     fpIdArray=new UInt_t[fArraySize];
204     fpTimeArray=new UInt_t[fArraySize];
205     fpCTimeArray=new UInt_t[fArraySize];
206     fpInputBlockCountArray=new UInt_t[fArraySize];
207     fpTotalInputSizeArray=new UInt_t[fArraySize];
208     fpNormalizedInputSizeArray=new UInt_t[fArraySize];
209     fpOutputBlockCountArray=new UInt_t[fArraySize];
210     fpTotalOutputSizeArray=new UInt_t[fArraySize];
211     fpInputOutputRatioArray=new UInt_t[fArraySize];
212     fpNormalizedInputOutputRatioArray=new UInt_t[fArraySize];
213     fpComponentCycleTimeArray=new UInt_t[fArraySize];
214     fpEventTypeArray=new UInt_t[fArraySize];
215     fpEventCountArray=new UInt_t[fArraySize];
216
217     fpStatTree=new TTree("CompStat", "HLT component statistics");
218     if (fpStatTree) {
219       fpStatTree->SetDirectory(0);
220       fpStatTree->Branch("cycleTime",        &fCycleTime, "cycleTime/F");
221       fpStatTree->Branch("nofSets",          &fNofSets, "nofSets/I");
222       fpStatTree->Branch("Level",            fpLevelArray, "Level[nofSets]/i");
223       fpStatTree->Branch("Specification",    fpSpecArray, "Specification[nofSets]/i");
224       fpStatTree->Branch("BlockNo",          fpBlockNoArray, "BlockNo[nofSets]/i");
225       fpStatTree->Branch("Id",               fpIdArray, "Id[nofSets]/i");
226       fpStatTree->Branch("Time",             fpTimeArray, "Time[nofSets]/i");
227       fpStatTree->Branch("CTime",            fpCTimeArray, "CTime[nofSets]/i");
228       fpStatTree->Branch("InputBlockCount",  fpInputBlockCountArray, "InputBlockCount[nofSets]/i");
229       fpStatTree->Branch("TotalInputSize",   fpTotalInputSizeArray, "TotalInputSize[nofSets]/i");
230       fpStatTree->Branch("NormalizedInputSize",   fpNormalizedInputSizeArray, "NormalizedInputSize[nofSets]/i");
231       fpStatTree->Branch("OutputBlockCount", fpOutputBlockCountArray, "OutputBlockCount[nofSets]/i");
232       fpStatTree->Branch("TotalOutputSize",  fpTotalOutputSizeArray, "TotalOutputSize[nofSets]/i");
233       fpStatTree->Branch("InputOutputRatio",  fpInputOutputRatioArray, "InputOutputRatio[nofSets]/i");
234       fpStatTree->Branch("NormalizedInputOutputRatio",  fpNormalizedInputOutputRatioArray, "NormalizedInputOutputRatio[nofSets]/i");
235       fpStatTree->Branch("ComponentCycleTime",fpComponentCycleTimeArray, "ComponentCycleTime[nofSets]/i");
236       fpStatTree->Branch("EventType",fpEventTypeArray, "EventType[nofSets]/i");
237       fpStatTree->Branch("EventCount",fpEventCountArray, "EventCount[nofSets]/i");
238     }
239   }
240
241   if (!fFileName.empty()) {
242     fFile=new TFile(fFileName.c_str(), "RECREATE");
243   }
244   return iResult;
245 }
246
247 int AliHLTCompStatCollector::DoDeinit( )
248 {
249   // see header file for class documentation
250   ClearAll();
251
252   if (fFile) {
253     fFile->Close();
254     delete fFile;
255     fFile=NULL;
256   }
257   return 0;
258 }
259
260 int AliHLTCompStatCollector::DoEvent( const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/)
261 {
262   // see header file for class documentation
263   int iResult=0;
264
265   AliHLTUInt32_t eventType=gkAliEventTypeUnknown;
266   IsDataEvent(&eventType);
267
268   ResetFillingVariables();
269   if (fpTimer) {
270     fCycleTime=fpTimer->RealTime()*1000000;
271   }
272
273   bool bEmbeddedTree=false;
274   bool bFolderCreated=false;
275   if ((bFolderCreated=(fpFolder==NULL))) {
276     fpFolder=new TFolder(HLTSTAT_FOLDER_NAME, HLTSTAT_FOLDER_DESC);
277     if (bEmbeddedTree) fpFolder->Add(fpStatTree);
278   }
279   if (!fpFolder) return -ENOMEM;
280   vector<TFolder*> newFolders;
281
282   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeComponentTable);
283        pBlock && iResult>=0;
284        pBlock=GetNextInputBlock()) {
285     string chainId, compId, compArgs;
286     vector<AliHLTUInt32_t> parents;
287     iResult=ExtractComponentTableEntry((const AliHLTUInt8_t*)pBlock->fPtr, pBlock->fSize,
288                                        chainId, compId, compArgs,
289                                        parents);
290     if (iResult>0) {
291       HLTDebug("%s(%s) 0x%08x", chainId.c_str(), compId.c_str(), pBlock->fSpecification);
292       TObject* pObj=NULL;
293       TFolder* pEntry=NULL;
294       if ((pObj=fpFolder->FindObjectAny(chainId.c_str()))!=NULL &&
295           (pEntry=dynamic_cast<TFolder*>(pObj))!=NULL ) {
296         
297       } else if (pObj) {
298         HLTError("entry %s exists in folder, but is not a sub-folder", chainId.c_str());
299       } else if (chainId.size()>0) {
300         pEntry=new TFolder(chainId.c_str(), chainId.c_str());
301         if (pEntry) {
302           pEntry->SetOwner();
303           TFolder* pProps=pEntry->AddFolder(HLTSTAT_ENTRY_PROPS_FOLDER_NAME, HLTSTAT_ENTRY_PROPS_FOLDER_DESC);
304           if (pProps) {
305             pProps->Add(new TObjString(compId.c_str()));
306             if (!compArgs.empty())
307               pProps->Add(new TObjString(compArgs.c_str()));
308             TNamed* pCRC=new TNamed(HLTSTAT_ENTRY_PROPS_IDOBJ_NAME, HLTSTAT_ENTRY_PROPS_IDOBJ_DESC);
309             if (pCRC) {
310               pCRC->SetUniqueID(pBlock->fSpecification);
311               pProps->Add(pCRC);
312             }
313           }
314           TFolder* pParents=pEntry->AddFolder(HLTSTAT_ENTRY_PARENT_FOLDER_NAME, HLTSTAT_ENTRY_PARENT_FOLDER_DESC);
315           if (pParents) {
316             for (vector<AliHLTUInt32_t>::iterator parent=parents.begin();
317                  parent!=parents.end(); parent++) {
318               TString name; name.Form("0x%08x", *parent);
319               pParents->Add(new TObjString(name));
320             }
321           }
322           if (parents.size()==0) {
323             newFolders.push_back(pEntry);
324           } else {
325             vector<TFolder*>::iterator iter=newFolders.begin();
326             vector<AliHLTUInt32_t>::iterator parent=parents.begin();
327             while (iter!=newFolders.end() && parent!=parents.end()) {
328               TObject* idobj=(*iter)->FindObjectAny(HLTSTAT_ENTRY_PROPS_IDOBJ_NAME);
329               AliHLTUInt32_t crcid=0;
330               if (idobj) crcid=idobj->GetUniqueID();
331               HLTDebug("check: %s 0x%08x", (*iter)->GetName(), crcid);
332               if (idobj && crcid==*parent) break;
333               if ((++parent!=parents.end())) continue;
334               parent=parents.begin();
335               iter++;
336             }
337             newFolders.insert(iter,pEntry);
338           }
339         }
340       } else {
341         HLTError("missing chain id for table entry 0x%08x (%p %d), skipping ...", pBlock->fSpecification, pBlock->fPtr, pBlock->fSize);
342       }
343     } else if (iResult!=0) {
344       HLTError("extraction of table entry 0x%08x (%p %d) failed with %d", pBlock->fSpecification, pBlock->fPtr, pBlock->fSize, iResult);
345     }
346     iResult=0;
347   }
348
349   if (newFolders.size()>0) {
350     vector<TFolder*> revert;
351     vector<TFolder*>::iterator iter=newFolders.begin();
352     while (iter!=newFolders.end()) {
353       revert.insert(revert.begin(), *iter);
354       HLTDebug("%s", (*iter)->GetName());
355       iter++;
356     }
357     newFolders.empty();
358     newFolders.assign(revert.begin(), revert.end());
359
360     vector<TFolder*>::iterator publisher=newFolders.begin();
361     while (publisher!=newFolders.end()) {
362       bool bRemove=false;
363       HLTDebug("checking %s for parents", (*publisher)->GetName());
364       TFolder* propsFolder=dynamic_cast<TFolder*>((*publisher)->FindObject(HLTSTAT_ENTRY_PROPS_FOLDER_NAME));
365       assert(propsFolder);
366       TObject* idobj=NULL;
367       if (propsFolder) idobj=propsFolder->FindObject(HLTSTAT_ENTRY_PROPS_IDOBJ_NAME);
368       assert(idobj);
369       AliHLTUInt32_t crcid=idobj->GetUniqueID();
370       TString idstr; idstr.Form("0x%08x", crcid);
371       if (idobj) {
372         for (vector<TFolder*>::iterator consumer=publisher+1;
373              consumer!=newFolders.end(); consumer++) {
374           HLTDebug("   checking %s", (*consumer)->GetName());
375           TFolder* parentFolder=dynamic_cast<TFolder*>((*consumer)->FindObject(HLTSTAT_ENTRY_PARENT_FOLDER_NAME));
376           assert(parentFolder);
377           if (parentFolder) {
378             TIter entries(parentFolder->GetListOfFolders());
379             while (TObject* entry=entries.Next())
380 #ifdef __DEBUG
381               if (entry) {
382                 HLTDebug("   searching %s in %s: %s", idstr.Data(), (*consumer)->GetName(), entry->GetName());
383               }
384 #endif
385             TObject* parent=parentFolder->FindObjectAny(idstr);
386             if (parent) {
387               parentFolder->Add(*publisher);
388               parentFolder->Remove(parent);
389               bRemove=true;
390             }
391           }
392         }
393       }
394       if (bRemove) publisher=newFolders.erase(publisher);
395       else publisher++;
396     }
397
398     for (publisher=newFolders.begin();
399          publisher!=newFolders.end(); publisher++) {
400       RemoveRecurrence(*publisher);
401       fpFolder->Add(*publisher);
402     }
403   }
404
405   int blockNo=0;
406   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeComponentStatistics);
407        pBlock && iResult>=0;
408        pBlock=GetNextInputBlock(), blockNo++) {
409     unsigned int current=fPosition;
410     iResult=FillVariablesSorted(pBlock->fPtr, pBlock->fSize, eventType);
411     for (; current<fPosition; current++) {
412       fpSpecArray[current]=pBlock->fSpecification;
413       fpBlockNoArray[current]=blockNo;
414     }
415     // indicate availability of component statistic block
416     iResult=1;
417   }
418
419   int totalOutputSize=0;
420   if (iResult>0 && eventType) {
421     fNofSets=fPosition;
422     fpStatTree->Fill();
423
424     // init the timer for the next cycle
425     if (!fpTimer)  fpTimer=new TStopwatch;
426     if (fpTimer) {
427       fpTimer->Reset();
428       fpTimer->Start();
429     }
430   }
431
432   if (eventType==gkAliEventTypeEndOfRun ||
433       (iResult>=0 && CheckPeriod())) {
434
435     // publish objects to component output
436     if ((fMode&kPublishObjects)!=0) {
437       if (!bEmbeddedTree) {
438         iResult=PushBack(fpStatTree, kAliHLTDataTypeTTree|kAliHLTDataOriginOut);
439         totalOutputSize+=GetLastObjectSize();
440       }
441       iResult=PushBack(fpFolder, kAliHLTDataTypeTObject|kAliHLTDataOriginOut);
442       totalOutputSize+=GetLastObjectSize();
443     }
444
445     // save objects to file
446     if ((fMode&kSaveObjects)!=0 && fFile!=NULL) {
447       HLTDebug("saving objects to file %s", fFileName.c_str());
448       fFile->cd();
449       if (!bEmbeddedTree) {
450         fpStatTree->Write("", TObject::kOverwrite);
451       }
452       fpFolder->Write("", TObject::kOverwrite);
453     }
454   }
455
456   if (iResult==-ENOSPC) {
457     fSizeEstimator+=totalOutputSize;
458   }
459
460   if (iResult>0) iResult=0;
461   return iResult;
462 }
463
464 void AliHLTCompStatCollector::ResetFillingVariables()
465 {
466   // see header file for class documentation
467   fCycleTime=0;
468   fNofSets=0;
469   fPosition=0;
470   memset(fpLevelArray, 0, sizeof(UInt_t)*fArraySize);
471   memset(fpSpecArray, 0, sizeof(UInt_t)*fArraySize);
472   memset(fpBlockNoArray, 0, sizeof(UInt_t)*fArraySize);
473   memset(fpIdArray, 0, sizeof(UInt_t)*fArraySize);
474   memset(fpTimeArray, 0, sizeof(UInt_t)*fArraySize);
475   memset(fpCTimeArray, 0, sizeof(UInt_t)*fArraySize);
476   memset(fpInputBlockCountArray, 0, sizeof(UInt_t)*fArraySize);
477   memset(fpTotalInputSizeArray, 0, sizeof(UInt_t)*fArraySize);
478   memset(fpNormalizedInputSizeArray, 0, sizeof(UInt_t)*fArraySize);
479   memset(fpOutputBlockCountArray, 0, sizeof(UInt_t)*fArraySize);
480   memset(fpTotalOutputSizeArray, 0, sizeof(UInt_t)*fArraySize);
481   memset(fpInputOutputRatioArray, 0, sizeof(UInt_t)*fArraySize);
482   memset(fpNormalizedInputOutputRatioArray, 0, sizeof(UInt_t)*fArraySize);
483   memset(fpComponentCycleTimeArray, 0, sizeof(UInt_t)*fArraySize);
484   memset(fpEventTypeArray, 0, sizeof(UInt_t)*fArraySize);
485   memset(fpEventCountArray, 0, sizeof(UInt_t)*fArraySize);
486 }
487
488 int AliHLTCompStatCollector::FillVariablesSorted(void* ptr, int size, AliHLTUInt32_t eventType)
489 {
490   // see header file for class documentation
491   int iResult=0;
492   if (size%sizeof(AliHLTComponentStatistics)) {
493     // older or invalid structure
494     HLTError("data block is not aligned to the size of the AliHLTComponentStatistics struct");
495     return -EINVAL;
496   }
497   
498   AliHLTComponentStatistics* pStat=reinterpret_cast<AliHLTComponentStatistics*>(ptr);
499   UInt_t nofStats=size/sizeof(AliHLTComponentStatistics);
500   vector<int> indexList;
501   UInt_t i=0;
502   for (i=0; i<nofStats; i++) {
503     vector<int>::iterator element=indexList.begin();
504     for (; element!=indexList.end(); element++) {
505       if (pStat[i].fLevel>pStat[*element].fLevel) {
506         break;
507       }
508     }
509     indexList.insert(element, i);
510   }
511
512   i=fPosition;
513   for (vector<int>::iterator element=indexList.begin();
514        element!=indexList.end();
515        element++, i++) {
516     if (i<fArraySize) {
517       fpLevelArray[i]=pStat[*element].fLevel;
518       fpIdArray[i]=pStat[*element].fId;
519       fpTimeArray[i]=pStat[*element].fTime;
520       fpCTimeArray[i]=pStat[*element].fCTime;
521       fpInputBlockCountArray[i]=pStat[*element].fInputBlockCount;
522       fpTotalInputSizeArray[i]=pStat[*element].fTotalInputSize;
523       fpNormalizedInputSizeArray[i]=pStat[*element].fTotalInputSize;
524       if (pStat[*element].fInputBlockCount>0) 
525         fpNormalizedInputSizeArray[i]/=pStat[*element].fInputBlockCount;
526       fpOutputBlockCountArray[i]=pStat[*element].fOutputBlockCount;
527       fpTotalOutputSizeArray[i]=pStat[*element].fTotalOutputSize;
528       if (pStat[*element].fTotalOutputSize>0)
529         fpInputOutputRatioArray[i]=pStat[*element].fTotalInputSize/pStat[*element].fTotalOutputSize;
530       if (pStat[*element].fInputBlockCount>0)
531         fpNormalizedInputOutputRatioArray[i]=fpInputOutputRatioArray[i]*pStat[*element].fOutputBlockCount/pStat[*element].fInputBlockCount;
532       fpComponentCycleTimeArray[i]=pStat[*element].fComponentCycleTime;
533       fpEventTypeArray[i]=eventType;
534       fpEventCountArray[i]=GetEventCount();
535     } else {
536       // TODO: dynamically grow arrays with placement new
537     }
538   }
539
540   if (i>=fArraySize) {
541     HLTWarning("too little space in branch variables to fill %d statistics blocks, available %d at position %d", i, fArraySize, fPosition);
542     fPosition=fArraySize;
543   } else {
544     fPosition=i;
545   }
546   
547   return iResult;
548 }
549
550 void AliHLTCompStatCollector::ClearAll()
551 {
552   // see header file for class documentation
553   if (fpTimer) delete fpTimer; fpTimer=NULL;
554   if (fpFolder) delete fpFolder; fpFolder=NULL;
555   if (fpStatTree) delete fpStatTree; fpStatTree=NULL;
556   if (fpLevelArray) delete fpLevelArray; fpLevelArray=NULL;
557   if (fpSpecArray) delete fpSpecArray; fpSpecArray=NULL;
558   if (fpBlockNoArray) delete fpBlockNoArray; fpBlockNoArray=NULL;
559   if (fpIdArray) delete fpIdArray; fpIdArray=NULL;
560   if (fpTimeArray) delete fpTimeArray; fpTimeArray=NULL;
561   if (fpCTimeArray) delete fpCTimeArray; fpCTimeArray=NULL;
562   if (fpInputBlockCountArray) delete fpInputBlockCountArray; fpInputBlockCountArray=NULL;
563   if (fpTotalInputSizeArray) delete fpTotalInputSizeArray; fpTotalInputSizeArray=NULL;
564   if (fpNormalizedInputSizeArray) delete fpNormalizedInputSizeArray; fpNormalizedInputSizeArray=NULL;
565   if (fpOutputBlockCountArray) delete fpOutputBlockCountArray; fpOutputBlockCountArray=NULL;
566   if (fpTotalOutputSizeArray) delete fpTotalOutputSizeArray; fpTotalOutputSizeArray=NULL;
567   if (fpInputOutputRatioArray) delete fpInputOutputRatioArray; fpInputOutputRatioArray=NULL;
568   if (fpNormalizedInputOutputRatioArray) delete fpNormalizedInputOutputRatioArray; fpNormalizedInputOutputRatioArray=NULL;
569   if (fpComponentCycleTimeArray) delete fpComponentCycleTimeArray; fpComponentCycleTimeArray=NULL;
570   if (fpEventTypeArray) delete fpEventTypeArray; fpEventTypeArray=NULL;
571   if (fpEventCountArray) delete fpEventCountArray; fpEventCountArray=NULL;
572 }
573
574 int AliHLTCompStatCollector::RemoveRecurrence(TFolder* pRoot) const
575 {
576   // see header file for class documentation
577   int iResult=0;
578   if (!pRoot) return -EINVAL;
579   TFolder* parentFolder=dynamic_cast<TFolder*>(pRoot->FindObject(HLTSTAT_ENTRY_PARENT_FOLDER_NAME));
580   assert(parentFolder);
581   vector<TFolder*> listRemove;
582   if (parentFolder) {
583     TIter entries(parentFolder->GetListOfFolders());
584     TFolder* entry=NULL;
585     TObject* obj=NULL;
586     while ((obj=entries.Next())!=NULL && (entry=dynamic_cast<TFolder*>(obj))!=NULL) {
587       TString name=entry->GetName();
588       HLTDebug("checking %s for recurrence", name.Data());
589       TIter tokens(parentFolder->GetListOfFolders());
590       TFolder* token=NULL;
591       while ((obj=tokens.Next())!=NULL && (token=dynamic_cast<TFolder*>(obj))!=NULL) {
592         if (name.CompareTo(token->GetName())==0) continue;
593         if ((obj=token->FindObjectAny(name))!=NULL) {
594           listRemove.push_back(entry);
595           HLTDebug("found recurrence in %s", token->GetName());
596           break;
597         } else {
598           HLTDebug("no recurrence found in %s", token->GetName());
599         }
600       }
601       RemoveRecurrence(entry);
602     }
603     for (vector<TFolder*>::iterator removeElement=listRemove.begin();
604          removeElement!=listRemove.end(); removeElement++) {
605       parentFolder->Remove(*removeElement);
606     }
607   }
608   
609   return iResult;  
610 }
611
612 bool AliHLTCompStatCollector::CheckPeriod(bool bUpdate)
613 {
614   // see header file for class documentation
615   bool result=true;
616   if (fEventModulo>0) {
617     if ((result=((GetEventCount()+1)%fEventModulo)==0)) {
618       return true;
619     }
620   }
621   if (fPeriod>0) {
622     if ((result=((difftime(time(NULL), fLastTime)>(double)fPeriod))) &&
623         bUpdate) {
624       fLastTime=time(NULL);
625     }
626   }
627   return result;
628 }