]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Cleaning the ESD: first step
authorbelikov <belikov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 1 Aug 2007 15:05:46 +0000 (15:05 +0000)
committerbelikov <belikov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 1 Aug 2007 15:05:46 +0000 (15:05 +0000)
STEER/AliESD.cxx
STEER/AliESD.h
STEER/AliESDEvent.cxx
STEER/AliESDEvent.h
STEER/AliReconstruction.cxx
STEER/AliReconstruction.h

index b54ad043afd401bdd94904bf7b0320b62de4422b..576e276954a8199ca1601d079f81c52219ff91e1 100644 (file)
@@ -213,6 +213,25 @@ void AliESD::Reset()
   fErrorLogs.Clear();
 }
 
+Bool_t  AliESD::RemoveTrack(Int_t /*i*/) {
+  // ---------------------------------------------------------
+  // Remove track
+  // ---------------------------------------------------------
+
+  // Check if this track comes from a reconstructed decay
+  // if (yes) return kFALSE
+
+  // Remap the indices of the daughters of recosntructed decays
+
+  // Remove the track
+  // delete fTracks->RemoveAt(i);
+
+  // Compress the array with tracks
+  // fTracks->Compress();
+
+  return kTRUE;
+}
+
 Int_t AliESD::AddV0(const AliESDv0 *v) {
   //
   // Add V0
index 13026e518a0e2cac6fbe9047ed6af9cbc98b46e6..1daea4764abb04f8b96b5facd8920ffea6c709e3 100644 (file)
@@ -75,6 +75,8 @@ public:
     return (AliESDTrdTrack *)fTrdTracks.UncheckedAt(i);
   }
 
+  Bool_t  RemoveTrack(Int_t i);
+
   Int_t  AddTrack(const AliESDtrack *t) {
     AliESDtrack * track = new(fTracks[fTracks.GetEntriesFast()]) AliESDtrack(*t);
     track->SetID(fTracks.GetEntriesFast()-1);
index cd4772aa38c3e66a6581377e65b0c3bca2db0db7..dbb22bd7b55ecb3a158e8051a891c1a16567112b 100644 (file)
@@ -372,6 +372,26 @@ void AliESDEvent::SetESDfriend(const AliESDfriend *ev) {
   }
 }
 
+Bool_t  AliESDEvent::RemoveTrack(Int_t /*i*/) {
+  // ---------------------------------------------------------
+  // Remove track
+  // ---------------------------------------------------------
+
+  // Check if this track comes from a reconstructed decay
+  // if (yes) return kFALSE
+
+  // Remap the indices of the daughters of recosntructed decays
+
+  // Remove the track
+  // delete fTracks->RemoveAt(i);
+
+  // Compress the array with tracks
+  // fTracks->Compress();
+
+  return kTRUE;
+}
+
+
 Int_t  AliESDEvent::AddTrack(const AliESDtrack *t) {
     // Add track
     TClonesArray &ftr = *fTracks;
index fe53fcf23f46b7d7386a91a7bbca56e394447470..c24398187fc4c404b89d0210f1b9a9ab844a4b1a 100644 (file)
@@ -173,7 +173,8 @@ public:
   const AliMultiplicity *GetMultiplicity() const {return fSPDMult;}
 
 
-  
+  Bool_t  RemoveTrack(Int_t i);
+
   AliESDtrack *GetTrack(Int_t i) const {
     return (AliESDtrack *)fTracks->UncheckedAt(i);
   }
index 60d045c1f38fea91451657118457d0752a799966..4938b0aebfdb0d6fd56442acc4203ff7ec072475 100644 (file)
@@ -187,6 +187,7 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename, const char* cdb
   fRunMuonTracking(kFALSE),
   fStopOnError(kFALSE),
   fWriteAlignmentData(kFALSE),
+  fCleanESD(kTRUE),
   fWriteESDfriend(kFALSE),
   fWriteAOD(kFALSE),
   fFillTriggerESD(kTRUE),
@@ -236,6 +237,7 @@ AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
   fRunMuonTracking(rec.fRunMuonTracking),
   fStopOnError(rec.fStopOnError),
   fWriteAlignmentData(rec.fWriteAlignmentData),
+  fCleanESD(rec.fCleanESD),
   fWriteESDfriend(rec.fWriteESDfriend),
   fWriteAOD(rec.fWriteAOD),
   fFillTriggerESD(rec.fFillTriggerESD),
@@ -773,6 +775,7 @@ Bool_t AliReconstruction::Run(const char* input)
     }
  
     // write ESD
