]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - HLT/VZERO/AliHLTVZERORecoComponent.cxx
Update master to aliroot
[u/mrichter/AliRoot.git] / HLT / VZERO / AliHLTVZERORecoComponent.cxx
index f5eca8b0b62fbdd6e40c51046d21b9b26bcc4cfb..30003de3d07786f0642b9e3743cdb899c02574c2 100644 (file)
@@ -1,5 +1,4 @@
-//-*- Mode: C++ -*-
-// $Id: AliHLTVZERORecoComponent.cxx $
+// $Id$
 /**************************************************************************
  * This file is property of and copyright by the ALICE HLT Project        * 
  * ALICE Experiment at CERN, All rights reserved.                         *
     @brief   VZERO reconstruction component
 */
 
-#if __GNUC__>= 3
-using namespace std;
-#endif
-
 #include "TTree.h"
 #include "TMap.h"
 #include "TObjString.h"
+#include "TDatime.h"
 
 #include "AliLog.h"
 #include "AliRunInfo.h"
@@ -42,6 +38,8 @@ using namespace std;
 #include "AliHLTDataTypes.h"
 #include "AliHLTVZERORecoComponent.h"
 
+using namespace std;
+
 /** ROOT macro for the implementation of ROOT specific class methods */
 ClassImp(AliHLTVZERORecoComponent)
 
@@ -54,8 +52,7 @@ ClassImp(AliHLTVZERORecoComponent)
 // #################################################################################
 AliHLTVZERORecoComponent::AliHLTVZERORecoComponent() :
   AliHLTProcessor(),
-  fRunInfo(NULL),
-  fDigitsTree(NULL),
+  fRunInfo(NULL),  
   fVZERORecoParam(NULL),
   fVZEROReconstructor(NULL),
   fRawReader(NULL) {
@@ -131,6 +128,8 @@ void AliHLTVZERORecoComponent::GetOCDBObjectDescription( TMap* const targetMap)
                 new TObjString("VZERO calibration object"));
   targetMap->Add(new TObjString("VZERO/Calib/TimeSlewing"),
                 new TObjString("VZERO calibration object"));
+  targetMap->Add(new TObjString("VZERO/Trigger/Data"),
+                new TObjString("VZERO calibration object"));
 
   return;
 }
@@ -212,12 +211,6 @@ Int_t AliHLTVZERORecoComponent::DoInit( Int_t argc, const Char_t** argv ) {
       break;
     }
 
-    fDigitsTree = new TTree("D", "Digits Tree");
-    if (!fDigitsTree) {
-      iResult=-ENOMEM;
-      break;
-    }
-
     fVZERORecoParam = new AliVZERORecoParam;
     if (!fVZERORecoParam) {
       iResult=-ENOMEM;
@@ -240,10 +233,6 @@ Int_t AliHLTVZERORecoComponent::DoInit( Int_t argc, const Char_t** argv ) {
       delete fRawReader;
     fRawReader = NULL;
 
-    if (!fDigitsTree)
-      delete fDigitsTree;
-    fDigitsTree = NULL;
-
     if (fVZERORecoParam)
       delete fVZERORecoParam;
     fVZERORecoParam = NULL;
@@ -293,10 +282,6 @@ Int_t AliHLTVZERORecoComponent::DoDeinit() {
     delete fRawReader;
   fRawReader = NULL;
   
-  if (!fDigitsTree)
-    delete fDigitsTree;
-  fDigitsTree = NULL;
-    
   if (fVZERORecoParam)
     delete fVZERORecoParam;
   fVZERORecoParam = NULL;
@@ -326,7 +311,7 @@ Int_t AliHLTVZERORecoComponent::DoEvent(const AliHLTComponentEventData& /*evtDat
   // -- Get VZERO raw dat a input block and set up the rawreader
   const AliHLTComponentBlockData* pBlock = GetFirstInputBlock(kAliHLTDataTypeDDLRaw|kAliHLTDataOriginVZERO);
   if (!pBlock) {
-    HLTInfo("No VZERO input block !!!");
+    ALIHLTERRORGUARD(1, "No VZERO input block at event %d", GetEventCount());
     return 0;
   }
   
@@ -337,24 +322,38 @@ Int_t AliHLTVZERORecoComponent::DoEvent(const AliHLTComponentEventData& /*evtDat
     iResult = -1;
   }
   
+  TTree *digitsTree = new TTree("D", "Digits Tree");
+  if (!digitsTree) {
+    iResult=-ENOMEM;
+  }
+
   if (iResult >= 0) {
 
     // -- Set VZERO EquipmentID
     fRawReader->SetEquipmentID(3584);
   
     // -- 1. step VZERO reconstruction
-    fVZEROReconstructor->ConvertDigits(fRawReader, fDigitsTree);
+    fVZEROReconstructor->ConvertDigits(fRawReader, digitsTree);
 
     // -- 2. step VZERO reconstruction -- fill AliESDVZERO object
-    fVZEROReconstructor->FillESD(fDigitsTree, NULL, NULL);
+    fVZEROReconstructor->FillESD(digitsTree, NULL, NULL);
+
+    AliESDVZERO *esdVZERO = fVZEROReconstructor->GetESDVZERO();
+    
+    // Send info every 10 s
+    const TDatime time;    
+    static UInt_t lastTime=0;
+    if (time.Get()-lastTime>10) {
+      lastTime=time.Get();
+      HLTInfo("VZERO Multiplicity A %f - C %f", esdVZERO->GetMTotV0A(), esdVZERO->GetMTotV0A() );
+    }
 
     // -- Send AliESDVZERO
-    PushBack(static_cast<TObject*>(fVZEROReconstructor->GetESDVZERO()),
-            kAliHLTDataTypeESDContent|kAliHLTDataOriginVZERO,0);
+    PushBack(esdVZERO, kAliHLTDataTypeESDContent|kAliHLTDataOriginVZERO,0);
   }
   
   // -- Clean up
-  fDigitsTree->Reset();
+  delete digitsTree;
   fRawReader->ClearBuffers();   
 
   return iResult;
@@ -382,6 +381,5 @@ Int_t AliHLTVZERORecoComponent::Reconfigure(const Char_t* cdbEntry, const Char_t
 // #################################################################################
 Int_t AliHLTVZERORecoComponent::ReadPreprocessorValues(const Char_t* /*modules*/) {
   // see header file for class documentation
-  ALIHLTERRORGUARD(5, "ReadPreProcessorValues not implemented for this component");
   return 0;
 }