1 #include "AliAnalysisTaskChargedJetsPA.h"
4 //TODO: FillHistogram can be done better with virtual TH1(?)
5 ClassImp(AliAnalysisTaskChargedJetsPA)
7 // ######################################################################################## DEFINE HISTOGRAMS
8 void AliAnalysisTaskChargedJetsPA::Init()
11 AliInfo("Creating histograms.");
13 // NOTE: Track & Cluster & QA histograms
17 AddHistogram1D<TH1D>("hNumberEvents", "Number of events (0 = before, 1 = after vertex cuts)", "", 2, 0, 2, "#Delta z(cm)","N^{Events}/cut");
19 AddHistogram1D<TH1D>("hAppliedEtaCorrectionFactor", "Applied #eta correction factor for the k_{T} background", "", 500, 0.5, 1.5, "Correction factor","dN^{Jets}/df");
20 AddHistogram1D<TH1D>("hAppliedEtaCorrectionFactor2", "Applied #eta correction factor for the k_{T} background 2", "", 500, 0.5, 1.5, "Correction factor","dN^{Jets}/df");
21 AddHistogram1D<TH1D>("hVertexZ", "Z distribution of the vertex", "", 400, -40., 40., "#Delta z(cm)","dN^{Events}/dz");
23 AddHistogram1D<TH1D>("hVertexR", "R distribution of the vertex", "", 100, 0., 1., "#Delta r(cm)","dN^{Events}/dr");
24 AddHistogram1D<TH1D>("hCentrality", "Centrality distribution", "", 5, 0., 100., "Centrality (classes)","dN^{Events}");
26 AddHistogram2D<TH2D>("hTrackCountAcc", "Number of tracks in acceptance vs. centrality", "LEGO2", 750, 0., 750., 5, 0, 100, "N tracks","Centrality", "dN^{Events}/dN^{Tracks}");
27 AddHistogram2D<TH2D>("hTrackPhiEta", "Track angular distribution", "LEGO2", 100, 0., 2*TMath::Pi(),100, -2.5, 2.5, "#phi","#eta","dN^{Tracks}/(d#phi d#eta)");
28 AddHistogram1D<TH1D>("hTrackPt", "Tracks p_{T} distribution", "", 20000, 0., 200., "p_{T} (GeV/c)","dN^{Tracks}/dp_{T}");
29 AddHistogram1D<TH1D>("hTrackCharge", "Charge", "", 11, -5, 5, "Charge (e)","dN^{Tracks}/dq");
30 AddHistogram1D<TH1D>("hTrackEta", "Track #eta distribution", "", 180, -fTrackEtaWindow, +fTrackEtaWindow, "#eta","dN^{Tracks}/d#eta");
32 AddHistogram2D<TH2D>("hClusterCountAcc", "Number of clusters in acceptance vs. centrality", "LEGO2", 750, 0., 750., 5, 0, 100, "N clusters","Centrality", "dN^{Events}/dN^{Clusters}");
33 AddHistogram1D<TH1D>("hClusterE", "Clusters energy distribution", "", 20000, 0., 200., "p_{T} (GeV/c)","dN^{Cluster}/dp_{T}");
36 // NOTE: Pythia histograms
39 AddHistogram1D<TH1D>("hPythiaPtHard", "Pythia p_{T} hard distribution", "", 2000, 0, 400, "p_{T} hard","dN^{Events}/dp_{T,hard}");
40 AddHistogram1D<TProfile>("hPythiaXSection", "Pythia cross section distribution", "", fNumPtHardBins+2, -1, fNumPtHardBins+1, "p_{T} hard bin","dN^{Events}/dp_{T,hard}");
41 AddHistogram1D<TH1D>("hPythiaNTrials", "Pythia trials (no correction for manual cuts)", "", fNumPtHardBins+2, -1, fNumPtHardBins+1, "p_{T} hard bin", "Trials");
44 // NOTE: Jet histograms
47 // ######## Jet spectra
48 AddHistogram1D<TH1D>("hJetPt", "Jets p_{T} distribution", "", 1000, 0., 200., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
49 AddHistogram1D<TH1D>("hJetPtBgrdSubtractedRC", "Jets p_{T} distribution, RC background subtracted", "", 500, -50., 200., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
50 AddHistogram1D<TH1D>("hJetPtBgrdSubtractedKT", "Jets p_{T} distribution, KT background subtracted, corrected for eta dependence)", "", 500, -50., 200., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
51 AddHistogram1D<TH1D>("hJetPtBgrdSubtractedKTNoEtaCorr", "Jets p_{T} distribution, KT background subtracted", "", 500, -50., 200., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
52 AddHistogram1D<TH1D>("hJetPtBgrdSubtractedKT2", "Jets p_{T} distribution, KT background 2 subtracted, corrected for eta dependence)", "", 500, -50., 200., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
53 AddHistogram1D<TH1D>("hJetPtBgrdSubtractedKT2NoEtaCorr", "Jets p_{T} distribution, KT background 2 subtracted", "", 500, -50., 200., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
55 AddHistogram1D<TH1D>("hJetPtBgrdSubtractedTR", "Jets p_{T} distribution, Track background subtracted", "", 500, -50., 200., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
58 AddHistogram1D<TH1D>("hJetArea", "Jets area distribution", "", 200, 0., 2., "Area","dN^{Jets}/dA");
59 AddHistogram2D<TH2D>("hJetPtArea", "Jets p_{T} distribution", "LEGO2", 1000, 0., 200.,100, 0., 1., "p_{T} (GeV/c)","Area","dN^{Jets}/(dp_{T}dA)");
60 AddHistogram1D<TH1D>("hJetDeltaPhi", "Jets combinatorial #Delta #phi", "", 250, 0., TMath::Pi(), "#Delta #phi","dN^{Jets}/d(#Delta #phi)");
61 AddHistogram2D<TH2D>("hJetDeltaPhiPt", "Jets combinatorial #Delta #phi vs. p_{T}", "LEGO2", 250, 0., TMath::Pi(), 20, 0.,100., "#Delta #phi","max(p_{T,1},p_{T,2}) (GeV/c)","dN^{Jets}/d(#Delta #phi)dp_{T}");
62 AddHistogram1D<TH1D>("hLeadingJetDeltaPhi", "1st and 2nd leading jet #Delta #phi", "", 250, 0., TMath::Pi(), "#Delta #phi","dN^{Jets}/d(#Delta #phi)");
63 AddHistogram2D<TH2D>("hLeadingJetDeltaPhiPt", "1st and 2nd leading jet #Delta #phi vs. p_{T}", "LEGO2", 250, 0., TMath::Pi(),20, 0.,100., "#Delta #phi","1st leading p_{T} (GeV/c)","dN^{Jets}/d(#Delta #phi)dp_{T}");
64 AddHistogram2D<TH2D>("hJetPtEta", "Jets p_{T} distribution", "LEGO2", 1000, 0., 200.,100, -0.6, 0.6, "p_{T} (GeV/c)","#eta","dN^{Jets}/(dp_{T}d#eta)");
65 AddHistogram2D<TH2D>("hJetPtPhi", "Jets p_{T} #phi distribution", "LEGO2", 1000, 0., 200.,100, 0.0, TMath::TwoPi(), "p_{T} (GeV/c)","#phi","dN^{Jets}/(dp_{T}d#phi)");
66 AddHistogram2D<TH2D>("hJetPtCentrality", "Jets p_{T} distribution", "LEGO2", 1000, 0., 200.,5, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Jets}/dp_{T}");
67 AddHistogram2D<TH2D>("hJetPhiEta", "Jets angular distribution", "LEGO2", 100, 0., 2*TMath::Pi(),100, -0.6, 0.6, "#phi","#eta","dN^{Jets}/(d#phi d#eta)");
68 AddHistogram1D<TH1D>("hJetCountAll", "Number of Jets", "", 200, 0., 200., "N jets","dN^{Events}/dN^{Jets}");
69 AddHistogram1D<TH1D>("hJetCountAccepted", "Number of accepted Jets", "", 200, 0., 200., "N jets","dN^{Events}/dN^{Jets}");
71 AddHistogram1D<TH1D>("hLeadingJetPt", "Leading jet p_{T}", "", 500, 0, 100, "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
72 AddHistogram1D<TH1D>("hSecondLeadingJetPt", "Second Leading jet p_{T}", "", 500, 0, 100, "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
74 AddHistogram1D<TH1D>("hDijetConstituentsPt", "Dijet constituents p_{T} distribution", "", 500, 0., 100., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
75 AddHistogram1D<TH1D>("hDijetLeadingJetPt", "Dijet leading jet p_{T} distribution", "", 500, 0., 100., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
76 AddHistogram2D<TH2D>("hDijetPtCorrelation", "Dijet constituents p_{T} correlation", "COLZ", 500, 5., 100., 500, 5., 100., "1st leading jet p_{T} (GeV/c)","2nd leading jet p_{T} (GeV/c)","dN^{Dijets}/d^{2}p_{T}");
78 AddHistogram1D<TH1D>("hDijet2ConstituentsPt", "Dijet2 constituents p_{T} distribution", "", 500, 0., 100., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
79 AddHistogram1D<TH1D>("hDijet2LeadingJetPt", "Dijet2 leading jet p_{T} distribution", "", 500, 0., 100., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
80 AddHistogram2D<TH2D>("hDijet2PtCorrelation", "Dijet2 constituents p_{T} correlation", "COLZ", 500, 5., 100., 500, 5., 100., "1st leading jet p_{T} (GeV/c)","2nd leading jet p_{T} (GeV/c)","dN^{Dijets}/d^{2}p_{T}");
83 // NOTE: Jet background histograms
84 if (fAnalyzeBackground)
86 // ########## Delta Pt
87 AddHistogram1D<TH1D>("hDeltaPtKT", "Background fluctuations #delta p_{T} (KT, 0 jets excluded)", "", 500, -20., 80., "#delta p_{T} (GeV/c)","dN^{Jets}/d#delta p_{T}");
88 AddHistogram1D<TH1D>("hDeltaPtKT1Excl", "Background fluctuations #delta p_{T} (KT, 1 jets excluded)", "", 500, -20., 80., "#delta p_{T} (GeV/c)","dN^{Jets}/d#delta p_{T}");
89 AddHistogram1D<TH1D>("hDeltaPtKT2Excl", "Background fluctuations #delta p_{T} (KT, 2 jets excluded)", "", 500, -20., 80., "#delta p_{T} (GeV/c)","dN^{Jets}/d#delta p_{T}");
91 Double_t dptEtaMin = -(fTrackEtaWindow-fRandConeRadius) + 2*(fTrackEtaWindow-fRandConeRadius)/fBackgroundEtaBins * fKTDeltaPtEtaBin;
92 Double_t dptEtaMax = -(fTrackEtaWindow-fRandConeRadius) + 2*(fTrackEtaWindow-fRandConeRadius)/fBackgroundEtaBins * (fKTDeltaPtEtaBin+1);
94 AddHistogram1D<TH1D>("hDeltaPtKTEta", Form("Background fluctuations #delta p_{T} (KT, 0 jets excluded, #eta=%1.3f to %1.3f)", dptEtaMin,dptEtaMax), "", 500, -20., 80., "#delta p_{T} (GeV/c)","dN^{Jets}/d#delta p_{T}");
95 AddHistogram1D<TH1D>("hDeltaPtKTEta1Excl", Form("Background fluctuations #delta p_{T} (KT, 1 jets excluded, #eta=%1.3f to %1.3f)", dptEtaMin,dptEtaMax), "", 500, -20., 80., "#delta p_{T} (GeV/c)","dN^{Jets}/d#delta p_{T}");
96 AddHistogram1D<TH1D>("hDeltaPtKTEta2Excl", Form("Background fluctuations #delta p_{T} (KT, 2 jets excluded, #eta=%1.3f to %1.3f)", dptEtaMin,dptEtaMax), "", 500, -20., 80., "#delta p_{T} (GeV/c)","dN^{Jets}/d#delta p_{T}");
98 AddHistogram1D<TH1D>("hDeltaPtKT2Eta2Excl", Form("Background fluctuations #delta p_{T} (KT 2, 2 jets excluded, #eta=%1.3f to %1.3f)", dptEtaMin,dptEtaMax), "", 500, -20., 80., "#delta p_{T} (GeV/c)","dN^{Jets}/d#delta p_{T}");
100 AddHistogram1D<TH1D>("hDeltaPtRC", "Background fluctuations #delta p_{T} (RC, 0 jets excluded)", "", 500, -20., 80., "#delta p_{T} (GeV/c)","dN^{Jets}/d#delta p_{T}");
101 AddHistogram1D<TH1D>("hDeltaPtRC1Excl", "Background fluctuations #delta p_{T} (RC, 1 jets excluded)", "", 500, -20., 80., "#delta p_{T} (GeV/c)","dN^{Jets}/d#delta p_{T}");
102 AddHistogram1D<TH1D>("hDeltaPtRC2Excl", "Background fluctuations #delta p_{T} (RC, 2 jets excluded)", "", 500, -20., 80., "#delta p_{T} (GeV/c)","dN^{Jets}/d#delta p_{T}");
104 AddHistogram1D<TH1D>("hDeltaPtTR", "Background fluctuations #delta p_{T} (TR, 0 jets excluded)", "", 500, -20., 80., "#delta p_{T} (GeV/c)","dN^{Jets}/d#delta p_{T}");
105 AddHistogram1D<TH1D>("hDeltaPtTR1Excl", "Background fluctuations #delta p_{T} (TR, 1 jets excluded)", "", 500, -20., 80., "#delta p_{T} (GeV/c)","dN^{Jets}/d#delta p_{T}");
106 AddHistogram1D<TH1D>("hDeltaPtTR2Excl", "Background fluctuations #delta p_{T} (TR, 2 jets excluded)", "", 500, -20., 80., "#delta p_{T} (GeV/c)","dN^{Jets}/d#delta p_{T}");
110 AddHistogram2D<TH2D>("hKTJetPhiEta", "KT Jets angular distribution", "LEGO2", 100, 0., 2*TMath::Pi(),100, -0.6, 0.6, "#phi","#eta","dN^{Jets}/(d#phi d#eta)");
111 AddHistogram2D<TH2D>("hKTLeadingJetPhiEta", "KT Leading jets angular distribution", "LEGO2", 100, 0., 2*TMath::Pi(),100, -0.6, 0.6, "#phi","#eta","dN^{Jets}/(d#phi d#eta)");
113 AddHistogram1D<TH1D>("hDijetBackground", "Background density (dijets excluded)", "", 400, 0., 40., "#rho (GeV/c)","dN^{Events}/d#rho");
114 AddHistogram1D<TH1D>("hDijetBackgroundMostCentral", "Background density (0-20 centrality, dijets excluded)", "", 400, 0., 40., "#rho (GeV/c)","dN^{Events}/d#rho");
115 AddHistogram2D<TH2D>("hDijetBackgroundVsCentrality", "Background density vs. centrality (dijets excluded)", "", 200, 0., 20., 5, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
117 AddHistogram1D<TH1D>("hDijetBackgroundPerpendicular", "Background density (dijets excluded)", "", 400, 0., 40., "#rho (GeV/c)","dN^{Events}/d#rho");
118 AddHistogram1D<TH1D>("hDijetBackgroundPerpendicularMostCentral", "Background density (0-20 centrality, dijets excluded)", "", 400, 0., 40., "#rho (GeV/c)","dN^{Events}/d#rho");
119 AddHistogram2D<TH2D>("hDijetBackgroundPerpendicularVsCentrality", "Background density vs. centrality (dijets excluded)", "", 200, 0., 20., 5, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
121 AddHistogram2D<TH2D>("hRCBackground", "RC background density (2 leading jets excluded, mean(8 RCs))", "LEGO2", fBackgroundEtaBins, -(fTrackEtaWindow-fRandConeRadius), +(fTrackEtaWindow-fRandConeRadius), 400, 0., 40., "#eta", "#rho (GeV/c)","dN^{Events}/d#rho");
123 AddHistogram1D<TH1D>("hAccConesInRCBackground", Form("Number of cones used for RC background (|#eta| < %1.1f)", fSignalJetEtaWindow), "", 8, 0, 8, "Used cones", "dN^{Events}/dN^{Cones}");
125 AddHistogram2D<TH2D>("hRCBackgroundMostCentral", "RC background density (0-20 centrality, 2 leading jets excluded)", "LEGO2", fBackgroundEtaBins, -(fTrackEtaWindow-fRandConeRadius), +(fTrackEtaWindow-fRandConeRadius), 400, 0., 40., "#eta", "#rho (GeV/c)","dN^{Events}/d#rho");
126 AddHistogram2D<TH2D>("hRCBackgroundMostPeripheral", "RC background density (80-100 centrality, 2 leading jets excluded)", "LEGO2", fBackgroundEtaBins, -(fTrackEtaWindow-fRandConeRadius), +(fTrackEtaWindow-fRandConeRadius), 400, 0., 40., "#eta", "#rho (GeV/c)","dN^{Events}/d#rho");
127 AddHistogram2D<TH2D>("hRCBackgroundVsCentrality", "RC background density vs centrality (2 leading jets excluded)", "LEGO2", 200, 0., 20., 5, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
130 AddHistogram2D<TH2D>("hKTBackground", "KT background density (2 leading jets excluded)", "LEGO2", fBackgroundEtaBins, -(fTrackEtaWindow-fRandConeRadius), +(fTrackEtaWindow-fRandConeRadius), 400, 0., 40., "#eta", "#rho (GeV/c)","dN^{Events}/d#rho");
131 AddHistogram2D<TH2D>("hKTBackgroundMostCentral", "KT background density (0-20 centrality, 2 leading jets excluded)", "LEGO2", fBackgroundEtaBins, -(fTrackEtaWindow-fRandConeRadius), +(fTrackEtaWindow-fRandConeRadius), 400, 0., 40., "#eta", "#rho (GeV/c)","dN^{Events}/d#rho");
132 AddHistogram2D<TH2D>("hKTBackgroundMostPeripheral", "KT background density (80-100 centrality, 2 leading jets excluded)", "LEGO2", fBackgroundEtaBins, -(fTrackEtaWindow-fRandConeRadius), +(fTrackEtaWindow-fRandConeRadius), 400, 0., 40., "#eta", "#rho (GeV/c)","dN^{Events}/d#rho");
133 AddHistogram2D<TH2D>("hKTBackgroundVsCentrality", "KT background density vs centrality (2 leading jets excluded)", "LEGO2", 200, 0., 20., 5, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
135 AddHistogram2D<TH2D>("hKTBackground2", "KT background 2 density (2 leading jets excluded)", "LEGO2", fBackgroundEtaBins, -(fTrackEtaWindow-fRandConeRadius), +(fTrackEtaWindow-fRandConeRadius), 400, 0., 40., "#eta", "#rho (GeV/c)","dN^{Events}/d#rho");
136 AddHistogram2D<TH2D>("hKTBackground2MostCentral", "KT background 2 density (0-20 centrality, 2 leading jets excluded)", "LEGO2", fBackgroundEtaBins, -(fTrackEtaWindow-fRandConeRadius), +(fTrackEtaWindow-fRandConeRadius), 400, 0., 40., "#eta", "#rho (GeV/c)","dN^{Events}/d#rho");
137 AddHistogram2D<TH2D>("hKTBackground2MostPeripheral", "KT background 2 density (80-100 centrality, 2 leading jets excluded)", "LEGO2", fBackgroundEtaBins, -(fTrackEtaWindow-fRandConeRadius), +(fTrackEtaWindow-fRandConeRadius), 400, 0., 40., "#eta", "#rho (GeV/c)","dN^{Events}/d#rho");
138 AddHistogram2D<TH2D>("hKTBackground2VsCentrality", "KT background 2 density vs centrality (2 leading jets excluded)", "LEGO2", 200, 0., 20., 5, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
141 AddHistogram2D<TH2D>("hTrackBackground", "Track background density (2 leading jets excluded)", "LEGO2", fBackgroundEtaBins, -(fTrackEtaWindow-fRandConeRadius), +(fTrackEtaWindow-fRandConeRadius), 400, 0., 40., "#eta", "#rho (GeV/c)","dN^{Events}/d#rho");
142 AddHistogram2D<TH2D>("hTrackBackgroundMostCentral", "Track background density (0-20 centrality, 2 leading jets excluded)", "LEGO2", fBackgroundEtaBins, -(fTrackEtaWindow-fRandConeRadius), +(fTrackEtaWindow-fRandConeRadius), 400, 0., 40., "#eta", "#rho (GeV/c)","dN^{Events}/d#rho");
143 AddHistogram2D<TH2D>("hTrackBackgroundMostPeripheral", "Track background density (80-100 centrality, 2 leading jets excluded)", "LEGO2", fBackgroundEtaBins, -(fTrackEtaWindow-fRandConeRadius), +(fTrackEtaWindow-fRandConeRadius), 400, 0., 40., "#eta", "#rho (GeV/c)","dN^{Events}/d#rho");
144 AddHistogram2D<TH2D>("hTrackBackgroundVsCentrality", "Track background density vs centrality (2 leading jets excluded)", "LEGO2", 200, 0., 20., 5, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
150 // register Histograms
151 for (Int_t i = 0; i < fHistCount; i++)
153 fOutputList->Add(fHistList->At(i));
156 PostData(1,fOutputList); // important for merging
160 //________________________________________________________________________
161 AliAnalysisTaskChargedJetsPA::AliAnalysisTaskChargedJetsPA(const char *name, const char* trackArrayName, const char* clusterArrayName, const char* jetArrayName, const char* backgroundJetArrayName) : AliAnalysisTaskSE(name), fOutputList(0), fAnalyzeQA(1), fAnalyzeJets(1), fAnalyzeBackground(1), fAnalyzePythia(0), fHasTracks(0), fHasClusters(0), fHasJets(0), fHasBackgroundJets(0), fIsMC(0), fJetArray(0), fTrackArray(0), fClusterArray(0), fBackgroundJetArray(0), fJetArrayName(0), fTrackArrayName(0), fClusterArrayName(0), fBackgroundJetArrayName(0), fNumPtHardBins(11), fRandConeRadius(0.4), fSignalJetRadius(0.4), fBackgroundJetRadius(0.4), fKTDeltaPtEtaBin(3), fTrackBackgroundConeRadius(0.4), fNumberRandCones(8), fNumberExcludedJets(2), fDijetMaxAngleDeviation(10.0), fBackgroundEtaBins(5), fJetBgrdCorrectionFactors(0), fSignalJetEtaWindow(0.5), fBackgroundJetEtaWindow(0.5), fTrackEtaWindow(0.9), fClusterEtaWindow(0.7), fVertexWindow(10.0), fVertexMaxR(1.0), fMinTrackPt(0.150), fMinClusterPt(0.300), fMinJetPt(1.0), fMinJetArea(0.4), fMinBackgroundJetPt(0.15), fMinDijetLeadingPt(10.0), fFirstLeadingJet(0), fSecondLeadingJet(0), fNumberSignalJets(0), fCrossSection(0.0), fTrials(0.0), fRandom(0), fInitialized(0), fTaskInstanceCounter(0), fHistList(0), fHistCount(0)
164 AliInfo("Calling constructor.");
168 // Define input and output slots here (never in the dummy constructor)
169 // Input slot #0 works with a TChain - it is connected to the default input container
170 // Output slot #1 writes into a TH1 container
173 // Every instance of this task gets his own number
174 static Int_t instance = 0;
175 fTaskInstanceCounter = instance;
178 fTrackArrayName = new TString(trackArrayName);
179 fClusterArrayName = new TString(clusterArrayName);
180 if (strcmp(fTrackArrayName->Data(),"") == 0)
185 if (fTrackArrayName->Contains("MCParticles")) //TODO: Hardcoded for now
189 fJetArrayName = new TString(jetArrayName);
190 if (strcmp(fJetArrayName->Data(),"") == 0)
191 fAnalyzeJets = kFALSE;
193 fAnalyzeJets = kTRUE;
195 fBackgroundJetArrayName = new TString(backgroundJetArrayName);
196 if (strcmp(fBackgroundJetArrayName->Data(),"") == 0)
197 fAnalyzeBackground = kFALSE;
199 fAnalyzeBackground = kTRUE;
201 DefineOutput(1, TList::Class());
203 fHistList = new TList();
206 AliInfo("Constructor done.");
211 //________________________________________________________________________
212 inline Double_t AliAnalysisTaskChargedJetsPA::GetConePt(Double_t eta, Double_t phi, Double_t radius)
214 Double_t tmpConePt = 0.0;
216 for (Int_t i = 0; i < fTrackArray->GetEntries(); i++)
218 AliVTrack* tmpTrack = static_cast<AliVTrack*>(fTrackArray->At(i));
219 if (IsTrackInAcceptance(tmpTrack))
220 if(IsTrackInCone(tmpTrack, eta, phi, radius))
221 tmpConePt = tmpConePt + tmpTrack->Pt();
226 //________________________________________________________________________
227 inline Double_t AliAnalysisTaskChargedJetsPA::GetPtHard()
229 Double_t tmpPtHard = -1.0;
232 AliError("MCEvent not accessible although demanded!");
235 AliGenPythiaEventHeader* pythiaHeader = dynamic_cast<AliGenPythiaEventHeader*>(MCEvent()->GenEventHeader());
237 AliError("Pythia Header not accessible!");
239 tmpPtHard = pythiaHeader->GetPtHard();
244 //________________________________________________________________________
245 inline Int_t AliAnalysisTaskChargedJetsPA::GetPtHardBin()
247 // ########## PT HARD BIN EDGES
248 const Int_t localkPtHardLowerEdges[] = { 0, 5,11,21,36,57, 84,117,152,191,234};
249 const Int_t localkPtHardHigherEdges[] = { 5,11,21,36,57,84,117,152,191,234,1000000};
251 Int_t tmpPtHardBin = 0;
252 Double_t tmpPtHard = GetPtHard();
254 for (tmpPtHardBin = 0; tmpPtHardBin <= fNumPtHardBins; tmpPtHardBin++)
255 if (tmpPtHard >= localkPtHardLowerEdges[tmpPtHardBin] && tmpPtHard < localkPtHardHigherEdges[tmpPtHardBin])
262 //________________________________________________________________________
263 inline Bool_t AliAnalysisTaskChargedJetsPA::IsTrackInCone(AliVTrack* track, Double_t eta, Double_t phi, Double_t radius)
265 // This is to use a full cone in phi even at the edges of phi (2pi -> 0) (0 -> 2pi)
266 Double_t trackPhi = 0.0;
267 if (track->Phi() > (TMath::TwoPi() - (radius-phi)))
268 trackPhi = track->Phi() - TMath::TwoPi();
269 else if (track->Phi() < (phi+radius - TMath::TwoPi()))
270 trackPhi = track->Phi() + TMath::TwoPi();
272 trackPhi = track->Phi();
274 if ( TMath::Abs(trackPhi-phi)*TMath::Abs(trackPhi-phi) + TMath::Abs(track->Eta()-eta)*TMath::Abs(track->Eta()-eta) <= radius*radius)
280 //________________________________________________________________________
281 inline Bool_t AliAnalysisTaskChargedJetsPA::IsTrackInAcceptance(AliVParticle* track)
284 if (TMath::Abs(track->Eta()) <= fTrackEtaWindow)
285 if (track->Pt() >= fMinTrackPt)
291 //________________________________________________________________________
292 inline Bool_t AliAnalysisTaskChargedJetsPA::IsClusterInAcceptance(AliVCluster* cluster)
295 // if (TMath::Abs(cluster->Eta()) <= fClusterEtaWindow)
296 // if (cluster->Phi() <= 187.0/360.0 * TMath::TwoPi());
297 // if (cluster->Phi() >= 80.0/360.0 * TMath::TwoPi());
298 if (cluster->E() >= fMinClusterPt)
305 //________________________________________________________________________
306 inline Bool_t AliAnalysisTaskChargedJetsPA::IsBackgroundJetInAcceptance(AliEmcalJet *jet)
309 if (TMath::Abs(jet->Eta()) <= fBackgroundJetEtaWindow)
310 if (jet->Pt() >= fMinBackgroundJetPt)
316 //________________________________________________________________________
317 inline Bool_t AliAnalysisTaskChargedJetsPA::IsSignalJetInAcceptance(AliEmcalJet *jet)
320 if (TMath::Abs(jet->Eta()) <= fSignalJetEtaWindow)
321 if (jet->Pt() >= fMinJetPt)
322 if (jet->Area() >= fMinJetArea)
328 //________________________________________________________________________
329 inline Bool_t AliAnalysisTaskChargedJetsPA::IsDijet(AliEmcalJet *jet1, AliEmcalJet *jet2)
331 // Output from GetDeltaPhi is < pi in any case
332 if ((jet1 != 0) && (jet2 != 0))
333 if((TMath::Pi() - GetDeltaPhi(jet1->Phi(),jet2->Phi())) < fDijetMaxAngleDeviation)
334 if ((jet1->Pt() > fMinDijetLeadingPt) && (jet2->Pt() > fMinDijetLeadingPt)) //TODO: Introduce recoil jet?
340 //________________________________________________________________________
341 void AliAnalysisTaskChargedJetsPA::ExecOnce()
344 AliInfo("Starting ExecOnce.");
346 fInitialized = kTRUE;
348 // Check for track array
349 if (strcmp(fTrackArrayName->Data(), "") != 0)
351 fTrackArray = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fTrackArrayName->Data()));
355 AliInfo(Form("%s: Could not retrieve tracks %s! This is OK, if tracks are not demanded.", GetName(), fTrackArrayName->Data()));
360 TClass *cl = fTrackArray->GetClass();
361 if (!cl->GetBaseClass("AliVParticle"))
363 AliError(Form("%s: Collection %s does not contain AliVParticle objects!", GetName(), fTrackArrayName->Data()));
369 // Check for cluster array
370 if (strcmp(fClusterArrayName->Data(), "") != 0)
372 fClusterArray = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fClusterArrayName->Data()));
373 fHasClusters = kTRUE;
376 AliInfo(Form("%s: Could not retrieve clusters %s! This is OK, if clusters are not demanded.", GetName(), fClusterArrayName->Data()));
377 fHasClusters = kFALSE;
381 TClass *cl = fClusterArray->GetClass();
382 if (!cl->GetBaseClass("AliVCluster"))
384 AliError(Form("%s: Collection %s does not contain AliVCluster objects!", GetName(), fClusterArrayName->Data()));
386 fHasClusters = kFALSE;
391 // Check for jet array
392 if (strcmp(fJetArrayName->Data(), "") != 0)
394 fJetArray = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fJetArrayName->Data()));
399 AliInfo(Form("%s: Could not retrieve jets %s! This is OK, if jets are not demanded.", GetName(), fJetArrayName->Data()));
404 if (!fJetArray->GetClass()->GetBaseClass("AliEmcalJet"))
406 AliError(Form("%s: Collection %s does not contain AliEmcalJet objects!", GetName(), fJetArrayName->Data()));
413 // Check for background object
414 if (strcmp(fBackgroundJetArrayName->Data(), "") != 0)
416 fBackgroundJetArray = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fBackgroundJetArrayName->Data()));
417 fHasBackgroundJets = kTRUE;
418 if (!fBackgroundJetArray)
420 AliInfo(Form("%s: Could not retrieve background jets %s! This is OK, if background is not demanded.", GetName(), fBackgroundJetArrayName->Data()));
421 fHasBackgroundJets = kFALSE;
425 // Look, if initialization is OK
426 if ((!fHasTracks && fAnalyzeQA) || (!fHasTracks && fAnalyzeBackground))
428 AliError(Form("%s: Tracks NOT successfully casted although demanded! Deactivating QA and background analysis",GetName()));
430 fAnalyzeBackground = kFALSE;
432 if ((!fHasJets && fAnalyzeJets) || (!fHasJets && fAnalyzeBackground))
434 AliError(Form("%s: Jets NOT successfully casted although demanded! Deactivating jet- and background analysis",GetName()));
435 fAnalyzeJets = kFALSE;
436 fAnalyzeBackground = kFALSE;
438 if (!fHasBackgroundJets && fAnalyzeBackground)
440 AliError(Form("%s: Background NOT successfully casted although demanded! Deactivating background analysis",GetName()));
441 fAnalyzeBackground = kFALSE;
448 AliInfo("ExecOnce done.");
453 //________________________________________________________________________
454 void AliAnalysisTaskChargedJetsPA::GetSignalJets()
457 fFirstLeadingJet = NULL;
458 fSecondLeadingJet = NULL;
459 fNumberSignalJets = 0;
461 Float_t maxJetPts[] = {0, 0};
462 Int_t jetIDArray[] = {-1, -1};
463 Int_t jetCount = fJetArray->GetEntries();
465 // Go through all jets and save signal jets and the two leading ones
466 for (Int_t i = 0; i < jetCount; i++)
468 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJetArray->At(i));
471 AliError(Form("%s: Could not receive jet %d", GetName(), i));
475 if (!IsSignalJetInAcceptance(jet)) continue;
477 if (jet->Pt() > maxJetPts[0])
479 maxJetPts[1] = maxJetPts[0];
480 jetIDArray[1] = jetIDArray[0];
481 maxJetPts[0] = jet->Pt();
484 else if (jet->Pt() > maxJetPts[1])
486 maxJetPts[1] = jet->Pt();
489 fSignalJets[fNumberSignalJets] = jet;
493 if (fNumberSignalJets > 0)
494 fFirstLeadingJet = static_cast<AliEmcalJet*>(fJetArray->At(jetIDArray[0]));
495 if (fNumberSignalJets > 1)
496 fSecondLeadingJet = static_cast<AliEmcalJet*>(fJetArray->At(jetIDArray[1]));
500 //________________________________________________________________________
501 Int_t AliAnalysisTaskChargedJetsPA::GetLeadingJets(TClonesArray* jetArray, Int_t* jetIDArray, Bool_t isSignalJets)
503 // Writes first two leading jets into already registered array jetIDArray
507 AliError("Could not get the jet array to get leading jets from it!");
511 Float_t maxJetPts[] = {0, 0};
515 Int_t jetCount = jetArray->GetEntries();
516 Int_t jetCountAccepted = 0;
518 for (Int_t i = 0; i < jetCount; i++)
520 AliEmcalJet* jet = static_cast<AliEmcalJet*>(jetArray->At(i));
523 AliError(Form("%s: Could not receive jet %d", GetName(), i));
529 if (!IsSignalJetInAcceptance(jet)) continue;
533 if (!IsBackgroundJetInAcceptance(jet)) continue;
536 if (jet->Pt() > maxJetPts[0])
538 maxJetPts[1] = maxJetPts[0];
539 jetIDArray[1] = jetIDArray[0];
540 maxJetPts[0] = jet->Pt();
543 else if (jet->Pt() > maxJetPts[1])
545 maxJetPts[1] = jet->Pt();
550 return jetCountAccepted;
552 //________________________________________________________________________
553 Double_t AliAnalysisTaskChargedJetsPA::GetJetBackgroundCorrFactor(Double_t eta, Double_t background)
555 Double_t tmpCorrFactor = 1.0;
557 if(fJetBgrdCorrectionFactors)
558 tmpCorrFactor = fJetBgrdCorrectionFactors->GetBinContent
560 fJetBgrdCorrectionFactors->GetXaxis()->FindBin(eta),
561 fJetBgrdCorrectionFactors->GetYaxis()->FindBin(background)
564 return tmpCorrFactor;
566 //________________________________________________________________________
567 Double_t AliAnalysisTaskChargedJetsPA::GetCorrectedJetPt(AliEmcalJet* jet, Double_t background, Bool_t useEtaCorrection)
570 AliInfo("Getting corrected jet spectra.");
575 AliError("Jet pointer passed to GetCorrectedJet() not valid!");
579 Double_t correctedPt = -1.0;
581 // Get correction factor from saved histo or similar in dependence of jet eta and background density
582 Double_t corrfactor = 1.0;
585 corrfactor = GetJetBackgroundCorrFactor(jet->Eta(), background);
588 // Get Eta corrected background
589 Double_t tmpCorrectedBackground = background * corrfactor;
591 // Subtract background
592 correctedPt = jet->Pt() - tmpCorrectedBackground * jet->Area();
595 AliInfo("Got corrected jet spectra.");
602 //________________________________________________________________________
603 void AliAnalysisTaskChargedJetsPA::GetDeltaPt(Double_t& deltaPt, Double_t rho, Int_t numberExcludeLeadingJets, Int_t usedEtaBin, Bool_t useEtaCorrection)
606 AliInfo("Getting Delta Pt.");
609 // Define the tmp delta pt
612 // Exclude UP TO numberExcludeLeadingJets
613 if (fNumberSignalJets < 2)
614 numberExcludeLeadingJets = fNumberSignalJets;
616 Double_t etaMin, etaMax;
619 etaMin = -(fTrackEtaWindow-fRandConeRadius);
620 etaMax = +(fTrackEtaWindow-fRandConeRadius);
624 etaMin = -(fTrackEtaWindow-fRandConeRadius) + 2*(fTrackEtaWindow-fRandConeRadius)/fBackgroundEtaBins * usedEtaBin;
625 etaMax = -(fTrackEtaWindow-fRandConeRadius) + 2*(fTrackEtaWindow-fRandConeRadius)/fBackgroundEtaBins * (usedEtaBin+1);
629 Double_t tmpRandConeEta = 0.0;
630 Double_t tmpRandConePhi = 0.0;
632 Bool_t coneValid = kTRUE;
635 tmpRandConeEta = etaMin + fRandom->Rndm()*(etaMax-etaMin);
636 tmpRandConePhi = fRandom->Rndm()*TMath::TwoPi();
638 // Apply eta correction on demand
640 rho = GetJetBackgroundCorrFactor(tmpRandConeEta, rho)*rho;
642 // Go through all excluded leading jets and check if there's an overlap
643 for(Int_t j=0;j<numberExcludeLeadingJets;j++)
645 AliEmcalJet* tmpJet = NULL;
648 tmpJet = fFirstLeadingJet;
650 tmpJet = fSecondLeadingJet;
652 AliFatal("Trying to exclude more than 2 jets for delta pt -- not implemented.");
654 Double_t excludedJetPhi = tmpJet->Phi();
655 Double_t excludedJetEta = tmpJet->Eta();
656 Double_t tmpDeltaPhi = GetDeltaPhi(tmpRandConePhi, excludedJetPhi);
658 // Check, if cone has overlap with jet
659 if ( tmpDeltaPhi*tmpDeltaPhi + TMath::Abs(tmpRandConeEta-excludedJetEta)*TMath::Abs(tmpRandConeEta-excludedJetEta) <= fRandConeRadius*fRandConeRadius)
666 // Get the cones' pt and calculate delta pt
668 deltaPt = GetConePt(tmpRandConeEta,tmpRandConePhi,fRandConeRadius) - (rho*fRandConeRadius*fRandConeRadius*TMath::Pi());
670 AliInfo("Got Delta Pt.");
674 //________________________________________________________________________
675 void AliAnalysisTaskChargedJetsPA::GetKTBackgroundDensity(Int_t numberExcludeLeadingJets, Double_t& rhoMedian, Double_t& areaMean, Double_t etaMin, Double_t etaMax)
678 AliInfo("Getting KT background density.");
681 // static declaration. Advantage: more speed. Disadvantage: Problematic for events with more than 1024 jets :)
682 static Double_t tmpRhos[1024];
683 static Double_t tmpAreas[1024];
684 Int_t maxJetIds[] = {-1, -1}; // Indices for excludes jets (up to two)
686 // Setting invalid values
690 // Exclude UP TO numberExcludeLeadingJets
691 Int_t numberBgrdJets = GetLeadingJets(fBackgroundJetArray, &maxJetIds[0], kFALSE);
692 if (numberBgrdJets < numberExcludeLeadingJets)
693 numberExcludeLeadingJets = numberBgrdJets;
694 if ((etaMin == 0) && (etaMax == 0))
696 etaMin = -fBackgroundJetEtaWindow;
697 etaMax = +fBackgroundJetEtaWindow;
700 Int_t jetCountAccepted = 0;
701 Int_t jetCount = fBackgroundJetArray->GetEntries();
703 for (Int_t i = 0; i < jetCount; i++)
705 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fBackgroundJetArray->At(i));
708 AliError(Form("%s: Could not receive jet %d", GetName(), i));
712 // exclude leading jets
713 if (numberExcludeLeadingJets > 0)
714 if (i == maxJetIds[0])
716 if (numberExcludeLeadingJets > 1)
717 if (i == maxJetIds[1])
722 if (!IsBackgroundJetInAcceptance(jet))
724 if (!((jet->Eta() >= etaMin) && (jet->Eta() < etaMax)))
728 tmpRhos[jetCountAccepted] = jet->Pt() / jet->Area();
729 tmpAreas[jetCountAccepted] = jet->Area();
733 if (jetCountAccepted > 0)
735 rhoMedian = TMath::Median(jetCountAccepted, tmpRhos);
736 areaMean = TMath::Mean(jetCountAccepted, tmpAreas);
739 AliInfo("Got KT background density.");
743 //________________________________________________________________________
744 void AliAnalysisTaskChargedJetsPA::GetKTBackground2Density(Int_t numberExcludeLeadingJets, Double_t& rhoMedian, Double_t& areaMean, Double_t etaMin, Double_t etaMax)
747 AliInfo("Getting KT background 2 density.");
750 // static declaration. Advantage: more speed. Disadvantage: Problematic for events with more than 1024 jets :)
751 static Double_t tmpRhos[1024];
752 static Double_t tmpAreas[1024];
754 // Setting invalid values
758 if ((etaMin == 0) && (etaMax == 0))
760 etaMin = -fBackgroundJetEtaWindow;
761 etaMax = +fBackgroundJetEtaWindow;
764 Int_t jetCountAccepted = 0;
765 Int_t jetCount = fBackgroundJetArray->GetEntries();
767 for (Int_t i = 0; i < jetCount; i++)
769 Bool_t jetValid = kTRUE;
770 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fBackgroundJetArray->At(i));
773 AliError(Form("%s: Could not receive jet %d", GetName(), i));
777 if (!((jet->Eta() >= etaMin) && (jet->Eta() < etaMax)))
779 if (!IsBackgroundJetInAcceptance(jet))
782 // Look, if theres an overlap of leading jets/ kT jet. If yes, exclude this jet
783 for(Int_t j=0;j<numberExcludeLeadingJets;j++)
785 AliEmcalJet* tmpLeadingJet = NULL;
788 tmpLeadingJet = fFirstLeadingJet;
790 tmpLeadingJet = fSecondLeadingJet;
792 AliFatal("Trying to exclude more than 2 jets in KT background 2 -- not implemented.");
796 Double_t tmpDeltaPhi = GetDeltaPhi(jet->Phi(), tmpLeadingJet->Phi());
797 if ( tmpDeltaPhi*tmpDeltaPhi + TMath::Abs(jet->Eta()-tmpLeadingJet->Eta())*TMath::Abs(jet->Eta()-tmpLeadingJet->Eta()) <= fBackgroundJetRadius*fBackgroundJetRadius)
808 tmpRhos[jetCountAccepted] = jet->Pt() / jet->Area();
809 tmpAreas[jetCountAccepted] = jet->Area();
813 if (jetCountAccepted > 0)
815 rhoMedian = TMath::Median(jetCountAccepted, tmpRhos);
816 areaMean = TMath::Mean(jetCountAccepted, tmpAreas);
819 AliInfo("Got KT background 2 density.");
824 //________________________________________________________________________
825 Int_t AliAnalysisTaskChargedJetsPA::GetRCBackgroundDensity(Int_t numberExcludeLeadingJets, Double_t& rhoMean, Double_t& rhoMedian, Double_t etaMin, Double_t etaMax, Int_t numberRandCones)
828 AliInfo("Getting RC background density.");
831 if(numberRandCones == 0)
832 numberRandCones = fNumberRandCones;
834 std::vector<AliEmcalJet> tmpCones(numberRandCones);
836 // Setting invalid values
840 // Exclude UP TO numberExcludeLeadingJets
841 if (fNumberSignalJets < 2)
842 numberExcludeLeadingJets = fNumberSignalJets;
844 // Search given amount of RCs
845 Int_t numAcceptedRCs = 0;
846 for(Int_t i=0;i<numberRandCones;i++)
848 Double_t tmpRandConeEta = 0.0;
849 Double_t tmpRandConePhi = 0.0;
850 Double_t excludedJetEta = 0.0;
851 Double_t excludedJetPhi = 0.0;
853 // Search random cone in acceptance with no overlap with already excluded jets (leading jets and random cones)
854 Bool_t coneValid = kTRUE;
856 // Set the random cone position
857 if ((etaMin == 0) && (etaMax == 0))
858 tmpRandConeEta = (fTrackEtaWindow-fRandConeRadius)*(2.0*fRandom->Rndm()-1.0); // full RC is in acceptance
860 tmpRandConeEta = etaMin + fRandom->Rndm()*(etaMax-etaMin);
862 tmpRandConePhi = fRandom->Rndm()*TMath::TwoPi();
864 // Go through all excluded leading jets and check if there's an overlap
866 for(Int_t j=0;j<numberExcludeLeadingJets;j++)
868 AliEmcalJet* tmpJet = NULL;
871 tmpJet = fFirstLeadingJet;
873 tmpJet = fSecondLeadingJet;
875 AliFatal("Trying to exclude more than 2 jets in RC background -- not implemented.");
877 excludedJetPhi = tmpJet->Phi();
878 excludedJetEta = tmpJet->Eta();
879 Double_t tmpDeltaPhi = GetDeltaPhi(tmpRandConePhi, excludedJetPhi);
881 if ( tmpDeltaPhi*tmpDeltaPhi + TMath::Abs(tmpRandConeEta-excludedJetEta)*TMath::Abs(tmpRandConeEta-excludedJetEta) <= fRandConeRadius*fRandConeRadius)
888 // RC is accepted, so save it
891 AliEmcalJet tmpJet(GetConePt(tmpRandConeEta, tmpRandConePhi, fRandConeRadius), tmpRandConeEta, tmpRandConePhi, 0.0);
892 tmpCones[numAcceptedRCs] = tmpJet;
897 // Calculate Rho and the mean from the RCs (no excluded jets are considered!)
898 if(numAcceptedRCs > 0)
900 std::vector<Double_t> tmpRho(numAcceptedRCs);
901 for (Int_t i=0; i<numAcceptedRCs;i++)
902 tmpRho[i] = tmpCones[i].Pt()/(fRandConeRadius*fRandConeRadius*TMath::Pi());
904 rhoMean = TMath::Mean(tmpRho.begin(), tmpRho.end());
905 rhoMedian = 0.0; // NOT IMPLEMENTED because TMath::Median is not working with iterators
909 AliInfo("Got RC background density.");
911 return numAcceptedRCs;
914 //________________________________________________________________________
915 void AliAnalysisTaskChargedJetsPA::GetTrackBackgroundDensity(Int_t numberExcludeLeadingJets, Double_t& rhoMean, Double_t& area, Double_t etaMin, Double_t etaMax)
918 AliInfo("Getting track background density.");
921 Double_t summedTracksPt = 0.0;
923 if ((etaMin == 0) && (etaMax == 0))
925 etaMin = -fTrackEtaWindow;
926 etaMax = +fTrackEtaWindow;
929 // Setting invalid values
932 // Exclude UP TO numberExcludeLeadingJets
933 if (fNumberSignalJets < 2)
934 numberExcludeLeadingJets = fNumberSignalJets;
937 Int_t trackCount = fTrackArray->GetEntries();
938 Int_t trackCountAccepted = 0;
939 for (Int_t i = 0; i < trackCount; i++)
941 Bool_t trackValid = kTRUE;
942 AliVTrack* tmpTrack = static_cast<AliVTrack*>(fTrackArray->At(i));
943 if (IsTrackInAcceptance(tmpTrack))
944 if ((tmpTrack->Eta() >= etaMin) && (tmpTrack->Eta() < etaMax))
946 for (Int_t j = 0; j < numberExcludeLeadingJets; j++)
948 AliEmcalJet* tmpJet = NULL;
950 tmpJet = fFirstLeadingJet;
952 tmpJet = fSecondLeadingJet;
954 AliFatal("Trying to exclude more than 2 jets in track background -- not implemented.");
956 if (IsTrackInCone(tmpTrack, tmpJet->Eta(), tmpJet->Phi(), fTrackBackgroundConeRadius))
964 // Add track pt to array
965 summedTracksPt = summedTracksPt + tmpTrack->Pt();
966 trackCountAccepted++;
971 if (trackCountAccepted > 0)
973 Double_t tmpArea = 0.0;
975 tmpArea = (2.0*fTrackEtaWindow) * TMath::TwoPi() * (etaMax-etaMin)/(2.0*fTrackEtaWindow); // area of the used eta strip
977 // Now: exclude the part of the leading jet that is in the strip
978 if (numberExcludeLeadingJets == 2)
979 tmpArea = tmpArea*(1.0-MCGetOverlapCircleRectancle(fFirstLeadingJet->Eta(), fFirstLeadingJet->Phi(), fTrackBackgroundConeRadius, etaMin, etaMax, 0., TMath::TwoPi()) -MCGetOverlapCircleRectancle(fSecondLeadingJet->Eta(), fSecondLeadingJet->Phi(), fTrackBackgroundConeRadius, etaMin, etaMax, 0., TMath::TwoPi()));
980 else if (numberExcludeLeadingJets == 1)
981 tmpArea = tmpArea*(1.0-MCGetOverlapCircleRectancle(fFirstLeadingJet->Eta(), fFirstLeadingJet->Phi(), fTrackBackgroundConeRadius, etaMin, etaMax, 0., TMath::TwoPi()));
983 rhoMean = summedTracksPt/tmpArea;
988 AliInfo("Got track background density.");
992 //________________________________________________________________________
993 void AliAnalysisTaskChargedJetsPA::GetTrackBackgroundDensity(Int_t numberExcludeLeadingJets, Double_t& rhoMean, Double_t& area, AliEmcalJet* excludeJet1, AliEmcalJet* excludeJet2, Bool_t doSearchPerpendicular)
996 AliInfo("Getting track background density.");
999 // Setting invalid values
1000 Double_t summedTracksPt = 0.0;
1004 Double_t tmpRadius = 0.0;
1005 if (doSearchPerpendicular)
1006 tmpRadius = 0.5*TMath::Pi(); // exclude 90 degrees around jets
1008 tmpRadius = fSignalJetRadius;
1010 numberExcludeLeadingJets = 2; // dijet is excluded here in any case
1014 if (!fTrackArray || !fJetArray)
1016 AliError("Could not get the track/jet array to calculate track rho!");
1020 Int_t trackCount = fTrackArray->GetEntries();
1021 Int_t trackCountAccepted = 0;
1022 for (Int_t i = 0; i < trackCount; i++)
1024 AliVTrack* tmpTrack = static_cast<AliVTrack*>(fTrackArray->At(i));
1025 if (IsTrackInAcceptance(tmpTrack))
1027 if (IsTrackInCone(tmpTrack, excludeJet1->Eta(), excludeJet1->Phi(), tmpRadius))
1030 if (numberExcludeLeadingJets > 1)
1031 if (IsTrackInCone(tmpTrack, excludeJet2->Eta(), excludeJet2->Phi(), tmpRadius))
1034 // Add track pt to array
1035 summedTracksPt = summedTracksPt + tmpTrack->Pt();
1036 trackCountAccepted++;
1040 if (trackCountAccepted > 0)
1042 Double_t tmpArea = 2.0*fTrackEtaWindow*TMath::TwoPi() - 2*(tmpRadius*tmpRadius * TMath::Pi()); //TPC area - excluding jet area
1043 rhoMean = summedTracksPt/tmpArea;
1048 AliInfo("Got track background density.");
1052 //________________________________________________________________________
1053 void AliAnalysisTaskChargedJetsPA::Calculate(AliVEvent* event)
1056 AliInfo("Starting Calculate().");
1058 ////////////////////// NOTE: initialization & casting
1061 AliError("??? Event pointer == 0 ???");
1066 ExecOnce(); // Get tracks, jets, background from arrays if not already given + Init Histos
1069 FillHistogram("hNumberEvents", 0.5); // number of events before manual cuts
1071 if ((event->GetPrimaryVertex()->GetNContributors() < 1) || (TMath::Abs(event->GetPrimaryVertex()->GetZ()) > fVertexWindow))
1074 if (TMath::Sqrt(event->GetPrimaryVertex()->GetX()*event->GetPrimaryVertex()->GetX() + event->GetPrimaryVertex()->GetY()*event->GetPrimaryVertex()->GetY()) > fVertexMaxR)
1077 FillHistogram("hNumberEvents", 1.5); // number of events after manual cuts
1080 AliInfo("Calculate()::Init done.");
1083 ////////////////////// NOTE: Get Centrality, (Leading)Signal jets and Background
1085 // Get centrality (V0A)
1086 AliCentrality* tmpCentrality = NULL;
1087 tmpCentrality = event->GetCentrality();
1088 Double_t centralityPercentile = 0.0;
1089 if (tmpCentrality != NULL)
1090 centralityPercentile = tmpCentrality->GetCentralityPercentile("V0A");
1093 if (fAnalyzeBackground || fAnalyzeJets)
1098 // Background with N excluded leading jets
1099 std::vector<Double_t> ktBackgroundRhoMedian(fBackgroundEtaBins+1);
1100 std::vector<Double_t> ktBackgroundAreaMean(fBackgroundEtaBins+1);
1101 std::vector<Double_t> ktBackground2RhoMedian(fBackgroundEtaBins+1);
1102 std::vector<Double_t> ktBackground2AreaMean(fBackgroundEtaBins+1);
1103 std::vector<Double_t> rcBackgroundRhoMean(fBackgroundEtaBins+1);
1104 std::vector<Double_t> rcBackgroundRhoMedian(fBackgroundEtaBins+1);
1105 std::vector<Double_t> trackBackgroundRhoMean(fBackgroundEtaBins+1);
1106 std::vector<Double_t> trackBackgroundArea(fBackgroundEtaBins+1);
1107 Double_t dijetBackground = -1.0; // calculation only done in events with dijets I!
1108 Double_t dijetBackgroundPerpendicular = -1.0; // calculation only done in events with dijets I!
1109 if (fAnalyzeBackground)
1112 // Get backgrounds in bins of eta
1113 for(Int_t i = 0; i<fBackgroundEtaBins; i++)
1115 // scheme: etaMin = RangeMin + l*binN; etaMax = RangeMin + l*(binN+1)
1117 Double_t etaMin = -(fTrackEtaWindow-fRandConeRadius) + 2*(fTrackEtaWindow-fRandConeRadius)/fBackgroundEtaBins * i;
1118 Double_t etaMax = -(fTrackEtaWindow-fRandConeRadius) + 2*(fTrackEtaWindow-fRandConeRadius)/fBackgroundEtaBins * (i+1);
1119 GetRCBackgroundDensity (fNumberExcludedJets, rcBackgroundRhoMean[i],rcBackgroundRhoMedian[i], etaMin, etaMax);
1120 GetTrackBackgroundDensity (fNumberExcludedJets, trackBackgroundRhoMean[i], trackBackgroundArea[i], etaMin, etaMax);
1121 GetKTBackgroundDensity (fNumberExcludedJets, ktBackgroundRhoMedian[i], ktBackgroundAreaMean[i], etaMin, etaMax);
1122 GetKTBackground2Density (fNumberExcludedJets, ktBackground2RhoMedian[i], ktBackground2AreaMean[i], etaMin, etaMax);
1127 // All eta in one bin
1128 tmpNRCs = GetRCBackgroundDensity (fNumberExcludedJets, rcBackgroundRhoMean[fBackgroundEtaBins], rcBackgroundRhoMedian[fBackgroundEtaBins]);
1129 FillHistogram("hAccConesInRCBackground", tmpNRCs);
1130 GetTrackBackgroundDensity (fNumberExcludedJets, trackBackgroundRhoMean[fBackgroundEtaBins], trackBackgroundArea[fBackgroundEtaBins]);
1131 GetKTBackgroundDensity (fNumberExcludedJets, ktBackgroundRhoMedian[fBackgroundEtaBins], ktBackgroundAreaMean[fBackgroundEtaBins]);
1132 GetKTBackground2Density (fNumberExcludedJets, ktBackground2RhoMedian[fBackgroundEtaBins], ktBackground2AreaMean[fBackgroundEtaBins]);
1137 AliInfo("Calculate()::Centrality&SignalJets&Background-Calculation done.");
1139 ////////////////////// NOTE: Pythia histograms
1142 FillHistogram("hPythiaPtHard", GetPtHard());
1143 FillHistogram("hPythiaNTrials", GetPtHardBin()-0.1, fTrials);
1144 FillHistogram("hPythiaXSection", GetPtHardBin()-0.1, fCrossSection);
1147 AliInfo("Calculate()::Pythia done.");
1151 ////////////////////// NOTE: Track & QA histograms
1155 FillHistogram("hVertexZ",event->GetPrimaryVertex()->GetZ());
1156 FillHistogram("hVertexR",TMath::Sqrt(event->GetPrimaryVertex()->GetX()*event->GetPrimaryVertex()->GetX() + event->GetPrimaryVertex()->GetY()*event->GetPrimaryVertex()->GetY()));
1157 FillHistogram("hCentrality",centralityPercentile);
1159 Int_t trackCountAcc = 0;
1160 Int_t nTracks = fTrackArray->GetEntries();
1161 for (Int_t i = 0; i < nTracks; i++)
1163 AliVTrack* track = static_cast<AliVTrack*>(fTrackArray->At(i));
1164 if (IsTrackInAcceptance(track))
1166 FillHistogram("hTrackPhiEta", track->Phi(),track->Eta(), 1);
1167 FillHistogram("hTrackPt", track->Pt());
1168 FillHistogram("hTrackEta", track->Eta());
1169 FillHistogram("hTrackCharge", track->Charge());
1173 FillHistogram("hTrackCountAcc", trackCountAcc, centralityPercentile);
1177 Int_t clusterCountAcc = 0;
1178 Int_t nClusters = fClusterArray->GetEntries();
1179 for (Int_t i = 0; i < nClusters; i++)
1181 AliVCluster* cluster = static_cast<AliVCluster*>(fClusterArray->At(i));
1182 if (IsClusterInAcceptance(cluster))
1184 FillHistogram("hClusterE", cluster->E());
1188 FillHistogram("hClusterCountAcc", clusterCountAcc, centralityPercentile);
1192 AliInfo("Calculate()::QA done.");
1195 ////////////////////// NOTE: Jet analysis and calculations
1199 FillHistogram("hJetCountAll", fJetArray->GetEntries());
1200 FillHistogram("hJetCountAccepted", fNumberSignalJets);
1201 if (fFirstLeadingJet)
1202 FillHistogram("hLeadingJetPt", fFirstLeadingJet->Pt());
1203 if (fSecondLeadingJet)
1204 FillHistogram("hSecondLeadingJetPt", fSecondLeadingJet->Pt());
1207 if(fNumberSignalJets >= 2)
1209 FillHistogram("hLeadingJetDeltaPhi", GetDeltaPhi(fFirstLeadingJet->Phi(), fSecondLeadingJet->Phi()));
1210 FillHistogram("hLeadingJetDeltaPhiPt", GetDeltaPhi(fFirstLeadingJet->Phi(), fSecondLeadingJet->Phi()), fFirstLeadingJet->Pt());
1212 if (IsDijet(fFirstLeadingJet, fSecondLeadingJet)) // Gettin' the money
1214 FillHistogram("hDijetConstituentsPt", fFirstLeadingJet->Pt()); FillHistogram("hDijetConstituentsPt", fSecondLeadingJet->Pt());
1215 FillHistogram("hDijetLeadingJetPt", fFirstLeadingJet->Pt());
1216 FillHistogram("hDijetPtCorrelation", fFirstLeadingJet->Pt(), fSecondLeadingJet->Pt());
1217 Double_t dummyArea = 0;
1218 GetTrackBackgroundDensity (2, dijetBackground, dummyArea, fFirstLeadingJet, fSecondLeadingJet, kFALSE);
1219 GetTrackBackgroundDensity (2, dijetBackgroundPerpendicular, dummyArea, fFirstLeadingJet, fSecondLeadingJet, kTRUE);
1223 // SIGNAL JET ANALYSIS
1224 for (Int_t i = 0; i<fNumberSignalJets; i++)
1226 AliEmcalJet* tmpJet = fSignalJets[i];
1228 FillHistogram("hJetPtArea", tmpJet->Pt(), tmpJet->Area());
1229 FillHistogram("hJetPtEta", tmpJet->Pt(), tmpJet->Eta());
1230 FillHistogram("hJetPtPhi", tmpJet->Pt(), tmpJet->Phi());
1231 FillHistogram("hJetPtCentrality", tmpJet->Pt(), centralityPercentile);
1232 FillHistogram("hJetArea", tmpJet->Area());
1233 FillHistogram("hJetPt", tmpJet->Pt());
1234 FillHistogram("hJetPhiEta", tmpJet->Phi(),tmpJet->Eta());
1236 // Background subtracted spectra
1238 FillHistogram("hJetPtBgrdSubtractedRC", GetCorrectedJetPt(tmpJet, rcBackgroundRhoMean[fBackgroundEtaBins]));
1239 FillHistogram("hJetPtBgrdSubtractedKT", GetCorrectedJetPt(tmpJet, ktBackgroundRhoMedian[fBackgroundEtaBins], kTRUE));
1240 FillHistogram("hJetPtBgrdSubtractedKT2", GetCorrectedJetPt(tmpJet, ktBackground2RhoMedian[fBackgroundEtaBins], kTRUE));
1241 FillHistogram("hJetPtBgrdSubtractedKTNoEtaCorr", GetCorrectedJetPt(tmpJet, ktBackgroundRhoMedian[fBackgroundEtaBins]));
1242 FillHistogram("hJetPtBgrdSubtractedKT2NoEtaCorr", GetCorrectedJetPt(tmpJet, ktBackground2RhoMedian[fBackgroundEtaBins]));
1243 FillHistogram("hJetPtBgrdSubtractedTR", GetCorrectedJetPt(tmpJet, trackBackgroundRhoMean[fBackgroundEtaBins]));
1246 Double_t tmpCorrFactor = GetJetBackgroundCorrFactor(tmpJet->Eta(), ktBackgroundRhoMedian[fBackgroundEtaBins]);
1247 FillHistogram("hAppliedEtaCorrectionFactor", tmpCorrFactor);
1248 tmpCorrFactor = GetJetBackgroundCorrFactor(tmpJet->Eta(), ktBackground2RhoMedian[fBackgroundEtaBins]);
1249 FillHistogram("hAppliedEtaCorrectionFactor2", tmpCorrFactor);
1251 // Signal jet vs. signal jet
1252 for (Int_t j = i+1; j<fNumberSignalJets; j++)
1254 AliEmcalJet* tmpJet2 = fSignalJets[j];
1255 FillHistogram("hJetDeltaPhi", GetDeltaPhi(tmpJet->Phi(), tmpJet2->Phi()));
1256 FillHistogram("hJetDeltaPhiPt", GetDeltaPhi(tmpJet->Phi(), tmpJet2->Phi()), max(tmpJet->Pt(), tmpJet2->Pt()));
1258 // ### Dijets II ###
1259 if (IsDijet(tmpJet, tmpJet2)) // Gettin' the money
1261 FillHistogram("hDijet2ConstituentsPt", tmpJet->Pt()); FillHistogram("hDijet2ConstituentsPt", tmpJet2->Pt());
1262 FillHistogram("hDijet2LeadingJetPt", fFirstLeadingJet->Pt());
1263 FillHistogram("hDijet2PtCorrelation", tmpJet->Pt(), tmpJet2->Pt());
1267 } //endif AnalyzeJets
1270 AliInfo("Calculate()::Jets done.");
1272 ////////////////////// NOTE: Background analysis
1274 if (fAnalyzeBackground)
1277 Int_t leadingJetIds[] = {-1, -1};
1278 GetLeadingJets(fBackgroundJetArray, &leadingJetIds[0], kFALSE);
1280 for (Int_t i = 0; i < fBackgroundJetArray->GetEntries(); i++)
1282 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fBackgroundJetArray->At(i));
1285 AliError(Form("%s: Could not receive kt jet %d", GetName(), i));
1288 if (!IsBackgroundJetInAcceptance(jet))
1290 if (!((jet->Eta() >= -fBackgroundJetEtaWindow) && (jet->Eta() < fBackgroundJetEtaWindow)))
1293 FillHistogram("hKTJetPhiEta", jet->Phi(),jet->Eta());
1294 if(i==leadingJetIds[0])
1295 FillHistogram("hKTLeadingJetPhiEta", jet->Phi(),jet->Eta());
1299 // ############# RC, Track, and KT background calculations
1300 Double_t etaMin = 0;
1301 for (Int_t i=0;i<fBackgroundEtaBins;i++)
1303 etaMin = -(fTrackEtaWindow-fRandConeRadius) + 2*(fTrackEtaWindow-fRandConeRadius)/fBackgroundEtaBins * (i+0.5);
1304 FillHistogram("hRCBackground", etaMin, rcBackgroundRhoMean[i]);
1305 FillHistogram("hTrackBackground", etaMin, trackBackgroundRhoMean[i]);
1306 FillHistogram("hKTBackground", etaMin, ktBackgroundRhoMedian[i]);
1307 FillHistogram("hKTBackground2", etaMin, ktBackground2RhoMedian[i]);
1308 if(centralityPercentile <= 20.)
1310 FillHistogram("hRCBackgroundMostCentral", etaMin, rcBackgroundRhoMean[i]);
1311 FillHistogram("hTrackBackgroundMostCentral", etaMin, trackBackgroundRhoMean[i]);
1312 FillHistogram("hKTBackgroundMostCentral", etaMin, ktBackgroundRhoMedian[i]);
1313 FillHistogram("hKTBackground2MostCentral", etaMin, ktBackground2RhoMedian[i]);
1315 else if(centralityPercentile >= 80.)
1317 FillHistogram("hRCBackgroundMostPeripheral", etaMin, rcBackgroundRhoMean[i]);
1318 FillHistogram("hTrackBackgroundMostPeripheral", etaMin, trackBackgroundRhoMean[i]);
1319 FillHistogram("hKTBackgroundMostPeripheral", etaMin, ktBackgroundRhoMedian[i]);
1320 FillHistogram("hKTBackground2MostPeripheral", etaMin, ktBackground2RhoMedian[i]);
1324 FillHistogram("hRCBackgroundVsCentrality", rcBackgroundRhoMean[fBackgroundEtaBins], centralityPercentile);
1325 FillHistogram("hTrackBackgroundVsCentrality", trackBackgroundRhoMean[fBackgroundEtaBins], centralityPercentile);
1326 FillHistogram("hKTBackgroundVsCentrality", ktBackgroundRhoMedian[fBackgroundEtaBins], centralityPercentile);
1327 FillHistogram("hKTBackground2VsCentrality", ktBackground2RhoMedian[fBackgroundEtaBins], centralityPercentile);
1329 if (dijetBackground >= 0)
1331 // Background in Dijet events
1332 FillHistogram("hDijetBackground", dijetBackground);
1333 if(centralityPercentile <= 20.)
1334 FillHistogram("hDijetBackgroundMostCentral", dijetBackground);
1335 FillHistogram("hDijetBackgroundVsCentrality", dijetBackground, centralityPercentile);
1337 if (dijetBackgroundPerpendicular >= 0)
1339 // Background in Dijet events
1340 FillHistogram("hDijetBackgroundPerpendicular", dijetBackgroundPerpendicular);
1341 if(centralityPercentile <= 20.)
1342 FillHistogram("hDijetBackgroundPerpendicularMostCentral", dijetBackgroundPerpendicular);
1343 FillHistogram("hDijetBackgroundPerpendicularVsCentrality", dijetBackgroundPerpendicular, centralityPercentile);
1346 // ########## Delta pT calculations (most central, kt is eta corrected)
1347 if (centralityPercentile <= 20.)
1349 Double_t tmpDeltaPtKT, tmpDeltaPtKT2Excl, tmpDeltaPtKT1Excl;
1350 Double_t tmpDeltaPtKTEta, tmpDeltaPtKTEta2Excl, tmpDeltaPtKTEta1Excl, tmpDeltaPtKT2Eta2Excl;
1351 Double_t tmpDeltaPtRC, tmpDeltaPtRC2Excl, tmpDeltaPtRC1Excl;
1352 Double_t tmpDeltaPtTR, tmpDeltaPtTR2Excl, tmpDeltaPtTR1Excl;
1354 GetDeltaPt(tmpDeltaPtKT, ktBackgroundRhoMedian[fBackgroundEtaBins], 0, -1, kTRUE);
1355 GetDeltaPt(tmpDeltaPtKTEta, ktBackgroundRhoMedian[fKTDeltaPtEtaBin], 0, fKTDeltaPtEtaBin);
1356 GetDeltaPt(tmpDeltaPtRC, rcBackgroundRhoMean[fBackgroundEtaBins], 0);
1357 GetDeltaPt(tmpDeltaPtTR, trackBackgroundRhoMean[fBackgroundEtaBins], 0);
1359 GetDeltaPt(tmpDeltaPtKT1Excl, ktBackgroundRhoMedian[fBackgroundEtaBins], 1, -1, kTRUE);
1360 GetDeltaPt(tmpDeltaPtKTEta1Excl, ktBackgroundRhoMedian[fKTDeltaPtEtaBin], 1, fKTDeltaPtEtaBin);
1361 GetDeltaPt(tmpDeltaPtRC1Excl, rcBackgroundRhoMean[fBackgroundEtaBins], 1);
1362 GetDeltaPt(tmpDeltaPtTR1Excl, trackBackgroundRhoMean[fBackgroundEtaBins], 1);
1364 GetDeltaPt(tmpDeltaPtKT2Excl, ktBackgroundRhoMedian[fBackgroundEtaBins], 2, -1, kTRUE);
1365 GetDeltaPt(tmpDeltaPtKTEta2Excl, ktBackgroundRhoMedian[fKTDeltaPtEtaBin], 2, fKTDeltaPtEtaBin);
1366 GetDeltaPt(tmpDeltaPtRC2Excl, rcBackgroundRhoMean[fBackgroundEtaBins], 2);
1367 GetDeltaPt(tmpDeltaPtTR2Excl, trackBackgroundRhoMean[fBackgroundEtaBins], 2);
1369 GetDeltaPt(tmpDeltaPtKT2Eta2Excl, ktBackground2RhoMedian[fKTDeltaPtEtaBin], 2, fKTDeltaPtEtaBin);
1372 if(tmpDeltaPtKT > -10000.0)
1373 FillHistogram("hDeltaPtKT", tmpDeltaPtKT);
1374 if(tmpDeltaPtKT1Excl > -10000.0)
1375 FillHistogram("hDeltaPtKT1Excl", tmpDeltaPtKT1Excl);
1376 if(tmpDeltaPtKT2Excl > -10000.0)
1377 FillHistogram("hDeltaPtKT2Excl", tmpDeltaPtKT2Excl);
1379 if(tmpDeltaPtKT > -10000.0)
1380 FillHistogram("hDeltaPtKTEta", tmpDeltaPtKTEta);
1381 if(tmpDeltaPtKTEta1Excl > -10000.0)
1382 FillHistogram("hDeltaPtKTEta1Excl", tmpDeltaPtKTEta1Excl);
1383 if(tmpDeltaPtKTEta2Excl > -10000.0)
1384 FillHistogram("hDeltaPtKTEta2Excl", tmpDeltaPtKTEta2Excl);
1385 if(tmpDeltaPtKT2Eta2Excl > -10000.0)
1386 FillHistogram("hDeltaPtKT2Eta2Excl", tmpDeltaPtKT2Eta2Excl);
1389 if(tmpDeltaPtRC > -10000.0)
1390 FillHistogram("hDeltaPtRC", tmpDeltaPtRC);
1391 if(tmpDeltaPtRC1Excl > -10000.0)
1392 FillHistogram("hDeltaPtRC1Excl", tmpDeltaPtRC1Excl);
1393 if(tmpDeltaPtRC2Excl > -10000.0)
1394 FillHistogram("hDeltaPtRC2Excl", tmpDeltaPtRC2Excl);
1396 if(tmpDeltaPtTR > -10000.0)
1397 FillHistogram("hDeltaPtTR", tmpDeltaPtTR);
1398 if(tmpDeltaPtTR1Excl > -10000.0)
1399 FillHistogram("hDeltaPtTR1Excl", tmpDeltaPtTR1Excl);
1400 if(tmpDeltaPtTR2Excl > -10000.0)
1401 FillHistogram("hDeltaPtTR2Excl", tmpDeltaPtTR2Excl);
1406 AliInfo("Calculate()::Background done.");
1410 //________________________________________________________________________
1411 Bool_t AliAnalysisTaskChargedJetsPA::Notify()
1413 // Implemented Notify() to read the cross sections
1414 // and number of trials from pyxsec.root
1417 AliInfo("Notify started.");
1422 TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
1423 TFile *currFile = tree->GetCurrentFile();
1425 TString file(currFile->GetName());
1427 if(file.Contains("root_archive.zip#")){
1428 Ssiz_t pos1 = file.Index("root_archive",12,TString::kExact);
1429 Ssiz_t pos = file.Index("#",1,pos1,TString::kExact);
1430 file.Replace(pos+1,20,"");
1433 // not an archive take the basename....
1434 file.ReplaceAll(gSystem->BaseName(file.Data()),"");
1437 TFile *fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec.root")); // problem that we cannot really test the existance of a file in a archive so we have to lvie with open error message from root
1439 // next trial fetch the histgram file
1440 fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec_hists.root"));
1442 // not a severe condition but inciate that we have no information
1446 // find the tlist we want to be independtent of the name so use the Tkey
1447 TKey* key = (TKey*)fxsec->GetListOfKeys()->At(0);
1452 TList *list = dynamic_cast<TList*>(key->ReadObj());
1457 fCrossSection = ((TProfile*)list->FindObject("h1Xsec"))->GetBinContent(1);
1458 fTrials = ((TH1F*)list->FindObject("h1Trials"))->GetBinContent(1);
1461 } // no tree pyxsec.root
1463 TTree *xtree = (TTree*)fxsec->Get("Xsection");
1469 Double_t xsection = 0;
1470 xtree->SetBranchAddress("xsection",&xsection);
1471 xtree->SetBranchAddress("ntrials",&ntrials);
1474 fCrossSection = xsection;
1478 AliInfo("Notify ended.");
1484 //________________________________________________________________________
1485 inline void AliAnalysisTaskChargedJetsPA::FillHistogram(const char * key, Double_t x)
1487 TH1* tmpHist = static_cast<TH1*>(fOutputList->FindObject(GetHistoName(key)));
1490 AliInfo(Form("Cannot find histogram <%s> ",key)) ;
1497 //________________________________________________________________________
1498 inline void AliAnalysisTaskChargedJetsPA::FillHistogram(const char * key, Double_t x, Double_t y)
1500 TH1* tmpHist = static_cast<TH1*>(fOutputList->FindObject(GetHistoName(key)));
1503 AliInfo(Form("Cannot find histogram <%s> ",key));
1507 if (tmpHist->IsA()->GetBaseClass("TH1"))
1508 static_cast<TH1*>(tmpHist)->Fill(x,y); // Fill x with y
1509 else if (tmpHist->IsA()->GetBaseClass("TH2"))
1510 static_cast<TH2*>(tmpHist)->Fill(x,y); // Fill x,y with 1
1513 //________________________________________________________________________
1514 inline void AliAnalysisTaskChargedJetsPA::FillHistogram(const char * key, Double_t x, Double_t y, Double_t add)
1516 TH2* tmpHist = static_cast<TH2*>(fOutputList->FindObject(GetHistoName(key)));
1519 AliInfo(Form("Cannot find histogram <%s> ",key));
1523 tmpHist->Fill(x,y,add);
1525 //________________________________________________________________________
1526 template <class T> T* AliAnalysisTaskChargedJetsPA::AddHistogram1D(const char* name, const char* title, const char* options, Int_t xBins, Double_t xMin, Double_t xMax, const char* xTitle, const char* yTitle)
1528 T* tmpHist = new T(GetHistoName(name), GetHistoName(title), xBins, xMin, xMax);
1530 tmpHist->GetXaxis()->SetTitle(xTitle);
1531 tmpHist->GetYaxis()->SetTitle(yTitle);
1532 tmpHist->SetOption(options);
1533 tmpHist->SetMarkerStyle(kFullCircle);
1536 fHistList->Add(tmpHist);
1542 //________________________________________________________________________
1543 template <class T> T* AliAnalysisTaskChargedJetsPA::AddHistogram2D(const char* name, const char* title, const char* options, Int_t xBins, Double_t xMin, Double_t xMax, Int_t yBins, Double_t yMin, Double_t yMax, const char* xTitle, const char* yTitle, const char* zTitle)
1545 T* tmpHist = new T(GetHistoName(name), GetHistoName(title), xBins, xMin, xMax, yBins, yMin, yMax);
1546 tmpHist->GetXaxis()->SetTitle(xTitle);
1547 tmpHist->GetYaxis()->SetTitle(yTitle);
1548 tmpHist->GetZaxis()->SetTitle(zTitle);
1549 tmpHist->SetOption(options);
1550 tmpHist->SetMarkerStyle(kFullCircle);
1553 fHistList->Add(tmpHist);
1559 //________________________________________________________________________
1560 void AliAnalysisTaskChargedJetsPA::Terminate(Option_t *)
1562 PostData(1, fOutputList);
1565 fOutputList = dynamic_cast<TList*> (GetOutputData(1)); // '1' refers to the output slot
1567 printf("ERROR: Output list not available\n");
1572 //________________________________________________________________________
1573 AliAnalysisTaskChargedJetsPA::~AliAnalysisTaskChargedJetsPA()
1575 // Destructor. Clean-up the output list, but not the histograms that are put inside
1576 // (the list is owner and will clean-up these histograms). Protect in PROOF case.
1577 if (fOutputList && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
1582 //________________________________________________________________________
1583 void AliAnalysisTaskChargedJetsPA::UserCreateOutputObjects()
1585 // called once to create user defined output objects like histograms, plots etc.
1586 // and to put it on the output list.
1587 // Note: Saving to file with e.g. OpenFile(0) is must be before creating other objects.
1589 fRandom = new TRandom3(0);
1591 fOutputList = new TList();
1592 fOutputList->SetOwner(); // otherwise it produces leaks in merging
1594 PostData(1, fOutputList);
1597 //________________________________________________________________________
1598 void AliAnalysisTaskChargedJetsPA::UserExec(Option_t *)
1601 AliInfo("UserExec() started.");
1604 Calculate(InputEvent());
1606 PostData(1, fOutputList);