]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/BASE/util/AliHLTCompStatCollector.cxx
corrected implementation of component description data blocks and extraction of the...
[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 "TStopwatch.h"
27 #include "TH1F.h"
28 #include "TH2F.h"
29 #include "TH2C.h"
30 #include "TTree.h"
31 #include "TFolder.h"
32 #include "TNamed.h"
33 #include "TString.h"
34 #include <cassert>
35
36 #define HLTSTAT_FOLDER_NAME               "HLTstat"
37 #define HLTSTAT_FOLDER_DESC               "ALICE HLT component statistics"
38 #define HLTSTAT_ENTRY_PARENT_FOLDER_NAME  "parents"
39 #define HLTSTAT_ENTRY_PARENT_FOLDER_DESC  "parent components"
40 #define HLTSTAT_ENTRY_PROPS_FOLDER_NAME   "props"
41 #define HLTSTAT_ENTRY_PROPS_FOLDER_DESC   "component properties"
42 #define HLTSTAT_ENTRY_PROPS_IDOBJ_NAME    "id"
43 #define HLTSTAT_ENTRY_PROPS_IDOBJ_DESC    "numerical id calculated from chain id"
44
45 /** ROOT macro for the implementation of ROOT specific class methods */
46 ClassImp(AliHLTCompStatCollector)
47
48 AliHLTCompStatCollector::AliHLTCompStatCollector()
49   :
50   AliHLTProcessor(),
51   fpTimer(NULL),
52   fpFolder(NULL),
53   fpStatTree(NULL),
54   fCycleTime(0),
55   fNofSets(0),
56   fArraySize(1000),
57   fPosition(0),
58   fpLevelArray(NULL),
59   fpSpecArray(NULL),
60   fpBlockNoArray(NULL),
61   fpIdArray(NULL),
62   fpTimeArray(NULL),
63   fpCTimeArray(NULL),
64   fpInputBlockCountArray(NULL),
65   fpTotalInputSizeArray(NULL),
66   fpOutputBlockCountArray(NULL),
67   fpTotalOutputSizeArray(NULL)
68 {
69   // see header file for class documentation
70   // or
71   // refer to README to build package
72   // or
73   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
74 }
75
76 AliHLTCompStatCollector::~AliHLTCompStatCollector()
77 {
78   // see header file for class documentation
79   ClearAll();
80 }
81
82 void AliHLTCompStatCollector::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
83 {
84   // see header file for class documentation
85   list.push_back(kAliHLTDataTypeComponentStatistics);
86 }
87
88 AliHLTComponentDataType AliHLTCompStatCollector::GetOutputDataType()
89 {
90   // see header file for class documentation
91   return kAliHLTMultipleDataType;
92 }
93
94 int AliHLTCompStatCollector::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
95 {
96   // see header file for class documentation
97   tgtList.clear();
98   tgtList.push_back(kAliHLTDataTypeHistogram);
99   tgtList.push_back(kAliHLTDataTypeTTree);
100   return tgtList.size();
101 }
102
103 void AliHLTCompStatCollector::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
104 {
105   // see header file for class documentation
106   constBase=1000;
107   inputMultiplier=100.0;
108 }
109
110 int AliHLTCompStatCollector::DoInit( int argc, const char** argv )
111 {
112   // see header file for class documentation
113   int iResult=0;
114   TString argument="";
115   int bMissingParam=0;
116   for (int i=0; i<argc && iResult>=0; i++) {
117     argument=argv[i];
118     if (argument.IsNull()) continue;
119
120     // -array-size
121     if (argument.CompareTo("-array-size")==0) {
122       if ((bMissingParam=(++i>=argc))) break;
123       TString param=argv[argc];
124       if (param.IsDigit()) {
125         fArraySize=param.Atoi();
126       } else {
127         HLTError("expecting number as parameter for option '-array-size'");
128         iResult=-EINVAL;
129       }
130
131     } else {
132       iResult=-EINVAL;
133     }
134   }
135   if (bMissingParam) {
136     HLTError("missing parameter for argument %s", argument.Data());
137     iResult=-EINVAL;
138   }
139
140   if (iResult>=0) {
141     fpLevelArray=new UInt_t[fArraySize];
142     fpSpecArray=new UInt_t[fArraySize];
143     fpBlockNoArray=new UInt_t[fArraySize];
144     fpIdArray=new UInt_t[fArraySize];
145     fpTimeArray=new UInt_t[fArraySize];
146     fpCTimeArray=new UInt_t[fArraySize];
147     fpInputBlockCountArray=new UInt_t[fArraySize];
148     fpTotalInputSizeArray=new UInt_t[fArraySize];
149     fpOutputBlockCountArray=new UInt_t[fArraySize];
150     fpTotalOutputSizeArray=new UInt_t[fArraySize];
151
152     fpStatTree=new TTree("CompStat", "HLT component statistics");
153     if (fpStatTree) {
154       fpStatTree->SetDirectory(0);
155       fpStatTree->Branch("cycleTime",        &fCycleTime, "cycleTime/F");
156       fpStatTree->Branch("nofSets",          &fNofSets, "nofSets/I");
157       fpStatTree->Branch("Level",            fpLevelArray, "Level[nofSets]/i");
158       fpStatTree->Branch("Specification",    fpSpecArray, "Specification[nofSets]/i");
159       fpStatTree->Branch("BlockNo",          fpBlockNoArray, "BlockNo[nofSets]/i");
160       fpStatTree->Branch("Id",               fpIdArray, "Id[nofSets]/i");
161       fpStatTree->Branch("Time",             fpTimeArray, "Time[nofSets]/i");
162       fpStatTree->Branch("CTime",            fpCTimeArray, "CTime[nofSets]/i");
163       fpStatTree->Branch("InputBlockCount",  fpInputBlockCountArray, "InputBlockCount[nofSets]/i");
164       fpStatTree->Branch("TotalInputSize",   fpTotalInputSizeArray, "TotalInputSize[nofSets]/i");
165       fpStatTree->Branch("OutputBlockCount", fpOutputBlockCountArray, "OutputBlockCount[nofSets]/i");
166       fpStatTree->Branch("TotalOutputSize",  fpTotalOutputSizeArray, "TotalOutputSize[nofSets]/i");
167     }
168   }
169   return iResult;
170 }
171
172 int AliHLTCompStatCollector::DoDeinit( )
173 {
174   // see header file for class documentation
175   ClearAll();
176
177   return 0;
178 }
179
180 int AliHLTCompStatCollector::DoEvent( const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/)
181 {
182   // see header file for class documentation
183   int iResult=0;
184
185   AliHLTUInt32_t eventType=gkAliEventTypeUnknown;
186   IsDataEvent(&eventType);
187
188   ResetFillingVariables();
189   if (fpTimer) {
190     fCycleTime=fpTimer->RealTime()*1000000;
191   }
192
193   bool bEmbeddedTree=false;
194   bool bFolderCreated=false;
195   if (bFolderCreated=(fpFolder==NULL)) {
196     fpFolder=new TFolder(HLTSTAT_FOLDER_NAME, HLTSTAT_FOLDER_DESC);
197     if (bEmbeddedTree) fpFolder->Add(fpStatTree);
198   }
199   if (!fpFolder) return -ENOMEM;
200   vector<TFolder*> newFolders;
201
202   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeComponentTable);
203        pBlock && iResult>=0;
204        pBlock=GetNextInputBlock()) {
205     string chainId, compId, compArgs;
206     vector<AliHLTUInt32_t> parents;
207     iResult=ExtractComponentTableEntry((const AliHLTUInt8_t*)pBlock->fPtr, pBlock->fSize,
208                                        chainId, compId, compArgs,
209                                        parents);
210     if (iResult>0) {
211       HLTDebug("%s(%s) 0x%08x", chainId.c_str(), compId.c_str(), pBlock->fSpecification);
212       TObject* pObj=NULL;
213       TFolder* pEntry=NULL;
214       if ((pObj=fpFolder->FindObjectAny(chainId.c_str()))!=NULL &&
215           (pEntry=dynamic_cast<TFolder*>(pObj))!=NULL ) {
216         
217       } else if (pObj) {
218         HLTError("entry %s exists in folder, but is not a sub-folder", chainId.c_str());
219       } else if (chainId.size()>0) {
220         pEntry=new TFolder(chainId.c_str(), chainId.c_str());
221         if (pEntry) {
222           pEntry->SetOwner();
223           TFolder* pProps=pEntry->AddFolder(HLTSTAT_ENTRY_PROPS_FOLDER_NAME, HLTSTAT_ENTRY_PROPS_FOLDER_DESC);
224           if (pProps) {
225             pProps->Add(new TObjString(compId.c_str()));
226             if (!compArgs.empty())
227               pProps->Add(new TObjString(compArgs.c_str()));
228             TNamed* pCRC=new TNamed(HLTSTAT_ENTRY_PROPS_IDOBJ_NAME, HLTSTAT_ENTRY_PROPS_IDOBJ_DESC);
229             if (pCRC) {
230               pCRC->SetUniqueID(pBlock->fSpecification);
231               pProps->Add(pCRC);
232             }
233           }
234           TFolder* pParents=pEntry->AddFolder(HLTSTAT_ENTRY_PARENT_FOLDER_NAME, HLTSTAT_ENTRY_PARENT_FOLDER_DESC);
235           if (pParents) {
236             for (vector<AliHLTUInt32_t>::iterator parent=parents.begin();
237                  parent!=parents.end(); parent++) {
238               TString name; name.Form("0x%08x", *parent);
239               pParents->Add(new TObjString(name));
240             }
241           }
242           if (parents.size()==0) {
243             newFolders.push_back(pEntry);
244           } else {
245             vector<TFolder*>::iterator iter=newFolders.begin();
246             vector<AliHLTUInt32_t>::iterator parent=parents.begin();
247             while (iter!=newFolders.end() && parent!=parents.end()) {
248               TObject* idobj=(*iter)->FindObjectAny(HLTSTAT_ENTRY_PROPS_IDOBJ_NAME);
249               AliHLTUInt32_t crcid=0;
250               if (idobj) crcid=idobj->GetUniqueID();
251               HLTDebug("check: %s 0x%08x", (*iter)->GetName(), crcid);
252               if (idobj && crcid==*parent) break;
253               if ((++parent!=parents.end())) continue;
254               parent=parents.begin();
255               iter++;
256             }
257             newFolders.insert(iter,pEntry);
258           }
259         }
260       } else {
261         HLTError("missing chain id for table entry 0x%08x (%p %d), skipping ...", pBlock->fSpecification, pBlock->fPtr, pBlock->fSize);
262       }
263     } else if (iResult!=0) {
264       HLTError("extraction of table entry 0x%08x (%p %d) failed with %d", pBlock->fSpecification, pBlock->fPtr, pBlock->fSize, iResult);
265     }
266     iResult=0;
267   }
268
269   if (newFolders.size()>0) {
270     vector<TFolder*> revert;
271     vector<TFolder*>::iterator iter=newFolders.begin();
272     while (iter!=newFolders.end()) {
273       revert.insert(revert.begin(), *iter);
274       HLTDebug("%s", (*iter)->GetName());
275       iter++;
276     }
277     newFolders.empty();
278     newFolders.assign(revert.begin(), revert.end());
279
280     vector<TFolder*>::iterator publisher=newFolders.begin();
281     while (publisher!=newFolders.end()) {
282       bool bRemove=false;
283       HLTDebug("checking %s for parents", (*publisher)->GetName());
284       TFolder* propsFolder=dynamic_cast<TFolder*>((*publisher)->FindObject(HLTSTAT_ENTRY_PROPS_FOLDER_NAME));
285       assert(propsFolder);
286       TObject* idobj=NULL;
287       if (propsFolder) idobj=propsFolder->FindObject(HLTSTAT_ENTRY_PROPS_IDOBJ_NAME);
288       assert(idobj);
289       AliHLTUInt32_t crcid=idobj->GetUniqueID();
290       TString idstr; idstr.Form("0x%08x", crcid);
291       if (idobj) {
292         for (vector<TFolder*>::iterator consumer=publisher+1;
293              consumer!=newFolders.end(); consumer++) {
294           HLTDebug("   checking %s", (*consumer)->GetName());
295           TFolder* parentFolder=dynamic_cast<TFolder*>((*consumer)->FindObject(HLTSTAT_ENTRY_PARENT_FOLDER_NAME));
296           assert(parentFolder);
297           if (parentFolder) {
298             TIter entries(parentFolder->GetListOfFolders());
299             while (TObject* entry=entries.Next())
300               HLTDebug("   searching %s in %s: %s", idstr.Data(), (*consumer)->GetName(), entry->GetName());
301
302             TObject* parent=parentFolder->FindObjectAny(idstr);
303             if (parent) {
304               parentFolder->Add(*publisher);
305               parentFolder->Remove(parent);
306               bRemove=true;
307             }
308           }
309         }
310       }
311       if (bRemove) publisher=newFolders.erase(publisher);
312       else publisher++;
313     }
314
315     for (publisher=newFolders.begin();
316          publisher!=newFolders.end(); publisher++) {
317       RemoveRecurrence(*publisher);
318       fpFolder->Add(*publisher);
319     }
320   }
321
322   int blockNo=0;
323   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeComponentStatistics);
324        pBlock && iResult>=0;
325        pBlock=GetNextInputBlock(), blockNo++) {
326     unsigned int current=fPosition;
327     iResult=FillVariablesSorted(pBlock->fPtr, pBlock->fSize);
328     for (; current<fPosition; current++) {
329       fpSpecArray[current]=pBlock->fSpecification;
330       fpBlockNoArray[current]=blockNo;
331     }
332     // indicate availability of component statistic block
333     iResult=1;
334   }
335
336   if (iResult>0) {
337     fNofSets=fPosition;
338     fpStatTree->Fill();
339     if (!bEmbeddedTree)
340       iResult=PushBack(fpStatTree, kAliHLTDataTypeTTree|kAliHLTDataOriginOut);
341
342     // init the timer for the next cycle
343     if (!fpTimer)  fpTimer=new TStopwatch;
344     if (fpTimer) {
345       fpTimer->Reset();
346       fpTimer->Start();
347     }
348   }
349
350   if (iResult>=0 /*&& eventType==gkAliEventTypeEndOfRun*/) {
351     PushBack(fpFolder, kAliHLTDataTypeTObject|kAliHLTDataOriginOut);
352   }
353
354   if (iResult>0) iResult=0;
355   return iResult;
356 }
357
358 void AliHLTCompStatCollector::ResetFillingVariables()
359 {
360   // see header file for class documentation
361   fCycleTime=0;
362   fNofSets=0;
363   fPosition=0;
364   memset(fpLevelArray, 0, sizeof(UInt_t)*fArraySize);
365   memset(fpSpecArray, 0, sizeof(UInt_t)*fArraySize);
366   memset(fpBlockNoArray, 0, sizeof(UInt_t)*fArraySize);
367   memset(fpIdArray, 0, sizeof(UInt_t)*fArraySize);
368   memset(fpTimeArray, 0, sizeof(UInt_t)*fArraySize);
369   memset(fpCTimeArray, 0, sizeof(UInt_t)*fArraySize);
370   memset(fpInputBlockCountArray, 0, sizeof(UInt_t)*fArraySize);
371   memset(fpTotalInputSizeArray, 0, sizeof(UInt_t)*fArraySize);
372   memset(fpOutputBlockCountArray, 0, sizeof(UInt_t)*fArraySize);
373   memset(fpTotalOutputSizeArray, 0, sizeof(UInt_t)*fArraySize);
374 }
375
376 int AliHLTCompStatCollector::FillVariablesSorted(void* ptr, int size)
377 {
378   // see header file for class documentation
379   int iResult=0;
380   if (size%sizeof(AliHLTComponentStatistics)) {
381     HLTError("data block is not aligned to the size of the AliHLTComponentStatistics struct");
382     return -EINVAL;
383   }
384   AliHLTComponentStatistics* pStat=reinterpret_cast<AliHLTComponentStatistics*>(ptr);
385   UInt_t nofStats=size/sizeof(AliHLTComponentStatistics);
386   vector<int> indexList;
387   UInt_t i=0;
388   for (i=0; i<nofStats; i++) {
389     vector<int>::iterator element=indexList.begin();
390     for (; element!=indexList.end(); element++) {
391       if (pStat[i].fLevel>pStat[*element].fLevel) {
392         break;
393       }
394     }
395     indexList.insert(element, i);
396   }
397
398   i=fPosition;
399   for (vector<int>::iterator element=indexList.begin();
400        element!=indexList.end();
401        element++, i++) {
402     if (i<fArraySize) {
403       fpLevelArray[i]=pStat[*element].fLevel;
404       fpIdArray[i]=pStat[*element].fId;
405       fpTimeArray[i]=pStat[*element].fTime;
406       fpCTimeArray[i]=pStat[*element].fCTime;
407       fpInputBlockCountArray[i]=pStat[*element].fInputBlockCount;
408       fpTotalInputSizeArray[i]=pStat[*element].fTotalInputSize;
409       fpOutputBlockCountArray[i]=pStat[*element].fOutputBlockCount;
410       fpTotalOutputSizeArray[i]=pStat[*element].fTotalOutputSize;
411     }
412   }
413
414   if (i>=fArraySize) {
415     HLTWarning("too little space in branch variables to fill %d statistics blocks, available %d at position %d", i, fArraySize, fPosition);
416     fPosition=fArraySize;
417   } else {
418     fPosition=i;
419   }
420   
421   return iResult;
422 }
423
424 void AliHLTCompStatCollector::ClearAll()
425 {
426   // see header file for class documentation
427   if (fpTimer) delete fpTimer; fpTimer=NULL;
428   if (fpFolder) delete fpFolder; fpFolder=NULL;
429   if (fpStatTree) delete fpStatTree; fpStatTree=NULL;
430   if (fpLevelArray) delete fpLevelArray; fpLevelArray=NULL;
431   if (fpSpecArray) delete fpSpecArray; fpSpecArray=NULL;
432   if (fpBlockNoArray) delete fpBlockNoArray; fpBlockNoArray=NULL;
433   if (fpIdArray) delete fpIdArray; fpIdArray=NULL;
434   if (fpTimeArray) delete fpTimeArray; fpTimeArray=NULL;
435   if (fpCTimeArray) delete fpCTimeArray; fpCTimeArray=NULL;
436   if (fpInputBlockCountArray) delete fpInputBlockCountArray; fpInputBlockCountArray=NULL;
437   if (fpTotalInputSizeArray) delete fpTotalInputSizeArray; fpTotalInputSizeArray=NULL;
438   if (fpOutputBlockCountArray) delete fpOutputBlockCountArray; fpOutputBlockCountArray=NULL;
439   if (fpTotalOutputSizeArray) delete fpTotalOutputSizeArray; fpTotalOutputSizeArray=NULL;
440 }
441
442 int AliHLTCompStatCollector::RemoveRecurrence(TFolder* pRoot) const
443 {
444   // see header file for class documentation
445   int iResult=0;
446   if (!pRoot) return -EINVAL;
447   TFolder* parentFolder=dynamic_cast<TFolder*>(pRoot->FindObject(HLTSTAT_ENTRY_PARENT_FOLDER_NAME));
448   assert(parentFolder);
449   vector<TFolder*> listRemove;
450   if (parentFolder) {
451     TIter entries(parentFolder->GetListOfFolders());
452     TFolder* entry=NULL;
453     TObject* obj=NULL;
454     while ((obj=entries.Next())!=NULL && (entry=dynamic_cast<TFolder*>(obj))!=NULL) {
455       TString name=entry->GetName();
456       HLTDebug("checking %s for recurrence", name.Data());
457       TIter tokens(parentFolder->GetListOfFolders());
458       TFolder* token=NULL;
459       while ((obj=tokens.Next())!=NULL && (token=dynamic_cast<TFolder*>(obj))!=NULL) {
460         if (name.CompareTo(token->GetName())==0) continue;
461         if ((obj=token->FindObjectAny(name))!=NULL) {
462           listRemove.push_back(entry);
463           HLTDebug("found recurrence in %s", token->GetName());
464           break;
465         } else {
466           HLTDebug("no recurrence found in %s", token->GetName());
467         }
468       }
469       RemoveRecurrence(entry);
470     }
471     for (vector<TFolder*>::iterator removeElement=listRemove.begin();
472          removeElement!=listRemove.end(); removeElement++) {
473       parentFolder->Remove(*removeElement);
474     }
475   }
476   
477   return iResult;  
478 }