]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Changes for using the functionality of AliPHOSRawDecoder
authorpolicheh <policheh@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 6 Feb 2007 11:02:00 +0000 (11:02 +0000)
committerpolicheh <policheh@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 6 Feb 2007 11:02:00 +0000 (11:02 +0000)
and AliPHOSRawDigiProducer

PHOS/AliPHOSCalibHistoProducer.cxx
PHOS/AliPHOSGetter.cxx

index a4396877141668264333f260abcc3ecd742b8a2b..0b8e6b0b3216526fb3b140fead4edd99c3a8f2ef 100644 (file)
 #include "TH1.h"
 #include "TFile.h"
 #include "AliRawReader.h"
-#include "AliCaloRawStream.h"
+#include "AliPHOSDigit.h"
+#include "AliPHOSRawDecoder.h"
+#include "AliPHOSRawDigiProducer.h"
+#include "AliPHOSGeometry.h"
+#include "TClonesArray.h"
 
 ClassImp(AliPHOSCalibHistoProducer)
 
@@ -79,89 +83,58 @@ void AliPHOSCalibHistoProducer::Run()
   // Reads raw data stream and fills amplitude histograms
   // The histograms are written to file every fUpdatingRate events
 
-  TH1F* gHighGain = 0;
-  TH1F* gLowGain = 0;
-  Int_t iBin = 0;
   Int_t iEvent = 0;
-  Int_t runNum = 0;
+  Int_t runNum = -111;
+  Int_t relId[4];
+  AliPHOSDigit* digit;
+  Double_t energy;
 
-  AliCaloRawStream in(fRawReader,"PHOS");
+  TClonesArray *digits = new TClonesArray("AliPHOSDigit",100);
+  AliPHOSGeometry* geo = AliPHOSGeometry::GetInstance("IHEP");
+
+  AliPHOSRawDecoder dc(fRawReader);
   if(fIsOldRCUFormat)
-    in.SetOldRCUFormat(kTRUE);
+    dc.SetOldRCUFormat(kTRUE);
 
   // Read raw data event by event
 
   while (fRawReader->NextEvent()) {
     runNum = fRawReader->GetRunNumber();
 
-    while ( in.Next() ) { 
-       
-      if(!gHighGain) gHighGain = new TH1F("gHighGain","High gain",
-                                         in.GetTimeLength(),0,in.GetTimeLength());
-      else
-       if(gHighGain->GetNbinsX() != in.GetTimeLength()) {
-         delete gHighGain;
-         gHighGain = new TH1F("gHighGain","High gain",in.GetTimeLength(),0,in.GetTimeLength());
-       }
-
-      if(!gLowGain)  gLowGain = new TH1F("gLowGain","Low gain",
-                                        in.GetTimeLength(),0,in.GetTimeLength());
-      else
-       if(gLowGain->GetNbinsX() != in.GetTimeLength()) {
-         delete gLowGain;
-         gLowGain = new TH1F("gLowGain","Low gain",in.GetTimeLength(),0,in.GetTimeLength());
-       }
-
-      Bool_t lowGainFlag = in.IsLowGain();
-      
-      if(lowGainFlag) 
-       gLowGain->SetBinContent(in.GetTimeLength()-iBin-1,in.GetSignal());
+    AliPHOSRawDigiProducer producer;
+    producer.MakeDigits(digits,&dc);
+    
+    for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
+      digit = (AliPHOSDigit*)digits->At(iDigit);
+      energy = digit->GetEnergy(); // no pedestal subtraction!
+      geo->AbsToRelNumbering(digit->GetId(),relId);        
+      Int_t mod = relId[0];
+      Int_t col = relId[3];
+      Int_t row = relId[2];
+
+      if(fAmpHisto[mod][col][row]) {
+       fAmpHisto[mod][col][row]->Fill(energy);
+      }
       else {
-       gHighGain->SetBinContent(in.GetTimeLength()-iBin-1,in.GetSignal());
+       char hname[128];
+       sprintf(hname,"mod%dcol%drow%d",mod,col,row);
+       fAmpHisto[mod][col][row] = new TH1F(hname,hname,100,0.,1000.);
       }
 
-      iBin++;
-
-      if(iBin==in.GetTimeLength()) {
-       iBin=0;
-
-       Double_t energy;
-
-       if(!lowGainFlag) {
-         energy = gHighGain->GetMaximum(); // no pedestal subtraction!
-       }
-       else {
-         energy = gLowGain->GetMaximum(); // no pedestal subtraction!
-       }
-           
-       Int_t mod = in.GetModule();
-       Int_t col = in.GetColumn();
-       Int_t row = in.GetRow();
-
-       if(fAmpHisto[mod][col][row]) {
-         if(!lowGainFlag) {
-           fAmpHisto[mod][col][row]->Fill(energy);
-         }
-       }
-       else {
-         char hname[128];
-         sprintf(hname,"mod%dcol%drow%d",mod,col,row);
-         fAmpHisto[mod][col][row] = new TH1F(hname,hname,100,0.,1000.);
-       }
-
-
-      }
+      // update histograms in local file every 100th event
+      if(iEvent%fUpdatingRate == 0) {
+       AliInfo(Form("Updating histo file, event %d, run %d\n",iEvent,runNum));
+       UpdateHistoFile();
+      } 
     }
 
-    // update histograms in local file every 100th event
-    if(iEvent%fUpdatingRate == 0) {
-      AliInfo(Form("Updating histo file, event %d, run %d\n",iEvent,runNum));
-      UpdateHistoFile();
-    } 
     iEvent++;
   }
 
