]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
calculation of fake rate and mc-reco-pt resolution
authorjgrosseo <jgrosseo@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 6 Sep 2012 07:19:41 +0000 (07:19 +0000)
committerjgrosseo <jgrosseo@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 6 Sep 2012 07:19:41 +0000 (07:19 +0000)
PWGCF/Correlations/Base/AliAnalyseLeadingTrackUE.cxx
PWGCF/Correlations/Base/AliAnalyseLeadingTrackUE.h
PWGCF/Correlations/Base/AliUEHist.cxx
PWGCF/Correlations/Base/AliUEHist.h
PWGCF/Correlations/Base/AliUEHistograms.cxx
PWGCF/Correlations/Base/AliUEHistograms.h
PWGCF/Correlations/DPhi/AliAnalysisTaskLeadingTrackUE.cxx
PWGCF/Correlations/DPhi/AliAnalysisTaskPhiCorrelations.cxx

index 34a8d3654dbd7ca7c32eb8f3ac6aed8664e28485..a8851ceac793cacefc0b4697a69730283d33e794 100644 (file)
@@ -327,6 +327,60 @@ TObjArray* AliAnalyseLeadingTrackUE::GetAcceptedParticles(TObject* obj, TObject*
   return tracks;
 }
 
+//-------------------------------------------------------------------
+TObjArray* AliAnalyseLeadingTrackUE::GetFakeParticles(TObject* obj, TObject* arrayMC, Bool_t onlyprimaries, Int_t particleSpecies, Bool_t useEtaPtCuts)
+{
+  // particleSpecies: -1 all particles are returned
+  //                  0 (pions) 1 (kaons) 2 (protons) 3 (others) particles
+
+  Int_t nTracks = NParticles(obj);
+  TObjArray* tracksReconstructed = new TObjArray;
+  TObjArray* tracksOriginal = new TObjArray;
+  TObjArray* tracksFake = new TObjArray;
+
+  // for TPC only tracks
+  Bool_t hasOwnership = kFALSE;
+  if ((fFilterBit == 128 || fFilterBit == 256 || fFilterBit == 512 || fFilterBit == 1024) && obj->InheritsFrom("AliESDEvent"))
+    hasOwnership = kTRUE;
+
+  tracksReconstructed->SetOwner(hasOwnership);
+  tracksFake->SetOwner(hasOwnership);
+
+  // Loop over tracks or jets
+  for (Int_t ipart=0; ipart<nTracks; ++ipart) {
+    AliVParticle* partReconstructed = ParticleWithCuts( obj, ipart, onlyprimaries, particleSpecies );
+    if (!partReconstructed) continue;
+
+    if (useEtaPtCuts)
+      if (TMath::Abs(partReconstructed->Eta()) > fTrackEtaCut || partReconstructed->Pt() < fTrackPtMin)
+      {
+        if (hasOwnership)
+          delete partReconstructed;
+        continue;
+      }
+
+    Int_t label = partReconstructed->GetLabel();
+    if (label == 0)
+    {
+      tracksFake->AddLast(partReconstructed);
+      continue;
+    }
+    if (hasOwnership)
+      delete partReconstructed;
+    AliVParticle* partOriginal = ParticleWithCuts(arrayMC, TMath::Abs(label),onlyprimaries, particleSpecies);
+    if (!partOriginal)continue;
+
+    tracksReconstructed->AddLast(partReconstructed);
+    tracksOriginal->AddLast(partOriginal);
+  }
+  TObjArray* pairs = new TObjArray;
+  pairs->SetOwner(kTRUE);
+  pairs->Add(tracksReconstructed);
+  pairs->Add(tracksOriginal);
+  pairs->Add(tracksFake);
+  return pairs;
+}
+
 //-------------------------------------------------------------------
 TObjArray* AliAnalyseLeadingTrackUE::GetMinMaxRegion(TList *transv1, TList *transv2)
 {
index 55de4ffbe568034ca02212b8e0ff184b7954c603..304cc50e9ad046168e37070ccdab261332abb212 100644 (file)
@@ -55,6 +55,7 @@ class AliAnalyseLeadingTrackUE : public TObject {
   void          QSortTracks(TObjArray &a, Int_t first, Int_t last);               // Sort by pT an array of AliVParticles 
   TObjArray*     SortRegions(const AliVParticle* leading, TObject* obj, TObject* arrayMC, Bool_t onlyprimaries = kTRUE); // Assign particles to towards, away or transverse regions
   TObjArray*     GetAcceptedParticles(TObject* obj, TObject* arrayMC, Bool_t onlyprimaries = kTRUE, Int_t particleSpecies = -1, Bool_t useEtaPtCuts = kFALSE); 
+  TObjArray*     GetFakeParticles(TObject* obj, TObject* arrayMC, Bool_t onlyprimaries, Int_t particleSpecies, Bool_t useEtaPtCuts);
   Bool_t         TriggerSelection(const TObject* obj);                                   // Select good triggers with AliPhysicsSelection class
   Bool_t         VertexSelection(const TObject* obj, Int_t ntracks, Double_t zed);       // Vertex selection: see implementation
   void                  RemoveInjectedSignals(TObjArray* tracks, TObject* arrayMC, Int_t maxLabel);
index d03cd5738460b6a0fd22723608369bdd4c3bf1fe..03363d646f85dd9e6734fbc20264c1383f17a4b7 100644 (file)
@@ -46,6 +46,7 @@ AliUEHist::AliUEHist(const char* reqHist) :
   fkRegions(4),
   fEventHist(0),
   fTrackHistEfficiency(0),
+  fFakePt(0),
   fEtaMin(0),
   fEtaMax(0),
   fPtMin(0),
@@ -299,7 +300,7 @@ AliUEHist::AliUEHist(const char* reqHist) :
   
   iTrackBin[2] = kNSpeciesBins;
 
-  fTrackHistEfficiency = new AliCFContainer("fTrackHistEfficiency", "Tracking efficiency", 3, 4, iTrackBin);
+  fTrackHistEfficiency = new AliCFContainer("fTrackHistEfficiency", "Tracking efficiency", 4, 4, iTrackBin);
   fTrackHistEfficiency->SetBinLimits(0, trackBins[0]);
   fTrackHistEfficiency->SetVarTitle(0, trackAxisTitle[0]);
   fTrackHistEfficiency->SetBinLimits(1, trackBins[1]);
@@ -308,6 +309,8 @@ AliUEHist::AliUEHist(const char* reqHist) :
   fTrackHistEfficiency->SetVarTitle(2, "particle species");
   fTrackHistEfficiency->SetBinLimits(3, trackBins[3]);
   fTrackHistEfficiency->SetVarTitle(3, trackAxisTitle[3]);
+
+  fFakePt = new TH3F("fFakePt","fFakePt;p_{T,rec};p_{T};centrality", 200, 0, 20, 200, 0, 20, 20, 0, 100);
 }
 
 //_____________________________________________________________________________
@@ -316,6 +319,7 @@ AliUEHist::AliUEHist(const AliUEHist &c) :
   fkRegions(4),
   fEventHist(0),
   fTrackHistEfficiency(0),
+  fFakePt(0),
   fEtaMin(0),
   fEtaMax(0),
   fPtMin(0),
@@ -373,6 +377,12 @@ AliUEHist::~AliUEHist()
     fTrackHistEfficiency = 0;
   }
 
