]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
add single inclusive TT + AllMatch histos
authormverweij <marta.verweij@cern.ch>
Tue, 2 Dec 2014 10:43:58 +0000 (11:43 +0100)
committermverweij <marta.verweij@cern.ch>
Tue, 2 Dec 2014 11:47:57 +0000 (12:47 +0100)
PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskEmcalHJetMass.cxx
PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskEmcalHJetMass.h

index 66c9af8b20fddddaa70c1f23ae0d79995b4736a8..c92384ad248a8babdb1bf30b6af527b2d581ca8f 100644 (file)
@@ -53,9 +53,11 @@ namespace EmcalHJetMassAnalysis {
     fh1PtHadron(0),
     fh3PtHPtJDPhi(0),
     fh3PtJet1VsMassVsHPtAllSel(0),
+    fh3PtJet1VsMassVsHPtAllSelMatch(0),
     fh3PtJet1VsMassVsHPtTagged(0),
     fh3PtJet1VsMassVsHPtTaggedMatch(0),
     fh3PtJet1VsRatVsHPtAllSel(0),
+    fh3PtJet1VsRatVsHPtAllSelMatch(0),
     fh3PtJet1VsRatVsHPtTagged(0),
     fh3PtJet1VsRatVsHPtTaggedMatch(0)
   {
@@ -64,9 +66,11 @@ namespace EmcalHJetMassAnalysis {
     fh1PtHadron                       = new TH1F*[fNcentBins];
     fh3PtHPtJDPhi                     = new TH3F*[fNcentBins];
     fh3PtJet1VsMassVsHPtAllSel        = new TH3F*[fNcentBins];
+    fh3PtJet1VsMassVsHPtAllSelMatch   = new TH3F*[fNcentBins];
     fh3PtJet1VsMassVsHPtTagged        = new TH3F*[fNcentBins];
     fh3PtJet1VsMassVsHPtTaggedMatch   = new TH3F*[fNcentBins];
     fh3PtJet1VsRatVsHPtAllSel         = new TH3F*[fNcentBins];
+    fh3PtJet1VsRatVsHPtAllSelMatch    = new TH3F*[fNcentBins];
     fh3PtJet1VsRatVsHPtTagged         = new TH3F*[fNcentBins];
     fh3PtJet1VsRatVsHPtTaggedMatch    = new TH3F*[fNcentBins];
 
@@ -74,9 +78,11 @@ namespace EmcalHJetMassAnalysis {
       fh1PtHadron[i]                       = 0;
       fh3PtHPtJDPhi[i]                     = 0;
       fh3PtJet1VsMassVsHPtAllSel[i]        = 0;
+      fh3PtJet1VsMassVsHPtAllSelMatch[i]   = 0;
       fh3PtJet1VsMassVsHPtTagged[i]        = 0;
       fh3PtJet1VsMassVsHPtTaggedMatch[i]   = 0;
       fh3PtJet1VsRatVsHPtAllSel[i]         = 0;
+      fh3PtJet1VsRatVsHPtAllSelMatch[i]    = 0;
       fh3PtJet1VsRatVsHPtTagged[i]         = 0;
       fh3PtJet1VsRatVsHPtTaggedMatch[i]    = 0;
     }
@@ -103,9 +109,11 @@ namespace EmcalHJetMassAnalysis {
     fh1PtHadron(0),
     fh3PtHPtJDPhi(0),
     fh3PtJet1VsMassVsHPtAllSel(0),
+    fh3PtJet1VsMassVsHPtAllSelMatch(0),
     fh3PtJet1VsMassVsHPtTagged(0),
     fh3PtJet1VsMassVsHPtTaggedMatch(0),
     fh3PtJet1VsRatVsHPtAllSel(0),
+    fh3PtJet1VsRatVsHPtAllSelMatch(0),
     fh3PtJet1VsRatVsHPtTagged(0),
     fh3PtJet1VsRatVsHPtTaggedMatch(0)
   {
@@ -114,9 +122,11 @@ namespace EmcalHJetMassAnalysis {
     fh1PtHadron                       = new TH1F*[fNcentBins];
     fh3PtHPtJDPhi                     = new TH3F*[fNcentBins];
     fh3PtJet1VsMassVsHPtAllSel        = new TH3F*[fNcentBins];
+    fh3PtJet1VsMassVsHPtAllSelMatch   = new TH3F*[fNcentBins];
     fh3PtJet1VsMassVsHPtTagged        = new TH3F*[fNcentBins];
     fh3PtJet1VsMassVsHPtTaggedMatch   = new TH3F*[fNcentBins];
     fh3PtJet1VsRatVsHPtAllSel         = new TH3F*[fNcentBins];
+    fh3PtJet1VsRatVsHPtAllSelMatch    = new TH3F*[fNcentBins];
     fh3PtJet1VsRatVsHPtTagged         = new TH3F*[fNcentBins];
     fh3PtJet1VsRatVsHPtTaggedMatch    = new TH3F*[fNcentBins];
  
@@ -124,9 +134,11 @@ namespace EmcalHJetMassAnalysis {
       fh1PtHadron[i]                       = 0;
       fh3PtHPtJDPhi[i]                     = 0;
       fh3PtJet1VsMassVsHPtAllSel[i]        = 0;
+      fh3PtJet1VsMassVsHPtAllSelMatch[i]   = 0;
       fh3PtJet1VsMassVsHPtTagged[i]        = 0;
       fh3PtJet1VsMassVsHPtTaggedMatch[i]   = 0;
       fh3PtJet1VsRatVsHPtAllSel[i]         = 0;
+      fh3PtJet1VsRatVsHPtAllSelMatch[i]    = 0;
       fh3PtJet1VsRatVsHPtTagged[i]         = 0;
       fh3PtJet1VsRatVsHPtTaggedMatch[i]    = 0;
     }
@@ -195,6 +207,11 @@ namespace EmcalHJetMassAnalysis {
       fh3PtJet1VsMassVsHPtAllSel[i] = new TH3F(histName.Data(),histTitle.Data(),nBinsPt,minPt,maxPt,nBinsM,minM,maxM,nBinsPtH,minPtH,maxPtH);
       fOutput->Add(fh3PtJet1VsMassVsHPtAllSel[i]);
 
+      histName = TString::Format("fh3PtJet1VsMassVsHPtAllSelMatch_%d",i);
+      histTitle = TString::Format("%s;#it{p}_{T,jet1};#it{M}_{jet1};#it{p}_{T,h}",histName.Data());
+      fh3PtJet1VsMassVsHPtAllSelMatch[i] = new TH3F(histName.Data(),histTitle.Data(),nBinsPt,minPt,maxPt,nBinsM,minM,maxM,nBinsPtH,minPtH,maxPtH);
+      fOutput->Add(fh3PtJet1VsMassVsHPtAllSelMatch[i]);
+
       histName = TString::Format("fh3PtJet1VsMassVsHPtTagged_%d",i);
       histTitle = TString::Format("%s;#it{p}_{T,jet1};#it{M}_{jet1};#it{p}_{T,h}",histName.Data());
       fh3PtJet1VsMassVsHPtTagged[i] =  new TH3F(histName.Data(),histTitle.Data(),nBinsPt,minPt,maxPt,nBinsM,minM,maxM,nBinsPtH,minPtH,maxPtH);
@@ -211,6 +228,11 @@ namespace EmcalHJetMassAnalysis {
       fh3PtJet1VsRatVsHPtAllSel[i] = new TH3F(histName.Data(),histTitle.Data(),nBinsPt,minPt,maxPt,nBinsR,minR,maxR,nBinsPtH,minPtH,maxPtH);
       fOutput->Add(fh3PtJet1VsRatVsHPtAllSel[i]);
 
+      histName = TString::Format("fh3PtJet1VsRatVsHPtAllSelMatch_%d",i);
+      histTitle = TString::Format("%s;#it{p}_{T,jet1};#it{M}_{jet1}/#it{p}_{T,jet1};#it{p}_{T,h}",histName.Data());
+      fh3PtJet1VsRatVsHPtAllSelMatch[i] = new TH3F(histName.Data(),histTitle.Data(),nBinsPt,minPt,maxPt,nBinsR,minR,maxR,nBinsPtH,minPtH,maxPtH);
+      fOutput->Add(fh3PtJet1VsRatVsHPtAllSelMatch[i]);
+
       histName = TString::Format("fh3PtJet1VsRatVsHPtTagged_%d",i);
       histTitle = TString::Format("%s;#it{p}_{T,jet1};#it{M}_{jet1}/#it{p}_{T,jet1};#it{p}_{T,h}",histName.Data());
       fh3PtJet1VsRatVsHPtTagged[i] =  new TH3F(histName.Data(),histTitle.Data(),nBinsPt,minPt,maxPt,nBinsR,minR,maxR,nBinsPtH,minPtH,maxPtH);
@@ -257,37 +279,45 @@ namespace EmcalHJetMassAnalysis {
       }
     }
     else if(fTriggerTrackType==kSingleInclusive) {
-      TArrayI arr; arr.Set(pCont->GetNParticles());
-      for(Int_t it = 0; it<fPtTTMin->GetSize(); it++) {
-        Int_t counter = -1;
-        arr.Reset();
-        pCont->ResetCurrentID();
-        while((vp = pCont->GetNextAcceptParticle())) {
-          if(vp->Pt()>fPtTTMin->At(it) && vp->Pt()<fPtTTMax->At(it) ) {
-            counter++;
-            arr.SetAt(pCont->GetCurrentID(),counter);
-          }
-        }
-        if(counter<0) continue;
-        //select trigger track randomly
-        fRandom->SetSeed(arr.At(0)); //random selection reproducible
-        Double_t rnd = fRandom->Uniform() * counter;
-        Int_t trigID = arr.At(TMath::FloorNint(rnd));
-        vp = pCont->GetParticle(trigID);
-        fh1PtHadron[fCentBin]->Fill(vp->Pt()); //all hadrons
-        AliEmcalJet* jet = NULL;
-        jCont->ResetCurrentID();
-        while((jet = jCont->GetNextAcceptJet())) {
-          Double_t dphi = GetDeltaPhi(vp,jet)-TMath::Pi();
-          fh3PtHPtJDPhi[fCentBin]->Fill(vp->Pt(),jet->Pt() - GetRhoVal(fContainerBase)*jet->Area(),dphi);
-          if(dphi>fDPhiHJetMax) continue;
-          FillHJetHistograms(vp->Pt(),jet);
-        }
-      }
+     for(Int_t it = 0; it<fPtTTMin->GetSize(); it++) {
+       vp = GetSingleInclusiveTT(pCont,fPtTTMin->At(it),fPtTTMax->At(it));
+       if(!vp) continue;
+       fh1PtHadron[fCentBin]->Fill(vp->Pt()); //all trigger tracks
+       AliEmcalJet* jet = NULL;
+       jCont->ResetCurrentID();
+       while((jet = jCont->GetNextAcceptJet())) {
+         Double_t dphi = GetDeltaPhi(vp,jet)-TMath::Pi();
+         fh3PtHPtJDPhi[fCentBin]->Fill(vp->Pt(),jet->Pt() - GetRhoVal(fContainerBase)*jet->Area(),dphi);
+         if(dphi>fDPhiHJetMax) continue;
+         FillHJetHistograms(vp->Pt(),jet);
+       }
+     }//trigger track types
     }
     return kTRUE;
   }
 
+  //________________________________________________________________________
+  AliVParticle* AliAnalysisTaskEmcalHJetMass::GetSingleInclusiveTT(AliParticleContainer *pCont, Double_t ptmin, Double_t ptmax) const {
+    AliVParticle *vp;
+    TArrayI arr; arr.Set(pCont->GetNParticles());
+    arr.Reset();
+    Int_t counter = -1;
+    pCont->ResetCurrentID();
+    while((vp = pCont->GetNextAcceptParticle())) {
+      if(vp->Pt()>=ptmin && vp->Pt()<ptmax ) {
+        counter++;
+        arr.SetAt(pCont->GetCurrentID(),counter);
+      }
+    }
+    if(counter<0) return NULL;
+    //select trigger track randomly
+    fRandom->SetSeed(arr.At(0)); //random selection reproducible
+    Double_t rnd = fRandom->Uniform() * counter;
+    Int_t trigID = arr.At(TMath::FloorNint(rnd));
+    vp = pCont->GetParticle(trigID);
+    return vp;
+  }
+
   //________________________________________________________________________
   Bool_t AliAnalysisTaskEmcalHJetMass::FillHJetHistograms(Double_t pt, const AliEmcalJet *jet)
   {
@@ -300,11 +330,6 @@ namespace EmcalHJetMassAnalysis {
     fh3PtJet1VsMassVsHPtAllSel[fCentBin]->Fill(ptJet,mJet,pt);
     fh3PtJet1VsRatVsHPtAllSel[fCentBin]->Fill(ptJet,rat,pt);
 
-    if(jet->GetTagStatus()<1 || !jet->GetTaggedJet())
-      return kFALSE;
-    fh3PtJet1VsMassVsHPtTagged[fCentBin]->Fill(ptJet,mJet,pt);
-    fh3PtJet1VsRatVsHPtTagged[fCentBin]->Fill(ptJet,rat,pt);
-
     Double_t fraction = -1.;
     if(fUseUnsubJet) {
       AliEmcalJet *jetUS = NULL;
@@ -327,6 +352,17 @@ namespace EmcalHJetMassAnalysis {
       fraction = jetCont->GetFractionSharedPt(jet);
     }
   
+    if(fMinFractionShared>0. && fraction>fMinFractionShared) {
+      fh3PtJet1VsMassVsHPtAllSelMatch[fCentBin]->Fill(ptJet,mJet,pt);
+      fh3PtJet1VsRatVsHPtAllSelMatch[fCentBin]->Fill(ptJet,rat,pt);
+    }
+
+    if(jet->GetTagStatus()<1 || !jet->GetTaggedJet())
+      return kFALSE;
+
+    fh3PtJet1VsMassVsHPtTagged[fCentBin]->Fill(ptJet,mJet,pt);
+    fh3PtJet1VsRatVsHPtTagged[fCentBin]->Fill(ptJet,rat,pt);
+
     if(fMinFractionShared>0. && fraction>fMinFractionShared) {
       fh3PtJet1VsMassVsHPtTaggedMatch[fCentBin]->Fill(ptJet,mJet,pt);
       fh3PtJet1VsRatVsHPtTaggedMatch[fCentBin]->Fill(ptJet,rat,pt);
index ec0dbfe857856e68c28cd56d5371f9b4ffa2e551..da4793e4f5e449b662ce479d1d34e8a265dfb612 100644 (file)
@@ -13,6 +13,7 @@ class TRandom3;
 class AliAnalysisManager;
 class AliJetContainer;
 class AliEmcalJet;
+class AliVParticle;
 
 #include "AliAnalysisTaskEmcalJet.h"
 
@@ -54,6 +55,7 @@ namespace EmcalHJetMassAnalysis {
     Double_t                            GetJetMass(const AliEmcalJet *jet) const;
     Double_t                            GetDeltaPhi(const AliVParticle *vp, const AliEmcalJet* jet) const;
     Double_t                            GetDeltaPhi(Double_t phi1,Double_t phi2) const; 
+    AliVParticle                       *GetSingleInclusiveTT(AliParticleContainer *pCont, Double_t ptmin, Double_t ptmax) const;
 
     Int_t                               fContainerBase;              // jets to be analyzed
     Int_t                               fContainerUnsub;             // unsubtracted jets
@@ -69,10 +71,12 @@ namespace EmcalHJetMassAnalysis {
     TH1F            **fh1PtHadron;                        //!pt of hadrons
     TH3F            **fh3PtHPtJDPhi;                      //!pt hadron vs pt jet vs delta phi
     TH3F            **fh3PtJet1VsMassVsHPtAllSel;         //!all jets after std selection pt vs mass vs track pt
+    TH3F            **fh3PtJet1VsMassVsHPtAllSelMatch;    //!all jets after std selection pt vs mass vs track pt matched to MC
     TH3F            **fh3PtJet1VsMassVsHPtTagged;         //!tagged jets pt vs mass vs track pt
     TH3F            **fh3PtJet1VsMassVsHPtTaggedMatch;    //!tagged jets pt vs mass vs track pt matched to MC
 
     TH3F            **fh3PtJet1VsRatVsHPtAllSel;          //!all jets after std selection pt vs mass/pt vs track pt
+    TH3F            **fh3PtJet1VsRatVsHPtAllSelMatch;     //!all jets after std selection pt vs mass/pt vs track pt matched to MC
     TH3F            **fh3PtJet1VsRatVsHPtTagged;          //!tagged jets pt vs mass/pt vs track pt
     TH3F            **fh3PtJet1VsRatVsHPtTaggedMatch;     //!tagged jets pt vs mas/pts vs track pt matched to MC
 
@@ -80,7 +84,7 @@ namespace EmcalHJetMassAnalysis {
     AliAnalysisTaskEmcalHJetMass(const AliAnalysisTaskEmcalHJetMass&);            // not implemented
     AliAnalysisTaskEmcalHJetMass &operator=(const AliAnalysisTaskEmcalHJetMass&); // not implemented
 
-    ClassDef(AliAnalysisTaskEmcalHJetMass, 3)
+    ClassDef(AliAnalysisTaskEmcalHJetMass, 4)
       };
 }
 #endif