]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
New FillESD() for raw data is added
authorkharlov <kharlov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 15 Nov 2006 16:05:03 +0000 (16:05 +0000)
committerkharlov <kharlov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 15 Nov 2006 16:05:03 +0000 (16:05 +0000)
PHOS/AliPHOSReconstructor.cxx
PHOS/AliPHOSReconstructor.h

index 79036e5bdfb7329b4de26f94196b4db73f1f429f..13d5c53f5e7c9e2d4179ce8d70f8a5b9713e2c10 100644 (file)
@@ -117,7 +117,9 @@ void AliPHOSReconstructor::FillESD(AliRunLoader* runLoader, AliESD* esd) const
   Int_t nOfRecParticles = recParticles->GetEntries();
 
   esd->SetNumberOfPHOSClusters(nOfRecParticles) ; 
-  esd->SetFirstPHOSCluster(esd->GetNumberOfCaloClusters()) ; 
+  esd->SetFirstPHOSCluster(esd->GetNumberOfCaloClusters()) ;
+  
+  AliDebug(2,Form("%d digits and %d rec. particles in event %d, option %s",gime->Digits()->GetEntries(),nOfRecParticles,eventNumber,GetOption()));
 
   for (Int_t recpart = 0 ; recpart < nOfRecParticles ; recpart++) {
     AliPHOSRecParticle * rp = dynamic_cast<AliPHOSRecParticle*>(recParticles->At(recpart));
@@ -132,7 +134,9 @@ void AliPHOSReconstructor::FillESD(AliRunLoader* runLoader, AliESD* esd) const
     Float_t xyz[3];
     for (Int_t ixyz=0; ixyz<3; ixyz++) 
       xyz[ixyz] = rp->GetPos()[ixyz];
-
+    
+    AliDebug(2,Form("Global position xyz=(%f,%f,%f)",xyz[0],xyz[1],xyz[2]));
+    
     Int_t  digitMult  = emcRP->GetDigitsMultiplicity();
     Int_t *digitsList = emcRP->GetDigitsList();
     UShort_t *amplList  = new UShort_t[digitMult];
@@ -164,8 +168,93 @@ void AliPHOSReconstructor::FillESD(AliRunLoader* runLoader, AliESD* esd) const
 
     // add the track to the esd object
     esd->AddCaloCluster(ec);
-    delete ec;
+    delete ec;    
+  }  
+}
+
+void AliPHOSReconstructor::FillESD(AliRunLoader* runLoader,
+                                  AliRawReader* rawReader, AliESD* esd) const
+{
+  //This function creates AliESDtracks from AliPHOSRecParticles 
+  //and writes them to the ESD in the case of raw data reconstruction.
+
+  Int_t eventNumber = runLoader->GetEventNumber() ;
+
+  if(eventNumber==0) {
+    rawReader->RewindEvents();
+    rawReader->NextEvent();
   }
+
+  AliPHOSGetter *gime = AliPHOSGetter::Instance();
+
+  Bool_t isOldRCUFormat = kFALSE;
+  TString opt = GetOption();
+  if(opt.Contains("OldRCUFormat"))
+    isOldRCUFormat = kTRUE;
+
+  gime->ReadRaw(rawReader,isOldRCUFormat) ;
+
+  TClonesArray *recParticles  = gime->RecParticles();
+  Int_t nOfRecParticles = recParticles->GetEntries();
+
+  esd->SetNumberOfPHOSClusters(nOfRecParticles) ; 
+  esd->SetFirstPHOSCluster(esd->GetNumberOfCaloClusters()) ;
+  
+  AliDebug(2,Form("%d digits and %d rec. particles in event %d, option %s",gime->Digits()->GetEntries(),nOfRecParticles,eventNumber,GetOption()));
+
+  for (Int_t recpart = 0 ; recpart < nOfRecParticles ; recpart++) {
+    AliPHOSRecParticle * rp = dynamic_cast<AliPHOSRecParticle*>(recParticles->At(recpart));
+
+    if(rp) {
+    Float_t xyz[3];
+    for (Int_t ixyz=0; ixyz<3; ixyz++) 
+      xyz[ixyz] = rp->GetPos()[ixyz];
+
+    AliDebug(2,Form("Global position xyz=(%f,%f,%f)",xyz[0],xyz[1],xyz[2]));
+    
+    AliPHOSTrackSegment *ts    = gime->TrackSegment(rp->GetPHOSTSIndex());
+    AliPHOSEmcRecPoint  *emcRP = gime->EmcRecPoint(ts->GetEmcIndex());
+    AliESDCaloCluster   *ec    = new AliESDCaloCluster() ; 
+
+    Int_t  digitMult  = emcRP->GetDigitsMultiplicity();
+    Int_t *digitsList = emcRP->GetDigitsList();
+    UShort_t *amplList  = new UShort_t[digitMult];
+    UShort_t *digiList  = new UShort_t[digitMult];
+
+    // Convert Float_t* and Int_t* to UShort_t* to save memory
+    for (Int_t iDigit=0; iDigit<digitMult; iDigit++) {
+      AliPHOSDigit *digit = gime->Digit(digitsList[iDigit]);
+      if(!digit) {
+       AliFatal(Form("Digit not found at the expected position %d!",iDigit));
+      }
+      else {
+       amplList[iDigit] = (UShort_t)(digit->GetEnergy()*500); // Energy in units of GeV/500
+       digiList[iDigit] = (UShort_t)(digit->GetId());
+      }
+    }
+
+    ec->SetGlobalPosition(xyz);                 //rec.point position in MARS
+    ec->SetClusterEnergy(rp->Energy());         //total particle energy
+    ec->SetClusterDisp(emcRP->GetDispersion()); //cluster dispersion
+    ec->SetPid          (rp->GetPID()) ;        //array of particle identification
+    ec->SetM02(emcRP->GetM2x()) ;               //second moment M2x
+    ec->SetM20(emcRP->GetM2z()) ;               //second moment M2z
+    ec->SetNExMax(emcRP->GetNExMax());          //number of local maxima
+    ec->SetNumberOfDigits(digitMult);           //digit multiplicity
+    ec->SetDigitAmplitude(amplList);            //energies in 1/500 of GeV
+    ec->SetDigitIndex(digiList);                //abs id of the cell
+    ec->SetEmcCpvDistance(-1);                  //not yet implemented
+    ec->SetClusterChi2(-1);                     //not yet implemented
+    ec->SetM11(-1) ;                            //not yet implemented
+
+    // add the track to the esd object
+    esd->AddCaloCluster(ec);
+    delete ec;    
+
+    }
+  }
+
+
 }
 
 AliTracker* AliPHOSReconstructor::CreateTracker(AliRunLoader* runLoader) const