+  if (fFakePt)
+  {
+    delete fFakePt;
+    fFakePt = 0;
+  }
+
   if (fCache)
   {
     delete fCache;
@@ -408,6 +418,9 @@ void AliUEHist::Copy(TObject& c) const
   if (fTrackHistEfficiency)
     target.fTrackHistEfficiency = dynamic_cast<AliCFContainer*> (fTrackHistEfficiency->Clone());
     
+  if (fFakePt)
+    target.fFakePt = dynamic_cast<TH3F*> (fFakePt->Clone());
+
   target.fEtaMin = fEtaMin;
   target.fEtaMax = fEtaMax;
   target.fPtMin = fPtMin;
@@ -441,7 +454,7 @@ Long64_t AliUEHist::Merge(TCollection* list)
   TObject* obj;
 
   // collections of objects
-  const UInt_t kMaxLists = fkRegions+2;
+  const UInt_t kMaxLists = fkRegions+3;
   TList** lists = new TList*[kMaxLists];
   
   for (UInt_t i=0; i<kMaxLists; i++)
@@ -460,6 +473,8 @@ Long64_t AliUEHist::Merge(TCollection* list)
     
     lists[fkRegions]->Add(entry->fEventHist);
     lists[fkRegions+1]->Add(entry->fTrackHistEfficiency);
+    if (entry->fFakePt)
+      lists[fkRegions+2]->Add(entry->fFakePt);
 
     count++;
   }
