]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
expanding compnent statistic structure to store cycle time
authorrichterm <richterm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 14 Jan 2009 13:41:28 +0000 (13:41 +0000)
committerrichterm <richterm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 14 Jan 2009 13:41:28 +0000 (13:41 +0000)
HLT/BASE/AliHLTComponent.cxx
HLT/BASE/AliHLTDataTypes.h
HLT/BASE/util/AliHLTCompStatCollector.cxx
HLT/BASE/util/AliHLTCompStatCollector.h

index 0daf1c7e81c60cecafb3b6298d11ac5de52ef32e..60fbe8107c18709674fbeba2ef80753d17ee39ff 100644 (file)
@@ -45,6 +45,7 @@ using namespace std;
  * default compression level for ROOT objects
  */
 #define ALIHLTCOMPONENT_DEFAULT_OBJECT_COMPRESSION 5
+#define ALIHLTCOMPONENT_STATTIME_SCALER 1000000
 
 /** ROOT macro for the implementation of ROOT specific class methods */
 ClassImp(AliHLTComponent);
@@ -240,9 +241,6 @@ int AliHLTComponent::Init(const AliHLTAnalysisEnvironment* comenv, void* environ
 #if defined(__DEBUG) || defined(HLT_COMPONENT_STATISTICS)
   // benchmarking stopwatch for the component statistics
   fpBenchmark=new TStopwatch;
-  if (fpBenchmark) {
-    fpBenchmark->Start();
-  }
 #endif
 
   return iResult;
@@ -1411,12 +1409,15 @@ int AliHLTComponent::ProcessEvent( const AliHLTComponentEventData& evtData,
 #if defined(__DEBUG) || defined(HLT_COMPONENT_STATISTICS)
   AliHLTComponentStatistics outputStat;
   memset(&outputStat, 0, sizeof(AliHLTComponentStatistics));
+  outputStat.fStructSize=sizeof(AliHLTComponentStatistics);
   outputStat.fId=fChainIdCrc;
-  compStats.push_back(outputStat);
   if (fpBenchmark) {
+    fpBenchmark->Stop();
+    outputStat.fComponentCycleTime=(AliHLTUInt32_t)(fpBenchmark->RealTime()*ALIHLTCOMPONENT_STATTIME_SCALER);
     fpBenchmark->Reset();
     fpBenchmark->Start();
   }
+  compStats.push_back(outputStat);
 #endif
 
   // data processing is skipped
@@ -1670,8 +1671,10 @@ int  AliHLTComponent::AddComponentStatistics(AliHLTComponentBlockDataList& block
   stats[0].fTotalOutputSize=offset;
   stats[0].fOutputBlockCount=blocks.size();
   if (fpBenchmark) {
-    stats[0].fTime=(AliHLTUInt32_t)(fpBenchmark->RealTime()*1000000);
-    stats[0].fCTime=(AliHLTUInt32_t)(fpBenchmark->CpuTime()*1000000);
+    fpBenchmark->Stop();
+    stats[0].fTime=(AliHLTUInt32_t)(fpBenchmark->RealTime()*ALIHLTCOMPONENT_STATTIME_SCALER);
+    stats[0].fCTime=(AliHLTUInt32_t)(fpBenchmark->CpuTime()*ALIHLTCOMPONENT_STATTIME_SCALER);
+    fpBenchmark->Continue();
   }
   if (offset+stats.size()*sizeof(AliHLTComponentStatistics)<=bufferSize) {
     AliHLTComponentBlockData bd;
index ef2410b8546d968a1ee334ad29dcafc69861ba90..f66434d6819b4e670c48bcf115fb7f60456836d2 100644 (file)
  *           kAliHLTMCObjectDataType added
  *           kAliHLTDataOriginOffline added
  *           kAliHLTDataOriginHLT added
+ *  11       extended AliHLTComponentStatistics: one more member to store the
+ *           cycle time between events per component.
  */
-#define ALIHLT_DATA_TYPES_VERSION 10
+#define ALIHLT_DATA_TYPES_VERSION 11
 
 //////////////////////////////////////////////////////////////////////////
 //
@@ -562,6 +564,9 @@ extern "C" {
    * fLevel is retrieved from incoming block statistics and incremented. 
    * Incoming block statistics are appended to the newly added one if
    * --enable-compstat=full has been chosen.
+   *
+   * ChangeLog:
+   *   2009-01-14 fComponentCycleTime added
    */
   struct AliHLTComponentStatistics
   {
@@ -574,6 +579,7 @@ extern "C" {
     AliHLTUInt32_t fTotalInputSize;
     AliHLTUInt32_t fOutputBlockCount;
     AliHLTUInt32_t fTotalOutputSize;
+    AliHLTUInt32_t fComponentCycleTime;
   };
 
   /**
index 2ec8876e8d54573d222278210781224da74b455c..bb2cd9641e44d26de6ebb4cb60290819f58f01d8 100644 (file)
@@ -66,6 +66,9 @@ AliHLTCompStatCollector::AliHLTCompStatCollector()
   fpTotalInputSizeArray(NULL),
   fpOutputBlockCountArray(NULL),
   fpTotalOutputSizeArray(NULL)
+  , fpComponentCycleTimeArray(NULL)
+  , fpEventTypeArray(NULL)
+  , fpEventCountArray(NULL)
   , fSizeEstimator(1000)
   , fMode(kPublishObjects)
   , fFileName()
@@ -189,6 +192,9 @@ int AliHLTCompStatCollector::DoInit( int argc, const char** argv )
     fpTotalInputSizeArray=new UInt_t[fArraySize];
     fpOutputBlockCountArray=new UInt_t[fArraySize];
     fpTotalOutputSizeArray=new UInt_t[fArraySize];
+    fpComponentCycleTimeArray=new UInt_t[fArraySize];
+    fpEventTypeArray=new UInt_t[fArraySize];
+    fpEventCountArray=new UInt_t[fArraySize];
 
     fpStatTree=new TTree("CompStat", "HLT component statistics");
     if (fpStatTree) {
@@ -205,6 +211,9 @@ int AliHLTCompStatCollector::DoInit( int argc, const char** argv )
       fpStatTree->Branch("TotalInputSize",   fpTotalInputSizeArray, "TotalInputSize[nofSets]/i");
       fpStatTree->Branch("OutputBlockCount", fpOutputBlockCountArray, "OutputBlockCount[nofSets]/i");
       fpStatTree->Branch("TotalOutputSize",  fpTotalOutputSizeArray, "TotalOutputSize[nofSets]/i");
+      fpStatTree->Branch("ComponentCycleTime",fpComponentCycleTimeArray, "ComponentCycleTime[nofSets]/i");
+      fpStatTree->Branch("EventType",fpEventTypeArray, "EventType[nofSets]/i");
+      fpStatTree->Branch("EventCount",fpEventCountArray, "EventCount[nofSets]/i");
     }
   }
 
@@ -376,7 +385,7 @@ int AliHLTCompStatCollector::DoEvent( const AliHLTComponentEventData& /*evtData*
        pBlock && iResult>=0;
        pBlock=GetNextInputBlock(), blockNo++) {
     unsigned int current=fPosition;
-    iResult=FillVariablesSorted(pBlock->fPtr, pBlock->fSize);
+    iResult=FillVariablesSorted(pBlock->fPtr, pBlock->fSize, eventType);
     for (; current<fPosition; current++) {
       fpSpecArray[current]=pBlock->fSpecification;
       fpBlockNoArray[current]=blockNo;
@@ -386,7 +395,7 @@ int AliHLTCompStatCollector::DoEvent( const AliHLTComponentEventData& /*evtData*
   }
 
   int totalOutputSize=0;
-  if (iResult>0) {
+  if (iResult>0 && eventType) {
     fNofSets=fPosition;
     fpStatTree->Fill();
 
@@ -446,16 +455,21 @@ void AliHLTCompStatCollector::ResetFillingVariables()
   memset(fpTotalInputSizeArray, 0, sizeof(UInt_t)*fArraySize);
   memset(fpOutputBlockCountArray, 0, sizeof(UInt_t)*fArraySize);
   memset(fpTotalOutputSizeArray, 0, sizeof(UInt_t)*fArraySize);
+  memset(fpComponentCycleTimeArray, 0, sizeof(UInt_t)*fArraySize);
+  memset(fpEventTypeArray, 0, sizeof(UInt_t)*fArraySize);
+  memset(fpEventCountArray, 0, sizeof(UInt_t)*fArraySize);
 }
 
-int AliHLTCompStatCollector::FillVariablesSorted(void* ptr, int size)
+int AliHLTCompStatCollector::FillVariablesSorted(void* ptr, int size, AliHLTUInt32_t eventType)
 {
   // see header file for class documentation
   int iResult=0;
   if (size%sizeof(AliHLTComponentStatistics)) {
+    // older or invalid structure
     HLTError("data block is not aligned to the size of the AliHLTComponentStatistics struct");
     return -EINVAL;
   }
+  
   AliHLTComponentStatistics* pStat=reinterpret_cast<AliHLTComponentStatistics*>(ptr);
   UInt_t nofStats=size/sizeof(AliHLTComponentStatistics);
   vector<int> indexList;
@@ -483,6 +497,9 @@ int AliHLTCompStatCollector::FillVariablesSorted(void* ptr, int size)
       fpTotalInputSizeArray[i]=pStat[*element].fTotalInputSize;
       fpOutputBlockCountArray[i]=pStat[*element].fOutputBlockCount;
       fpTotalOutputSizeArray[i]=pStat[*element].fTotalOutputSize;
+      fpComponentCycleTimeArray[i]=pStat[*element].fComponentCycleTime;
+      fpEventTypeArray[i]=eventType;
+      fpEventCountArray[i]=GetEventCount();
     } else {
       // TODO: dynamically grow arrays with placement new
     }
@@ -514,6 +531,9 @@ void AliHLTCompStatCollector::ClearAll()
   if (fpTotalInputSizeArray) delete fpTotalInputSizeArray; fpTotalInputSizeArray=NULL;
   if (fpOutputBlockCountArray) delete fpOutputBlockCountArray; fpOutputBlockCountArray=NULL;
   if (fpTotalOutputSizeArray) delete fpTotalOutputSizeArray; fpTotalOutputSizeArray=NULL;
+  if (fpComponentCycleTimeArray) delete fpComponentCycleTimeArray; fpComponentCycleTimeArray=NULL;
+  if (fpEventTypeArray) delete fpEventTypeArray; fpEventTypeArray=NULL;
+  if (fpEventCountArray) delete fpEventCountArray; fpEventCountArray=NULL;
 }
 
 int AliHLTCompStatCollector::RemoveRecurrence(TFolder* pRoot) const
index ca510b65a82673b1a17ec22b4de00b42d3186ed8..74fba1251ec9338bdcbe8891e0883f3e12547f1f 100644 (file)
@@ -124,7 +124,7 @@ class AliHLTCompStatCollector : public AliHLTProcessor
   /**
    * Fill the lists from the component statistics block.
    */
-  int FillVariablesSorted(void* ptr, int size);
+  int FillVariablesSorted(void* ptr, int size, AliHLTUInt32_t eventType);
 
   /** delete all internal objects */
   void ClearAll();
@@ -181,6 +181,12 @@ class AliHLTCompStatCollector : public AliHLTProcessor
   UInt_t* fpOutputBlockCountArray; //!transient
   /** branch filling variable */
   UInt_t* fpTotalOutputSizeArray; //!transient
+  /** branch filling variable */
+  UInt_t* fpComponentCycleTimeArray; //!transient
+  /** branch filling variable */
+  UInt_t* fpEventTypeArray; //!transient
+  /** branch filling variable */
+  UInt_t* fpEventCountArray; //!transient
 
   /** const base of GetOutputSize, updated on error in DoEvent */
   int fSizeEstimator; //! transient
@@ -203,6 +209,6 @@ class AliHLTCompStatCollector : public AliHLTProcessor
   /** event modulo to save/publish onjects */
   unsigned int fEventModulo; //! transient
 
-  ClassDef(AliHLTCompStatCollector, 3)
+  ClassDef(AliHLTCompStatCollector, 4)
 };
 #endif