index f1bef35d1bc49bcdb1f48f557492014ecf969de4..2d5ac8497aac3301406308afad12d1d74d952e25 100644 (file)
@@ -8,6 +8,9 @@
 /* History of cvs commits:
  *
  * $Log$
+ * Revision 1.8  2005/05/28 14:19:04  schutz
+ * Compilation warnings fixed by T.P.
+ *
  */
 
 //_________________________________________________________________________
@@ -27,7 +30,7 @@ class AliPHOSTrackSegmentMaker ;
 class AliPHOSPID ;
 class AliPHOSSDigitizer ;
 class AliESD ;
-class AliRawReaderFile 
+class AliRawReader; 
 
 // --- Standard library ---
 
@@ -50,6 +53,7 @@ public:
   AliTracker *CreateTracker(AliRunLoader* runLoader) const;
   using AliReconstructor::FillESD;
   virtual void               FillESD(AliRunLoader* runLoader, AliESD* esd) const ;
+  virtual void FillESD(AliRunLoader* runLoader,AliRawReader* rawReader,AliESD* esd) const;
   using AliReconstructor::Reconstruct;
   virtual void               Reconstruct(AliRunLoader* runLoader) const ;
   virtual void               Reconstruct(AliRunLoader* runLoader, AliRawReader * rawreader) const ;