@@ -469,6 +484,8 @@ Long64_t AliUEHist::Merge(TCollection* list)
   
   fEventHist->Merge(lists[fkRegions]);
   fTrackHistEfficiency->Merge(lists[fkRegions+1]);
+  if (fFakePt)
+    fFakePt->Merge(lists[fkRegions+2]);
 
   for (UInt_t i=0; i<kMaxLists; i++)
     delete lists[i];
@@ -850,6 +867,7 @@ TH2* AliUEHist::GetSumOfRatios2(AliUEHist* mixed, AliUEHist::CFStep step, AliUEH
   TH2* eventMixedAll = 0;
   
   Int_t totalEvents = 0;
+  Int_t nCorrelationFunctions = 0;
   
   GetHistsZVtxMult(step, region, ptLeadMin, ptLeadMax, &trackSameAll, &eventSameAll);
   mixed->GetHistsZVtxMult(step, region, ptLeadMin, ptLeadMax, &trackMixedAll, &eventMixedAll);
@@ -959,11 +977,13 @@ TH2* AliUEHist::GetSumOfRatios2(AliUEHist* mixed, AliUEHist::CFStep step, AliUEH
 
       delete tracksSame;
       delete tracksMixed;
+      
+      nCorrelationFunctions++;
     }
   }
 
   if (totalTracks) {
-    Printf("Dividing %f tracks by %d events", totalTracks->Integral(), totalEvents);
+    Printf("Dividing %f tracks by %d events (%d correlation function(s))", totalTracks->Integral(), totalEvents, nCorrelationFunctions);
     if (totalEvents > 0)
       totalTracks->Scale(1.0 / totalEvents);
   
@@ -2141,6 +2161,12 @@ TH2D* AliUEHist::GetTrackingEfficiency()
   return dynamic_cast<TH2D*> (GetTrackEfficiency(kCFStepAnaTopology, kCFStepTrackedOnlyPrim, 0, 1));
 }
   