-  UpdateHistoFile(); 
+  UpdateHistoFile();
+  digits->Delete();
+  delete digits;
+  delete geo;
   AliInfo(Form("%d events of run %d processed.",iEvent,runNum));
 }
 
index a6c31a61aa83b5c5e9c5a7043e63a96e2fd8104f..5ac90f3e29deddfac0ca026832032443d1639db7 100644 (file)
@@ -69,6 +69,7 @@
 #include "AliCDBLocal.h"
 #include "AliCDBStorage.h"
 #include "AliCDBManager.h"
+#include "AliPHOSRawDigiProducer.h"
 
 ClassImp(AliPHOSGetter)
   
@@ -744,73 +745,9 @@ Int_t AliPHOSGetter::ReadRaw(AliRawReader *rawReader,Bool_t isOldRCUFormat)
   dc.SetOldRCUFormat(isOldRCUFormat);
   dc.SubtractPedestals(kTRUE);
 
-  Int_t relId[4], absId =0;
-  Bool_t lowGainFlag = kFALSE ; 
-
   TClonesArray * digits = Digits() ;
-  digits->Clear() ;
-
-  Int_t    iDigit   = 0 ; 
-  Double_t time     = 0. ;
-  Int_t    iOldDigit;
-  Bool_t   seen;    
-
-  while (dc.NextDigit()) {  
-
-    lowGainFlag = dc.IsLowGain();
-    time = dc.GetTime();
-
-    relId[0] = dc.GetModule();
-    relId[1] = 0;
-    relId[2] = dc.GetRow();
-    relId[3] = dc.GetColumn();
-
-    if(!PHOSGeometry()) AliFatal("Couldn't find PHOSGeometry!");
-    PHOSGeometry()->RelToAbsNumbering(relId, absId);
-    AliDebug(2,Form("relId=(mod,row,col)=(%d,%d,%d), absId=%d\n\n",
-                   relId[0],relId[2],relId[3],absId));
-
-    // Add low gain digit only 
-    //if the high gain digit does not exist in the digits array
-
-    seen = kFALSE;
-    
-    if(lowGainFlag) {
-      for (iOldDigit=iDigit-1; iOldDigit>=0; iOldDigit--) {
-       if ((dynamic_cast<AliPHOSDigit*>(digits->At(iOldDigit)))->GetId() == absId) {
-         seen = kTRUE;
-         break;
-       }
-      }
-      if (!seen) {
-       new((*digits)[iDigit]) AliPHOSDigit(-1,absId,(Float_t)dc.GetEnergy(),time);
-       iDigit++;
-      }
-    }
-
-    // Add high gain digit only if it is not saturated;
-    // replace low gain digit by a high gain one
-    else {
-      if (dc.GetEnergy() >= 1023) continue;
-      for (iOldDigit=iDigit-1; iOldDigit>=0; iOldDigit--) {
-       if ((dynamic_cast<AliPHOSDigit*>(digits->At(iOldDigit)))->GetId() == absId) {
-         digits->RemoveAt(iOldDigit);
-         new((*digits)[iOldDigit]) AliPHOSDigit(-1,absId,(Float_t)dc.GetEnergy(),time);
-         seen = kTRUE;
-         break;
-       }
-      }
-      if (!seen) {
-       new((*digits)[iDigit]) AliPHOSDigit(-1,absId,(Float_t)dc.GetEnergy(),time);
-       iDigit++;
-      }
-    }
-    
-    
-  }
-  
-  digits->Compress() ;
-  digits->Sort() ;
+  AliPHOSRawDigiProducer pr;
+  pr.MakeDigits(digits,&dc);
 
   //!!!!for debug!!!
   Int_t modMax=-111;
@@ -819,7 +756,8 @@ Int_t AliPHOSGetter::ReadRaw(AliRawReader *rawReader,Bool_t isOldRCUFormat)
   Float_t eMax=-333;
   //!!!for debug!!!
 
-  for(iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
+  Int_t relId[4];
+  for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
     AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(iDigit);
     if(digit->GetEnergy()>eMax) {
       PHOSGeometry()->AbsToRelNumbering(digit->GetId(),relId);