]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
added corrections for delta phi UE event distribution
authorjgrosseo <jgrosseo@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 21 Oct 2010 08:14:01 +0000 (08:14 +0000)
committerjgrosseo <jgrosseo@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 21 Oct 2010 08:14:01 +0000 (08:14 +0000)
tracking efficiency is stored in settings tree for each worker

PWG4/JetTasks/AliAnalysisTaskLeadingTrackUE.cxx
PWG4/JetTasks/AliUEHist.cxx
PWG4/JetTasks/AliUEHist.h
PWG4/JetTasks/AliUEHistograms.cxx
PWG4/JetTasks/AliUEHistograms.h
PWG4/JetTasks/inputFiles/ue_trackingefficiency.root

index 864c737b471496c3a7312801dfcadc000a2a49d2..d44f7d599bc9a69cb8120030d39002401590a69c 100644 (file)
@@ -295,6 +295,7 @@ void  AliAnalysisTaskLeadingTrackUE::AddSettingsTree()
   settingsTree->Branch("fSelectBit", &fSelectBit,"EventSelectionBit/I");
   settingsTree->Branch("fLeadingTrackEtaCut", &fLeadingTrackEtaCut, "LeadingTrackEtaCut/D");
   settingsTree->Branch("fUseChargeHadrons", &fUseChargeHadrons,"UseChHadrons/O");
+  settingsTree->Branch("fkTrackingEfficiency", "TH1D", &fkTrackingEfficiency);
   settingsTree->Fill();
   fListOfHistos->Add(settingsTree);
 }  
index 4480b0337f978a85d0fe37837532f5a03d6317c8..3776d73f61627250c7a69042ef6ff15e6bf57276 100644 (file)
@@ -89,10 +89,8 @@ AliUEHist::AliUEHist(const char* reqHist) :
     leadingpTBins[i] = 0.5 * i;
   
   // pT,lead binning 2
-  const Int_t kNLeadingpTBins2 = 10;
-  Double_t leadingpTBins2[kNLeadingpTBins2+1];
-  for (Int_t i=0; i<=kNLeadingpTBins2; i++)
-    leadingpTBins2[i] = 5.0 * i;
+  const Int_t kNLeadingpTBins2 = 13;
+  Double_t leadingpTBins2[] = { 0.0, 0.5, 1.0, 2.0, 4.0, 6.0, 8.0, 10.0, 15.0, 20.0, 30.0, 40.0, 50.0, 100.0 };
   
   // phi,lead
   const Int_t kNLeadingPhiBins = 40;