+//____________________________________________________________________
+TH2D* AliUEHist::GetFakeRate()
+{
+  return dynamic_cast<TH2D*> (GetTrackEfficiency(kCFStepTracked, kCFStepReconstructed, 0, 1));
+}
+
 //____________________________________________________________________
 TH2D* AliUEHist::GetTrackingEfficiencyCentrality()
 {
index 3ebf2ab359e5d4cf065d10c0a1f9261afc983f09..98c7d8898fb02a11464d661e242acfc0c4676fb5 100644 (file)
@@ -15,6 +15,7 @@ class AliCFContainer;
 class TH1;
 class TH1F;
 class TH3;
+class TH3F;
 class TH1D;
 class TH2;
 class TH2D;
@@ -41,7 +42,8 @@ class AliUEHist : public TObject
   AliCFContainer* GetTrackHist(Region region) { return fTrackHist[region]; }
   AliCFContainer* GetEventHist() { return fEventHist; }
   AliCFContainer* GetTrackHistEfficiency()     { return fTrackHistEfficiency; }
-  
+  TH3F* GetMCRecoPtCorrelation() { return fFakePt; } 
   void SetTrackHist(Region region, AliCFContainer* hist) { fTrackHist[region] = hist; }
   void SetEventHist(AliCFContainer* hist) { fEventHist = hist; }
   void SetTrackHistEfficiency(AliCFContainer* hist) { fTrackHistEfficiency = hist; }
@@ -66,6 +68,8 @@ class AliUEHist : public TObject
   TH2D* GetTrackingEfficiency();
   TH2D* GetTrackingEfficiencyCentrality();
   
+  TH2D* GetFakeRate();
+
   TH1D* GetTrackingContamination(Int_t axis);
   TH2D* GetTrackingContamination();
   TH2D* GetTrackingContaminationCentrality();
@@ -130,7 +134,8 @@ protected:
   AliCFContainer* fTrackHist[4];      // container for track level distributions in four regions (toward, away, min, max) and at four analysis steps
   AliCFContainer* fEventHist;         // container for event level distribution at four analysis steps
   AliCFContainer* fTrackHistEfficiency; // container for tracking efficiency and contamination (all particles filled including leading one): axes: eta, pT, particle species
-  
+  TH3F* fFakePt;
   Float_t fEtaMin;                    // eta min for projections
   Float_t fEtaMax;                    // eta max for projections
   Float_t fPtMin;                     // pT min for projections (for track pT, not pT,lead)
index a2eb16eaf1173daa1311995cffb554888155567a..7162523d89e8143e5298d1f4bf3e9eb092393f2c 100644 (file)
@@ -662,7 +662,7 @@ void AliUEHistograms::FillCorrelations(Double_t centrality, Float_t zVtx, AliUEH
 }
   
 //____________________________________________________________________
-void AliUEHistograms::FillTrackingEfficiency(TObjArray* mc, TObjArray* recoPrim, TObjArray* recoAll, Int_t particleType, Double_t centrality)
+void AliUEHistograms::FillTrackingEfficiency(TObjArray* mc, TObjArray* recoPrim, TObjArray* recoAll, TObjArray* fake, Int_t particleType, Double_t centrality)
 {
   // fills the tracking efficiency objects
   //
@@ -670,15 +670,20 @@ void AliUEHistograms::FillTrackingEfficiency(TObjArray* mc, TObjArray* recoPrim,
   // recoPrim: reconstructed primaries (again MC particles)
   // recoAll: reconstructed (again MC particles)
   // particleType is: 0 for pion, 1 for kaon, 2 for proton, 3 for others
-  
-  for (Int_t step=0; step<3; step++)
+  for (Int_t step=0; step<4; step++)
   {
     TObjArray* list = mc;
     if (step == 1)
       list = recoPrim;
     else if (step == 2)
       list = recoAll;
-      
+    else if (step == 3)
+      list = fake;
+    
+    if (!list)
+      continue;
+
     for (Int_t i=0; i<list->GetEntriesFast(); i++)
     {
       AliVParticle* particle = (AliVParticle*) list->At(i);
@@ -695,6 +700,26 @@ void AliUEHistograms::FillTrackingEfficiency(TObjArray* mc, TObjArray* recoPrim,
   }
 }
 
+//____________________________________________________________________
+void AliUEHistograms::FillFakePt(TObjArray* fake, Double_t centrality)
+{
+  TObjArray* tracksReco = (TObjArray*) fake->At(0);
+  TObjArray* tracksMC = (TObjArray*) fake->At(1);
+  
+  for (Int_t i=0; i<tracksReco->GetEntriesFast(); i++)
+  {
+    AliVParticle* particle1 = (AliVParticle*) tracksReco->At(i);
+    AliVParticle* particle2 = (AliVParticle*) tracksMC->At(i);
+    Double_t vars[3];
+    vars[0] = particle1->Pt();
+    vars[1] = particle2->Pt();
+    vars[2] = centrality;
+    for (Int_t j=0; j<fgkUEHists; j++)
+      if (GetUEHist(j))
+        GetUEHist(j)->GetMCRecoPtCorrelation()->Fill(vars[0],vars[1],vars[2]);
+  }
+}
+
 //____________________________________________________________________
 void AliUEHistograms::FillEvent(Int_t eventType, Int_t step)
 {
index a7ee286471c82eb6f3ec92ad2b5895eeddec11d0..bfda7256c1897e3f167bba1fea0ea81c4523f66d 100644 (file)
@@ -32,8 +32,9 @@ class AliUEHistograms : public TNamed
   void Fill(AliVParticle* leadingMC, AliVParticle* leadingReco);
   void FillEvent(Int_t eventType, Int_t step);
   void FillEvent(Double_t centrality, Int_t step);
-  void FillTrackingEfficiency(TObjArray* mc, TObjArray* recoPrim, TObjArray* recoAll, Int_t particleType, Double_t centrality = 0);
-  
+  void FillTrackingEfficiency(TObjArray* mc, TObjArray* recoPrim, TObjArray* recoAll, TObjArray* fake, Int_t particleType, Double_t centrality = 0);
+  void FillFakePt(TObjArray* fake, Double_t centrality);
   void CopyReconstructedData(AliUEHistograms* from);
   void DeepCopy(AliUEHistograms* from);
   
index 565abde89be4d5b543136145310759c8ddc07802..4762ea4ce867fa667f507d7522af7f78a7a29a2b 100644 (file)
@@ -400,8 +400,8 @@ void  AliAnalysisTaskLeadingTrackUE::AnalyseCorrectionMode()
                TObjArray* primMCParticles = fAnalyseUE->GetAcceptedParticles(fArrayMC, 0x0, kTRUE, particleSpecies);
                TObjArray* primRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(fAOD, fArrayMC, kTRUE, particleSpecies);
                TObjArray* allRecoTracksMatched = fAnalyseUE->GetAcceptedParticles(fAOD, fArrayMC, kFALSE, particleSpecies);
-              
-               fHistosUE->FillTrackingEfficiency(primMCParticles, primRecoTracksMatched, allRecoTracksMatched, particleSpecies);
+
+               fHistosUE->FillTrackingEfficiency(primMCParticles, primRecoTracksMatched, allRecoTracksMatched, 0, particleSpecies);
               
                delete primMCParticles;
                delete primRecoTracksMatched;
index b508b396bb4106b65f949010bf7b1cd08a3efd87..c5bdd0d8610c6855419950ea0c1ea41a2495ffa6 100644 (file)
@@ -519,7 +519,7 @@ void  AliAnalysisTaskPhiCorrelations::AnalyseCorrectionMode()
          fAnalyseUE->RemoveInjectedSignals(allRecoTracksMatched, mc, skipParticlesAbove);
        }
       
-        fHistos->FillTrackingEfficiency(primMCParticles, primRecoTracksMatched, allRecoTracksMatched, particleSpecies, centrality);
+        fHistos->FillTrackingEfficiency(primMCParticles, primRecoTracksMatched, allRecoTracksMatched, 0, particleSpecies, centrality);
         
 //     Printf("%d --> %d %d %d", particleSpecies, primMCParticles->GetEntries(), primRecoTracksMatched->GetEntries(), allRecoTracksMatched->GetEntries());
 
@@ -527,6 +527,10 @@ void  AliAnalysisTaskPhiCorrelations::AnalyseCorrectionMode()
         delete primRecoTracksMatched;
         delete allRecoTracksMatched;
       }
+      TObjArray* fakeParticles = fAnalyseUE->GetFakeParticles(inputEvent, mc, kFALSE, -1, kTRUE);
+      fHistos->FillTrackingEfficiency(0, 0, 0, (TObjArray*) fakeParticles->At(2), -1, centrality);
+      fHistos->FillFakePt(fakeParticles, centrality);
+      delete fakeParticles;
     
       // (MC-true all particles)
       // STEP 2