+    if (fCleanESD) CleanESD(esd);
     if (fWriteESDfriend) {
       new (esdf) AliESDfriend(); // Reset...
       esd->GetESDfriend(esdf);
@@ -1245,6 +1248,32 @@ Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
   return kTRUE;
 }
 
+//_____________________________________________________________________________
+Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
+  //
+  // Remove the data which are not needed for the physics analysis.
+  //
+
+  AliInfo("Cleaning the ESD...");
+
+  const AliESDVertex *vertex=esd->GetVertex();
+  Double_t vz=vertex->GetZv();
+  
+  Int_t nTracks=esd->GetNumberOfTracks();
+  for (Int_t i=0; i<nTracks; i++) {
+    AliESDtrack *track=esd->GetTrack(i);
+
+    Float_t xy,z; track->GetImpactParameters(xy,z);
+    if (TMath::Abs(xy) < 50.)    continue;  
+    if (vertex->GetStatus())
+      if (TMath::Abs(vz-z) < 5.) continue;  
+
+    esd->RemoveTrack(i);
+  }
+
+  return kTRUE;
+}
+
 //_____________________________________________________________________________
 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
 {
index 615534dd4a42d91bb01d21e99994e171b1d62a20..1af41fe507c17c26fbb00200f6f172c8f9983474 100644 (file)
@@ -74,6 +74,7 @@ public:
   void SetRunHLTTracking(Bool_t flag=kTRUE) {fRunHLTTracking=flag;};
   void SetStopOnError(Bool_t flag=kTRUE) {fStopOnError=flag;}
   void SetWriteAlignmentData(Bool_t flag=kTRUE){fWriteAlignmentData=flag;}
+  void SetCleanESD(Bool_t flag=kTRUE){fCleanESD=flag;}
   void SetWriteESDfriend(Bool_t flag=kTRUE){fWriteESDfriend=flag;}
   void SetWriteAOD(Bool_t flag=kTRUE){fWriteAOD=flag;}
   void SetFillTriggerESD(Bool_t flag=kTRUE){fFillTriggerESD=flag;}
@@ -104,6 +105,7 @@ private:
   Bool_t         RunHLTTracking(AliESDEvent*& esd);
   Bool_t         RunMuonTracking(AliESDEvent*& esd);
   Bool_t         RunTracking(AliESDEvent*& esd);
+  Bool_t         CleanESD(AliESDEvent *esd);
   Bool_t         FillESD(AliESDEvent*& esd, const TString& detectors);
   Bool_t         FillTriggerESD(AliESDEvent*& esd);
   Bool_t         FillRawEventHeaderESD(AliESDEvent*& esd);
@@ -132,6 +134,7 @@ private:
   Bool_t         fRunMuonTracking;     // run the HLT tracking
   Bool_t         fStopOnError;        // stop or continue on errors
   Bool_t         fWriteAlignmentData; // write track space-points flag
+  Bool_t         fCleanESD;           // clean ESD flag
   Bool_t         fWriteESDfriend;     // write ESD friend flag
   Bool_t         fWriteAOD;           // write AOD flag
   Bool_t         fFillTriggerESD;     // fill trigger info into ESD
@@ -168,7 +171,7 @@ private:
   TString       fCDBUri;             // Uri of the default CDB storage
   TObjArray      fSpecCDBUri;         // Array with detector specific CDB storages
 
-  ClassDef(AliReconstruction, 11)      // class for running the reconstruction
+  ClassDef(AliReconstruction, 12)      // class for running the reconstruction
 };
 
 #endif