]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - VZERO/AliVZEROReconstructor.cxx
bugfix: corrected defines to use right default algorithms
[u/mrichter/AliRoot.git] / VZERO / AliVZEROReconstructor.cxx
index c4da8bff44fc5731e808ec2cd40571ea3405eb71..1d345a8cc8c7da7acc1de938c187c683ba0005c9 100644 (file)
@@ -25,7 +25,7 @@
 #include "AliRawReader.h"
 #include "AliVZEROReconstructor.h"
 #include "AliVZERORawStream.h"
-#include "AliESD.h"
+#include "AliESDEvent.h"
 
 ClassImp(AliVZEROReconstructor)
 
@@ -33,7 +33,6 @@ ClassImp(AliVZEROReconstructor)
 AliVZEROReconstructor:: AliVZEROReconstructor(): AliReconstructor(),
    fESDVZERO(0x0),
    fESD(0x0),
-   fRunLoader(0x0),
    fCalibData(GetCalibData())
 {
   // Default constructor  
@@ -62,17 +61,18 @@ AliVZEROReconstructor::~AliVZEROReconstructor()
 }
 
 //_____________________________________________________________________________
-void AliVZEROReconstructor::Init(AliRunLoader* runLoader)
+void AliVZEROReconstructor::Init()
 {
-/// initializer
+// initializer
 
-  fRunLoader = runLoader;
   fESDVZERO  = new AliESDVZERO;
 }
 
 //______________________________________________________________________
 void AliVZEROReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
 {
+// converts to digits
+
   if (!digitsTree) {
     AliError("No digits tree!");
     return;
@@ -83,12 +83,13 @@ void AliVZEROReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digits
 
   rawReader->Reset();
   AliVZERORawStream rawStream(rawReader);
-  while (rawStream.Next()) {
-    Int_t pmNumber = rawStream.GetCell();
-    Int_t adc = rawStream.GetADC();  
-    Int_t time = rawStream.GetTime();
+  if (rawStream.Next()) {
+    for(Int_t iChannel = 0; iChannel < 64; iChannel++) {
+    Int_t adc = rawStream.GetADC(iChannel);  
+    Int_t time = rawStream.GetTime(iChannel);
     new ((*digitsArray)[digitsArray->GetEntriesFast()])
-      AliVZEROdigit(pmNumber,adc,time);
+      AliVZEROdigit(iChannel,adc,time);
+    }
   }
 
   digitsTree->Fill();
@@ -96,8 +97,10 @@ void AliVZEROReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digits
 
 //______________________________________________________________________
 void AliVZEROReconstructor::FillESD(TTree* digitsTree, TTree* /*clustersTree*/,
-                                   AliESD* esd) const
+                                   AliESDEvent* esd) const
 {
+// fills multiplicities to the ESD
+
   if (!digitsTree) {
     AliError("No digits tree!");
     return;
@@ -107,27 +110,27 @@ void AliVZEROReconstructor::FillESD(TTree* digitsTree, TTree* /*clustersTree*/,
   TBranch* digitBranch = digitsTree->GetBranch("VZERODigit");
   digitBranch->SetAddress(&digitsArray);
 
-  Int_t   NbPMV0A = 0;
-  Int_t   NbPMV0C = 0;
-  Int_t   MTotV0A = 0;
-  Int_t   MTotV0C = 0;
-  Float_t ADCV0A  = 0.0;
-  Float_t ADCV0C  = 0.0;
-  Float_t MultV0A[4];
-  Float_t MultV0C[4];
-  Int_t   MRingV0A[4];
-  Int_t   MRingV0C[4];
+  Int_t   nbPMV0A = 0;
+  Int_t   nbPMV0C = 0;
+  Int_t   mTotV0A = 0;
+  Int_t   mTotV0C = 0;
+  Float_t adcV0A  = 0.0;
+  Float_t adcV0C  = 0.0;
+  Float_t multV0A[4];
+  Float_t multV0C[4];
+  Int_t   mRingV0A[4];
+  Int_t   mRingV0C[4];
   
-  Int_t   ADC[64]; 
-  Float_t MIP[64];
+  Int_t   adc[64]; 
+  Float_t mip[64];
   for (Int_t i=0; i<64; i++){
-       ADC[i] = 0;
-       MIP[i] = 110.0;}
+       adc[i] = 0;
+       mip[i] = 110.0;}
   for (Int_t j=0; j<4; j++){
-       MultV0A[j]  = 0.0;
-       MultV0C[j]  = 0.0;
-       MRingV0A[j] = 0;
-       MRingV0C[j] = 0;}
+       multV0A[j]  = 0.0;
+       multV0C[j]  = 0.0;
+       mRingV0A[j] = 0;
+       mRingV0C[j] = 0;}
      
   // loop over VZERO entries
   Int_t nEntries = (Int_t)digitsTree->GetEntries();
@@ -138,43 +141,46 @@ void AliVZEROReconstructor::FillESD(TTree* digitsTree, TTree* /*clustersTree*/,
     
     for (Int_t d=0; d<nDigits; d++) {    
       AliVZEROdigit* digit = (AliVZEROdigit*)digitsArray->At(d);      
-      Int_t  PMNumber      = digit->PMNumber();  
-      ADC[PMNumber] = digit->ADC();  
-      if (PMNumber<=31) {
-        if (PMNumber<=7) MultV0C[0]=MultV0C[0]+ float(ADC[PMNumber])/MIP[PMNumber];
-       if (PMNumber>=8  && PMNumber<=15) MultV0C[1]=MultV0C[1]+ float(ADC[PMNumber])/MIP[PMNumber];
-       if (PMNumber>=16 && PMNumber<=23) MultV0C[2]=MultV0C[2]+ float(ADC[PMNumber])/MIP[PMNumber];
-       if (PMNumber>=24 && PMNumber<=31) MultV0C[3]=MultV0C[3]+ float(ADC[PMNumber])/MIP[PMNumber];
-        ADCV0C = ADCV0C + float(ADC[PMNumber])/MIP[PMNumber];
-       if(ADC[PMNumber] > 4) NbPMV0C++;
-      }        
-      if (PMNumber>=32) {
-        if (PMNumber>=32 && PMNumber<=39) MultV0A[0]=MultV0A[0]+ float(ADC[PMNumber])/MIP[PMNumber];
-       if (PMNumber>=40 && PMNumber<=47) MultV0A[1]=MultV0A[1]+ float(ADC[PMNumber])/MIP[PMNumber];
-       if (PMNumber>=48 && PMNumber<=55) MultV0A[2]=MultV0A[2]+ float(ADC[PMNumber])/MIP[PMNumber];
-       if (PMNumber>=56 && PMNumber<=63) MultV0A[3]=MultV0A[3]+ float(ADC[PMNumber])/MIP[PMNumber];
-        ADCV0A = ADCV0A + float(ADC[PMNumber])/MIP[PMNumber];
-       if(ADC[PMNumber] > 4) NbPMV0A++;
+      Int_t  pmNumber      = digit->PMNumber();  
+      adc[pmNumber] = digit->ADC(); 
+      // cut of ADC at MIP/2
+      if  (adc[pmNumber] > (mip[pmNumber]/2)) { 
+        if (pmNumber<=31) {
+          if (pmNumber<=7) multV0C[0]=multV0C[0]+ float(adc[pmNumber])/mip[pmNumber];
+         if (pmNumber>=8  && pmNumber<=15) multV0C[1]=multV0C[1]+ float(adc[pmNumber])/mip[pmNumber];
+         if (pmNumber>=16 && pmNumber<=23) multV0C[2]=multV0C[2]+ float(adc[pmNumber])/mip[pmNumber];
+         if (pmNumber>=24 && pmNumber<=31) multV0C[3]=multV0C[3]+ float(adc[pmNumber])/mip[pmNumber];
+          adcV0C = adcV0C + float(adc[pmNumber])/mip[pmNumber];
+         nbPMV0C++;
+        }      
+        if (pmNumber>=32 ) {
+          if (pmNumber>=32 && pmNumber<=39) multV0A[0]=multV0A[0]+ float(adc[pmNumber])/mip[pmNumber];
+         if (pmNumber>=40 && pmNumber<=47) multV0A[1]=multV0A[1]+ float(adc[pmNumber])/mip[pmNumber];
+         if (pmNumber>=48 && pmNumber<=55) multV0A[2]=multV0A[2]+ float(adc[pmNumber])/mip[pmNumber];
+         if (pmNumber>=56 && pmNumber<=63) multV0A[3]=multV0A[3]+ float(adc[pmNumber])/mip[pmNumber];
+          adcV0A = adcV0A + float(adc[pmNumber])/mip[pmNumber];
+         nbPMV0A++;
+        }
       }
     } // end of loop over digits
     
   } // end of loop over events in digits tree
   
-  MTotV0A = int(ADCV0A + 0.5);
-  MTotV0C = int(ADCV0C + 0.5);
+  mTotV0A = int(adcV0A + 0.5);
+  mTotV0C = int(adcV0C + 0.5);
   for (Int_t j=0; j<4; j++){       
-       MRingV0A[j] = int(MultV0A[j] + 0.5);
-       MRingV0C[j] = int(MultV0C[j] + 0.5);}
+       mRingV0A[j] = int(multV0A[j] + 0.5);
+       mRingV0C[j] = int(multV0C[j] + 0.5);}
      
-  AliDebug(1,Form("VZERO multiplicities : %d (V0A) %d (V0C)", MTotV0A, MTotV0C));
-  AliDebug(1,Form("Number of PMs fired  : %d (V0A) %d (V0C)", NbPMV0A, NbPMV0C));
-
-  fESDVZERO->SetNbPMV0A(NbPMV0A);
-  fESDVZERO->SetNbPMV0C(NbPMV0C);
-  fESDVZERO->SetMTotV0A(MTotV0A);
-  fESDVZERO->SetMTotV0C(MTotV0C);
-  fESDVZERO->SetMRingV0A(MRingV0A);
-  fESDVZERO->SetMRingV0C(MRingV0C);
+  AliDebug(1,Form("VZERO multiplicities : %d (V0A) %d (V0C)", mTotV0A, mTotV0C));
+  AliDebug(1,Form("Number of PMs fired  : %d (V0A) %d (V0C)", nbPMV0A, nbPMV0C));
+
+  fESDVZERO->SetNbPMV0A(nbPMV0A);
+  fESDVZERO->SetNbPMV0C(nbPMV0C);
+  fESDVZERO->SetMTotV0A(mTotV0A);
+  fESDVZERO->SetMTotV0C(mTotV0C);
+  fESDVZERO->SetMRingV0A(mRingV0A);
+  fESDVZERO->SetMRingV0C(mRingV0C);
   
   if (esd) { 
     AliDebug(1, Form("Writing VZERO data to ESD tree"));
@@ -185,7 +191,8 @@ void AliVZEROReconstructor::FillESD(TTree* digitsTree, TTree* /*clustersTree*/,
 //_____________________________________________________________________________
 AliCDBStorage* AliVZEROReconstructor::SetStorage(const char *uri) 
 {
-  
+// Sets the storage  
+
   Bool_t deleteManager = kFALSE;
   
   AliCDBManager *manager = AliCDBManager::Instance();
@@ -210,8 +217,7 @@ AliCDBStorage* AliVZEROReconstructor::SetStorage(const char *uri)
 //_____________________________________________________________________________
 AliVZEROCalibData* AliVZEROReconstructor::GetCalibData() const
 {
-
-  // Getting calibration object for VZERO set
+  // Gets calibration object for VZERO set
 
   AliCDBManager *man = AliCDBManager::Instance();
 
@@ -219,22 +225,21 @@ AliVZEROCalibData* AliVZEROReconstructor::GetCalibData() const
 
   entry = man->Get("VZERO/Calib/Data");
 
-  if(!entry){
-    AliWarning("Load of calibration data from default storage failed!");
-    AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
-    Int_t runNumber = man->GetRun();
-    entry = man->GetStorage("local://$ALICE_ROOT")
-      ->Get("VZERO/Calib/Data",runNumber);
-       
-  }
+//   if(!entry){
+//     AliWarning("Load of calibration data from default storage failed!");
+//     AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
+//     Int_t runNumber = man->GetRun();
+//     entry = man->GetStorage("local://$ALICE_ROOT")
+//       ->Get("VZERO/Calib/Data",runNumber);
+//     
+//   }
 
   // Retrieval of data in directory VZERO/Calib/Data:
 
-
   AliVZEROCalibData *calibdata = 0;
 
   if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
-  if (!calibdata)  AliError("No calibration data from calibration database !");
+  if (!calibdata)  AliFatal("No calibration data from calibration database !");
 
   return calibdata;
 }