@@ -484,6 +482,7 @@ TH1D* AliUEHist::GetUEHist(AliUEHist::CFStep step, AliUEHist::Region region, Flo
   {
     fTrackHist[region]->GetGrid(step)->SetRangeUser(2, ptLeadMin, ptLeadMax);
     tracks = fTrackHist[region]->GetGrid(step)->Project(4);
+    Printf("Calculated histogram in %.2f <= pT <= %.2f --> %f entries", ptLeadMin, ptLeadMax, tracks->Integral());
     fTrackHist[region]->GetGrid(step)->SetRangeUser(2, 0, -1);
     
     // normalize to get a density (deta dphi)
@@ -625,13 +624,26 @@ void AliUEHist::Correct(AliUEHist* corrections)
   {
     if (!fTrackHist[region])
       continue;
-      
-    //TH1D* leadingBias = (TH1D*) corrections->GetBias(kCFStepReconstructed, kCFStepTracked, region, "z"); // from MC
-    TH1D* leadingBias = (TH1D*) GetBias(kCFStepBiasStudy, kCFStepReconstructed, region, "z");          // from data
-    CorrectTracks(kCFStepReconstructed, kCFStepTracked, region, leadingBias, 2);
+
+    const char* projAxis = "z";
+    Int_t secondBin = -1;
+
+    if (fTrackHist[region]->GetNVar() == 5)
+    {
+      projAxis = "yz";
+      secondBin = 4;
+    }
+    
+    #if 0
+      TH1* leadingBias = (TH1*) corrections->GetBias(kCFStepReconstructed, kCFStepTracked, region, projAxis); // from MC
+      Printf("WARNING: Using MC bias correction");
+    #else
+      TH1* leadingBias = (TH1*) GetBias(kCFStepBiasStudy, kCFStepReconstructed, region, projAxis);          // from data
+    #endif
+    CorrectTracks(kCFStepReconstructed, kCFStepTracked, region, leadingBias, 2, secondBin);
     if (region == kMin && fCombineMinMax)
     {
-      CorrectTracks(kCFStepReconstructed, kCFStepTracked, kMax, leadingBias, 2);
+      CorrectTracks(kCFStepReconstructed, kCFStepTracked, kMax, leadingBias, 2, secondBin);
       delete leadingBias;
       break;
     }
@@ -672,7 +684,7 @@ void AliUEHist::Correct(AliUEHist* corrections)
   vertexCorrectionObs->Reset();
   
   TF1* func = new TF1("func", "[1]+[0]/(x-[2])");
-  vertexCorrection->Fit(func, "0", "", 0.8, 4);
+  vertexCorrection->Fit(func, "0", "", 0, 3);
 
   for (Int_t i=1; i<=vertexCorrectionObs->GetNbinsX(); i++)
   {
@@ -688,12 +700,15 @@ void AliUEHist::Correct(AliUEHist* corrections)
   vertexCorrection->DrawCopy();
   vertexCorrectionObs->SetLineColor(2);
   vertexCorrectionObs->DrawCopy("same");
+  func->SetRange(0, 4);
+  func->DrawClone("same");
   #endif
   
   CorrectTracks(kCFStepVertex, kCFStepTriggered, vertexCorrectionObs, 3);
   CorrectEvents(kCFStepVertex, kCFStepTriggered, vertexCorrectionObs, 1);
   delete vertexCorrectionObs;
   delete vertexCorrection;
+  delete func;
   
   // copy 
   CorrectTracks(kCFStepTriggered, kCFStepAll, 0, -1);
@@ -955,11 +970,11 @@ TH1* AliUEHist::GetBias(CFStep step1, CFStep step2, Int_t region, const char* ax
   SetBinLimits(tmp->GetGrid(step2));
   
   TH1D* events1 = fEventHist->Project(0, step1);
-  TH3D* hist1 = tmp->Project(0, 1, 2, step1);
+  TH3D* hist1 = tmp->Project(0, tmp->GetNVar()-1, 2, step1);
   WeightHistogram(hist1, events1);
   
   TH1D* events2 = fEventHist->Project(0, step2);
-  TH3D* hist2 = tmp->Project(0, 1, 2, step2);
+  TH3D* hist2 = tmp->Project(0, tmp->GetNVar()-1, 2, step2);
   WeightHistogram(hist2, events2);
   
   TH1* generated = hist1;
@@ -978,12 +993,7 @@ TH1* AliUEHist::GetBias(CFStep step1, CFStep step2, Int_t region, const char* ax
       hist1->GetXaxis()->SetRangeUser(fEtaMin, fEtaMax);
       hist2->GetXaxis()->SetRangeUser(fEtaMin, fEtaMax);
     }
-    if (fPtMax > fPtMin && !TString(axis).Contains("y"))
-    {
-      hist1->GetYaxis()->SetRangeUser(fPtMin, fPtMax);
-      hist2->GetYaxis()->SetRangeUser(fPtMin, fPtMax);
-    }
-    
+   
     generated = hist1->Project3D(axis);
     measured  = hist2->Project3D(axis);
     
index f853fa5728a9f1e98545ff18971feb93d5049c47..879766ef39aa865f8400a28841342debac2b3e32 100644 (file)
@@ -40,6 +40,7 @@ class AliUEHist : public TObject
   
   void SetTrackHist(Region region, AliCFContainer* hist) { fTrackHist[region] = hist; }
   void SetEventHist(AliCFContainer* hist) { fEventHist = hist; }
+  void SetTrackHistEfficiency(AliCFContainer* hist) { fTrackHistEfficiency = hist; }
   
   void CopyReconstructedData(AliUEHist* from);
   
@@ -47,7 +48,7 @@ class AliUEHist : public TObject
   
   TH1* GetTrackEfficiency(CFStep step1, CFStep step2, Int_t axis1, Int_t axis2 = -1, Int_t source = 1);
   TH1* GetEventEfficiency(CFStep step1, CFStep step2, Int_t axis1, Int_t axis2 = -1, Float_t ptLeadMin = -1, Float_t ptLeadMax = -1);
-  TH1* GetBias(CFStep step1, CFStep step2, Int_t region, const char* axis, Float_t leadPtMin = 0, Float_t leadPtMax = 0);
+  TH1* GetBias(CFStep step1, CFStep step2, Int_t region, const char* axis, Float_t leadPtMin = 0, Float_t leadPtMax = -1);
   
   TH1D* GetTrackingEfficiency(Int_t axis);
   TH2D* GetTrackingEfficiency();
@@ -75,6 +76,9 @@ class AliUEHist : public TObject
   
   void CountEmptyBins(AliUEHist::CFStep step, Float_t ptLeadMin, Float_t ptLeadMax);
   
+  void SetBinLimits(AliCFGridSparse* grid);
+  void ResetBinLimits(AliCFGridSparse* grid);
+  
   AliUEHist(const AliUEHist &c);
   AliUEHist& operator=(const AliUEHist& corr);
   virtual void Copy(TObject& c) const;
@@ -84,8 +88,6 @@ class AliUEHist : public TObject
 protected:
   void SetStepNames(AliCFContainer* container);
   void WeightHistogram(TH3* hist1, TH1* hist2);
-  void SetBinLimits(AliCFGridSparse* grid);
-  void ResetBinLimits(AliCFGridSparse* grid);
 
   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
index 79c3e8b0da09658109712b3410683e3787d7c684..3c4347b1a058c3dd8593515ce3468673046a4bbe 100644 (file)
@@ -183,6 +183,20 @@ AliUEHistograms::~AliUEHistograms()
   }
 }
 
+AliUEHist* AliUEHistograms::GetUEHist(Int_t id)
+{
+  // returns AliUEHist object, useful for loops
+  
+  switch (id)
+  {
+    case 0: return fNumberDensitypT; break;
+    case 1: return fSumpT; break;
+    case 2: return fNumberDensityPhi; break;
+  }
+  
+  return 0;
+}
+
 //____________________________________________________________________
 Int_t AliUEHistograms::CountParticles(TList* list, Float_t ptMin)
 {
@@ -391,7 +405,7 @@ void AliUEHistograms::Correct(AliUEHistograms* corrections)
   
   fNumberDensitypT->Correct(corrections->fNumberDensitypT);
   fSumpT->Correct(corrections->fSumpT);
-  //fNumberDensityPhi->Correct(corrections->fNumberDensityPhi);
+  fNumberDensityPhi->Correct(corrections->fNumberDensityPhi);
 }
 
 //____________________________________________________________________
index 1877e2c7804ea47881e70a87c05f9b04674d2879..14a3abf6d9be22426fffc82452789eb76e080582 100644 (file)
@@ -31,6 +31,8 @@ class AliUEHistograms : public TObject
   
   void CopyReconstructedData(AliUEHistograms* from);
   
+  AliUEHist* GetUEHist(Int_t id);
+  
   AliUEHist* GetNumberDensitypT() { return fNumberDensitypT; }
   AliUEHist* GetSumpT() { return fSumpT; }
   AliUEHist* GetNumberDensityPhi() { return fNumberDensityPhi; }
index 2c76b9ccf1f793680b667591552775ca373e9a5c..26bad5ef166ca509be7edf4a63133ea6fe6a7e78 100644 (file)
Binary files a/PWG4/JetTasks/inputFiles/ue_trackingefficiency.root and b/PWG4/JetTasks/inputFiles/ue_trackingefficiency.root differ