]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskChargedJetsPA.cxx
01b78154c45dc588e8825553a0f09b23dfc29138
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / UserTasks / AliAnalysisTaskChargedJetsPA.cxx
1 #include "AliAnalysisTaskChargedJetsPA.h"
2
3
4 //TODO: FillHistogram can be done better with virtual TH1(?)
5 ClassImp(AliAnalysisTaskChargedJetsPA)
6
7 // ######################################################################################## DEFINE HISTOGRAMS
8 void AliAnalysisTaskChargedJetsPA::Init()
9 {
10   #ifdef DEBUGMODE
11     AliInfo("Creating histograms.");
12   #endif
13   // NOTE: Track & Cluster & QA histograms
14   if (fAnalyzeQA)
15   {
16   
17     AddHistogram1D<TH1D>("hNumberEvents", "Number of events (0 = before, 1 = after vertex cuts)", "", 2, 0, 2, "#Delta z(cm)","N^{Events}/cut");
18       
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");
22
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}");
25
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");
31
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}");
34   }
35
36   // NOTE: Pythia histograms
37   if (fAnalyzePythia)
38   {
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");
42   }
43
44   // NOTE: Jet histograms
45   if (fAnalyzeJets)
46   {
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}");
54
55     AddHistogram1D<TH1D>("hJetPtBgrdSubtractedTR", "Jets p_{T} distribution, Track background subtracted", "", 500, -50., 200., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
56
57
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}");
70
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}");
73
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}");
77
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}");
81
82   }
83   // NOTE: Jet background histograms
84   if (fAnalyzeBackground)
85   {
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}");
90
91     Double_t dptEtaMin = -(fTrackEtaWindow-fRandConeRadius) + 2*(fTrackEtaWindow-fRandConeRadius)/fBackgroundEtaBins *  fKTDeltaPtEtaBin;
92     Double_t dptEtaMax = -(fTrackEtaWindow-fRandConeRadius) + 2*(fTrackEtaWindow-fRandConeRadius)/fBackgroundEtaBins * (fKTDeltaPtEtaBin+1);
93
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}");
97
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}");
99
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}");
103
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}");
107
108
109
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)");
112
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");
116
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");
120   
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");
122
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}");
124
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");
128
129     
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");
134   
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");
139
140
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");
145
146   }
147
148
149
150   // register Histograms
151   for (Int_t i = 0; i < fHistCount; i++)
152   {
153     fOutputList->Add(fHistList->At(i));
154   }
155   
156   PostData(1,fOutputList); // important for merging
157
158 }
159
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)
162 {
163   #ifdef DEBUGMODE
164     AliInfo("Calling constructor.");
165   #endif
166
167   // 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
171   // Constructor
172
173   // Every instance of this task gets his own number
174   static Int_t instance = 0;
175   fTaskInstanceCounter = instance;
176   instance++;
177
178   fTrackArrayName = new TString(trackArrayName);
179   fClusterArrayName = new TString(clusterArrayName);
180   if (strcmp(fTrackArrayName->Data(),"") == 0)
181     fAnalyzeQA = kFALSE;
182   else
183   {
184     fAnalyzeQA = kTRUE;
185     if (fTrackArrayName->Contains("MCParticles")) //TODO: Hardcoded for now
186       fIsMC = kTRUE;
187   }
188
189   fJetArrayName = new TString(jetArrayName);
190   if (strcmp(fJetArrayName->Data(),"") == 0)
191     fAnalyzeJets = kFALSE;
192   else
193     fAnalyzeJets = kTRUE;
194     
195   fBackgroundJetArrayName = new TString(backgroundJetArrayName);
196   if (strcmp(fBackgroundJetArrayName->Data(),"") == 0)
197     fAnalyzeBackground = kFALSE;
198   else
199     fAnalyzeBackground = kTRUE;
200
201   DefineOutput(1, TList::Class());
202  
203   fHistList = new TList();
204
205   #ifdef DEBUGMODE
206     AliInfo("Constructor done.");
207   #endif
208   
209 }
210
211 //________________________________________________________________________
212 inline Double_t AliAnalysisTaskChargedJetsPA::GetConePt(Double_t eta, Double_t phi, Double_t radius)
213 {
214   Double_t tmpConePt = 0.0;
215
216   for (Int_t i = 0; i < fTrackArray->GetEntries(); i++)
217   {
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();
222   }
223   return tmpConePt;
224 }
225
226 //________________________________________________________________________
227 inline Double_t AliAnalysisTaskChargedJetsPA::GetPtHard()
228 {
229   Double_t tmpPtHard = -1.0;
230
231   if (!MCEvent())
232     AliError("MCEvent not accessible although demanded!");
233   else
234   {
235     AliGenPythiaEventHeader* pythiaHeader = dynamic_cast<AliGenPythiaEventHeader*>(MCEvent()->GenEventHeader());
236     if (!pythiaHeader)
237       AliError("Pythia Header not accessible!");
238     else
239       tmpPtHard = pythiaHeader->GetPtHard();
240   }
241   return tmpPtHard;
242 }
243
244 //________________________________________________________________________
245 inline Int_t AliAnalysisTaskChargedJetsPA::GetPtHardBin()
246 {
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};
250
251   Int_t tmpPtHardBin = 0;
252   Double_t tmpPtHard = GetPtHard();
253  
254   for (tmpPtHardBin = 0; tmpPtHardBin <= fNumPtHardBins; tmpPtHardBin++)
255     if (tmpPtHard >= localkPtHardLowerEdges[tmpPtHardBin] && tmpPtHard < localkPtHardHigherEdges[tmpPtHardBin])
256       break;
257
258   return tmpPtHardBin;
259 }
260
261
262 //________________________________________________________________________
263 inline Bool_t AliAnalysisTaskChargedJetsPA::IsTrackInCone(AliVTrack* track, Double_t eta, Double_t phi, Double_t radius)
264 {
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();
271   else
272     trackPhi = track->Phi();
273   
274   if ( TMath::Abs(trackPhi-phi)*TMath::Abs(trackPhi-phi) + TMath::Abs(track->Eta()-eta)*TMath::Abs(track->Eta()-eta) <= radius*radius)
275     return kTRUE;
276   
277   return kFALSE;
278 }
279
280 //________________________________________________________________________
281 inline Bool_t AliAnalysisTaskChargedJetsPA::IsTrackInAcceptance(AliVParticle* track)
282 {
283   if (track != 0)
284     if (TMath::Abs(track->Eta()) <= fTrackEtaWindow)
285       if (track->Pt() >= fMinTrackPt)
286         return kTRUE;
287
288   return kFALSE;
289 }
290
291 //________________________________________________________________________
292 inline Bool_t AliAnalysisTaskChargedJetsPA::IsClusterInAcceptance(AliVCluster* cluster)
293 {
294   if (cluster != 0)
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)
299             return kTRUE;
300
301   return kFALSE;
302 }
303
304
305 //________________________________________________________________________
306 inline Bool_t AliAnalysisTaskChargedJetsPA::IsBackgroundJetInAcceptance(AliEmcalJet *jet)
307 {   
308   if (jet != 0)
309     if (TMath::Abs(jet->Eta()) <= fBackgroundJetEtaWindow)
310       if (jet->Pt() >= fMinBackgroundJetPt)
311         return kTRUE;
312
313   return kFALSE;
314 }
315
316 //________________________________________________________________________
317 inline Bool_t AliAnalysisTaskChargedJetsPA::IsSignalJetInAcceptance(AliEmcalJet *jet)
318 {   
319   if (jet != 0)
320     if (TMath::Abs(jet->Eta()) <= fSignalJetEtaWindow)
321       if (jet->Pt() >= fMinJetPt)
322         if (jet->Area() >= fMinJetArea)
323           return kTRUE;
324
325   return kFALSE;
326 }
327
328 //________________________________________________________________________
329 inline Bool_t AliAnalysisTaskChargedJetsPA::IsDijet(AliEmcalJet *jet1, AliEmcalJet *jet2)
330 {   
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?
335         return kTRUE;
336
337   return kFALSE;
338 }
339
340 //________________________________________________________________________
341 void AliAnalysisTaskChargedJetsPA::ExecOnce()
342 {
343   #ifdef DEBUGMODE
344     AliInfo("Starting ExecOnce.");
345   #endif
346   fInitialized = kTRUE;
347
348   // Check for track array
349   if (strcmp(fTrackArrayName->Data(), "") != 0)
350   {
351     fTrackArray = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fTrackArrayName->Data()));
352     fHasTracks = kTRUE;
353     if (!fTrackArray) 
354     {
355       AliInfo(Form("%s: Could not retrieve tracks %s! This is OK, if tracks are not demanded.", GetName(), fTrackArrayName->Data())); 
356       fHasTracks = kFALSE;
357     } 
358     else
359     {
360       TClass *cl = fTrackArray->GetClass();
361       if (!cl->GetBaseClass("AliVParticle"))
362       {
363         AliError(Form("%s: Collection %s does not contain AliVParticle objects!", GetName(), fTrackArrayName->Data())); 
364         fTrackArray = 0;
365         fHasTracks = kFALSE;
366       }
367     }
368   }
369   // Check for cluster array
370   if (strcmp(fClusterArrayName->Data(), "") != 0)
371   {
372     fClusterArray = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fClusterArrayName->Data()));
373     fHasClusters = kTRUE;
374     if (!fClusterArray) 
375     {
376       AliInfo(Form("%s: Could not retrieve clusters %s! This is OK, if clusters are not demanded.", GetName(), fClusterArrayName->Data())); 
377       fHasClusters = kFALSE;
378     } 
379     else
380     {
381       TClass *cl = fClusterArray->GetClass();
382       if (!cl->GetBaseClass("AliVCluster"))
383       {
384         AliError(Form("%s: Collection %s does not contain AliVCluster objects!", GetName(), fClusterArrayName->Data())); 
385         fClusterArray = 0;
386         fHasClusters = kFALSE;
387       }
388     }
389   }
390
391   // Check for jet array
392   if (strcmp(fJetArrayName->Data(), "") != 0)
393   {
394     fJetArray = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fJetArrayName->Data()));
395     fHasJets = kTRUE;
396
397     if (!fJetArray) 
398     {
399       AliInfo(Form("%s: Could not retrieve jets %s! This is OK, if jets are not demanded.", GetName(), fJetArrayName->Data())); 
400       fHasJets = kFALSE;
401     } 
402     else
403     {
404       if (!fJetArray->GetClass()->GetBaseClass("AliEmcalJet")) 
405       {
406         AliError(Form("%s: Collection %s does not contain AliEmcalJet objects!", GetName(), fJetArrayName->Data())); 
407         fJetArray = 0;
408         fHasJets = kFALSE;
409       }
410     }
411   }
412
413   // Check for background object
414   if (strcmp(fBackgroundJetArrayName->Data(), "") != 0)
415   {
416     fBackgroundJetArray = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fBackgroundJetArrayName->Data()));
417     fHasBackgroundJets = kTRUE;
418     if (!fBackgroundJetArray)
419     {
420       AliInfo(Form("%s: Could not retrieve background jets %s! This is OK, if background is not demanded.", GetName(), fBackgroundJetArrayName->Data())); 
421       fHasBackgroundJets = kFALSE;
422     }
423   }
424
425   // Look, if initialization is OK
426   if ((!fHasTracks && fAnalyzeQA) || (!fHasTracks && fAnalyzeBackground))
427   {
428     AliError(Form("%s: Tracks NOT successfully casted although demanded! Deactivating QA and background analysis",GetName()));
429     fAnalyzeQA = kFALSE;
430     fAnalyzeBackground = kFALSE;
431   }
432   if ((!fHasJets && fAnalyzeJets) || (!fHasJets && fAnalyzeBackground))
433   {
434     AliError(Form("%s: Jets NOT successfully casted although demanded!  Deactivating jet- and background analysis",GetName()));
435     fAnalyzeJets = kFALSE;
436     fAnalyzeBackground = kFALSE;
437   }
438   if (!fHasBackgroundJets && fAnalyzeBackground)
439   {
440     AliError(Form("%s: Background NOT successfully casted although demanded!  Deactivating background analysis",GetName()));
441     fAnalyzeBackground = kFALSE;
442   }
443
444   // Histogram init
445   Init();
446
447   #ifdef DEBUGMODE
448     AliInfo("ExecOnce done.");
449   #endif
450
451 }
452
453 //________________________________________________________________________
454 void AliAnalysisTaskChargedJetsPA::GetSignalJets()
455 {
456   // Reset vars
457   fFirstLeadingJet = NULL;
458   fSecondLeadingJet = NULL;
459   fNumberSignalJets = 0;
460
461   Float_t maxJetPts[] = {0, 0};
462   Int_t jetIDArray[]   = {-1, -1};
463   Int_t jetCount = fJetArray->GetEntries();
464
465   // Go through all jets and save signal jets and the two leading ones
466   for (Int_t i = 0; i < jetCount; i++)
467   {
468     AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJetArray->At(i));
469     if (!jet)
470     {
471       AliError(Form("%s: Could not receive jet %d", GetName(), i));
472       continue;
473     }
474
475     if (!IsSignalJetInAcceptance(jet)) continue;
476     
477     if (jet->Pt() > maxJetPts[0]) 
478     {
479       maxJetPts[1] = maxJetPts[0];
480       jetIDArray[1] = jetIDArray[0];
481       maxJetPts[0] = jet->Pt();
482       jetIDArray[0] = i;
483     }
484     else if (jet->Pt() > maxJetPts[1]) 
485     {
486       maxJetPts[1] = jet->Pt();
487       jetIDArray[1] = i;
488     }
489     fSignalJets[fNumberSignalJets] = jet;
490     fNumberSignalJets++;
491   }
492   
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]));
497
498 }
499
500 //________________________________________________________________________
501 Int_t AliAnalysisTaskChargedJetsPA::GetLeadingJets(TClonesArray* jetArray, Int_t* jetIDArray, Bool_t isSignalJets)
502 {
503 // Writes first two leading jets into already registered array jetIDArray
504
505   if (!jetArray)
506   {
507     AliError("Could not get the jet array to get leading jets from it!");
508     return 0;
509   }
510
511   Float_t maxJetPts[] = {0, 0};
512   jetIDArray[0] = -1;
513   jetIDArray[1] = -1;
514
515   Int_t jetCount = jetArray->GetEntries();
516   Int_t jetCountAccepted = 0;
517
518   for (Int_t i = 0; i < jetCount; i++)
519   {
520     AliEmcalJet* jet = static_cast<AliEmcalJet*>(jetArray->At(i));
521     if (!jet) 
522     {
523       AliError(Form("%s: Could not receive jet %d", GetName(), i));
524       continue;
525     }
526
527     if(isSignalJets)
528     {
529       if (!IsSignalJetInAcceptance(jet)) continue;
530     }
531     else
532     {
533       if (!IsBackgroundJetInAcceptance(jet)) continue;
534     }    
535
536     if (jet->Pt() > maxJetPts[0]) 
537     {
538       maxJetPts[1] = maxJetPts[0];
539       jetIDArray[1] = jetIDArray[0];
540       maxJetPts[0] = jet->Pt();
541       jetIDArray[0] = i;
542     }
543     else if (jet->Pt() > maxJetPts[1]) 
544     {
545       maxJetPts[1] = jet->Pt();
546       jetIDArray[1] = i;
547     }
548     jetCountAccepted++;
549   }
550   return jetCountAccepted;
551 }
552 //________________________________________________________________________
553 Double_t AliAnalysisTaskChargedJetsPA::GetJetBackgroundCorrFactor(Double_t eta, Double_t background)
554 {
555   Double_t tmpCorrFactor = 1.0;
556
557   if(fJetBgrdCorrectionFactors)
558     tmpCorrFactor = fJetBgrdCorrectionFactors->GetBinContent
559                           (
560                             fJetBgrdCorrectionFactors->GetXaxis()->FindBin(eta),
561                             fJetBgrdCorrectionFactors->GetYaxis()->FindBin(background)
562                           );
563
564   return tmpCorrFactor;
565 }
566 //________________________________________________________________________
567 Double_t AliAnalysisTaskChargedJetsPA::GetCorrectedJetPt(AliEmcalJet* jet, Double_t background, Bool_t useEtaCorrection)
568 {
569   #ifdef DEBUGMODE
570     AliInfo("Getting corrected jet spectra.");
571   #endif
572
573   if(!jet)
574   {
575     AliError("Jet pointer passed to GetCorrectedJet() not valid!");
576     return -1.0;
577   }
578
579   Double_t correctedPt = -1.0;
580
581   // Get correction factor from saved histo or similar in dependence of jet eta and background density
582   Double_t corrfactor = 1.0;
583   if(useEtaCorrection)
584   {
585     corrfactor = GetJetBackgroundCorrFactor(jet->Eta(), background);
586   }
587
588   // Get Eta corrected background
589   Double_t tmpCorrectedBackground = background * corrfactor;
590
591   // Subtract background
592   correctedPt = jet->Pt() - tmpCorrectedBackground * jet->Area();
593
594   #ifdef DEBUGMODE
595     AliInfo("Got corrected jet spectra.");
596   #endif 
597
598   return correctedPt;
599 }
600
601
602 //________________________________________________________________________
603 void AliAnalysisTaskChargedJetsPA::GetDeltaPt(Double_t& deltaPt, Double_t rho, Int_t numberExcludeLeadingJets, Int_t usedEtaBin, Bool_t useEtaCorrection)
604 {
605   #ifdef DEBUGMODE
606     AliInfo("Getting Delta Pt.");
607   #endif
608
609   // Define the tmp delta pt
610   deltaPt = -10000.0;
611
612   // Exclude UP TO numberExcludeLeadingJets
613   if (fNumberSignalJets < 2)
614     numberExcludeLeadingJets = fNumberSignalJets;
615   
616   Double_t etaMin, etaMax;
617   if (usedEtaBin==-1)
618   {
619     etaMin = -(fTrackEtaWindow-fRandConeRadius);
620     etaMax = +(fTrackEtaWindow-fRandConeRadius);
621   }
622   else
623   {
624     etaMin = -(fTrackEtaWindow-fRandConeRadius) + 2*(fTrackEtaWindow-fRandConeRadius)/fBackgroundEtaBins *  usedEtaBin;
625     etaMax = -(fTrackEtaWindow-fRandConeRadius) + 2*(fTrackEtaWindow-fRandConeRadius)/fBackgroundEtaBins * (usedEtaBin+1);
626   }  
627
628
629   Double_t tmpRandConeEta = 0.0;
630   Double_t tmpRandConePhi = 0.0;
631
632   Bool_t coneValid = kTRUE;
633
634
635   tmpRandConeEta = etaMin + fRandom->Rndm()*(etaMax-etaMin);
636   tmpRandConePhi = fRandom->Rndm()*TMath::TwoPi();
637
638   // Apply eta correction on demand
639   if(useEtaCorrection)
640     rho = GetJetBackgroundCorrFactor(tmpRandConeEta, rho)*rho;
641
642   // Go through all excluded leading jets and check if there's an overlap
643   for(Int_t j=0;j<numberExcludeLeadingJets;j++)
644   {
645     AliEmcalJet* tmpJet = NULL;
646
647     if (j==0)
648       tmpJet = fFirstLeadingJet;
649     else if (j==1)
650       tmpJet = fSecondLeadingJet;
651     else
652       AliFatal("Trying to exclude more than 2 jets for delta pt -- not implemented.");
653
654     Double_t excludedJetPhi = tmpJet->Phi();
655     Double_t excludedJetEta = tmpJet->Eta();
656     Double_t tmpDeltaPhi = GetDeltaPhi(tmpRandConePhi, excludedJetPhi);
657
658     // Check, if cone has overlap with jet
659     if ( tmpDeltaPhi*tmpDeltaPhi + TMath::Abs(tmpRandConeEta-excludedJetEta)*TMath::Abs(tmpRandConeEta-excludedJetEta) <= fRandConeRadius*fRandConeRadius)
660     {
661       coneValid = kFALSE;
662       break;
663     }
664   }
665
666   // Get the cones' pt and calculate delta pt
667   if (coneValid)
668     deltaPt = GetConePt(tmpRandConeEta,tmpRandConePhi,fRandConeRadius) - (rho*fRandConeRadius*fRandConeRadius*TMath::Pi());
669   #ifdef DEBUGMODE
670     AliInfo("Got Delta Pt.");
671   #endif
672 }
673
674 //________________________________________________________________________
675 void AliAnalysisTaskChargedJetsPA::GetKTBackgroundDensity(Int_t numberExcludeLeadingJets, Double_t& rhoMedian, Double_t& areaMean, Double_t etaMin, Double_t etaMax)
676 {
677   #ifdef DEBUGMODE
678     AliInfo("Getting KT background density.");
679   #endif
680
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)
685
686   // Setting invalid values
687   rhoMedian = -1.0;
688   areaMean= -1.0;
689
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))
695   {
696     etaMin = -fBackgroundJetEtaWindow;
697     etaMax = +fBackgroundJetEtaWindow;
698   }
699
700   Int_t jetCountAccepted = 0;
701   Int_t jetCount = fBackgroundJetArray->GetEntries();
702
703   for (Int_t i = 0; i < jetCount; i++)
704   {
705     AliEmcalJet* jet = static_cast<AliEmcalJet*>(fBackgroundJetArray->At(i));
706     if (!jet)
707     {
708       AliError(Form("%s: Could not receive jet %d", GetName(), i));
709       continue;
710     } 
711
712     // exclude leading jets
713     if (numberExcludeLeadingJets > 0)
714       if (i == maxJetIds[0])
715         continue;
716     if (numberExcludeLeadingJets > 1)
717       if (i == maxJetIds[1])
718         continue;
719       
720
721
722     if (!IsBackgroundJetInAcceptance(jet))
723       continue;
724     if (!((jet->Eta() >= etaMin) && (jet->Eta() < etaMax)))
725       continue;
726
727     
728     tmpRhos[jetCountAccepted] = jet->Pt() / jet->Area();
729     tmpAreas[jetCountAccepted] = jet->Area();
730     jetCountAccepted++;
731   }
732
733   if (jetCountAccepted > 0)
734   {
735     rhoMedian = TMath::Median(jetCountAccepted, tmpRhos);
736     areaMean   = TMath::Mean(jetCountAccepted, tmpAreas);
737   }
738   #ifdef DEBUGMODE
739     AliInfo("Got KT background density.");
740   #endif
741 }
742
743 //________________________________________________________________________
744 void AliAnalysisTaskChargedJetsPA::GetKTBackground2Density(Int_t numberExcludeLeadingJets, Double_t& rhoMedian, Double_t& areaMean, Double_t etaMin, Double_t etaMax)
745 {
746   #ifdef DEBUGMODE
747     AliInfo("Getting KT background 2 density.");
748   #endif
749
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];
753
754   // Setting invalid values
755   rhoMedian = -1.0;
756   areaMean= -1.0;
757
758   if ((etaMin == 0) && (etaMax == 0))
759   {
760     etaMin = -fBackgroundJetEtaWindow;
761     etaMax = +fBackgroundJetEtaWindow;
762   }
763
764   Int_t jetCountAccepted = 0;
765   Int_t jetCount = fBackgroundJetArray->GetEntries();
766
767   for (Int_t i = 0; i < jetCount; i++)
768   {
769     Bool_t jetValid = kTRUE;
770     AliEmcalJet* jet = static_cast<AliEmcalJet*>(fBackgroundJetArray->At(i));
771     if (!jet)
772     {
773       AliError(Form("%s: Could not receive jet %d", GetName(), i));
774       continue;
775     } 
776
777     if (!((jet->Eta() >= etaMin) && (jet->Eta() < etaMax)))
778       continue;
779     if (!IsBackgroundJetInAcceptance(jet))
780       continue;
781
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++)
784     {
785       AliEmcalJet* tmpLeadingJet = NULL;
786
787       if (j==0)
788         tmpLeadingJet = fFirstLeadingJet;
789       else if (j==1)
790         tmpLeadingJet = fSecondLeadingJet;
791       else
792         AliFatal("Trying to exclude more than 2 jets in KT background 2 -- not implemented.");
793
794       if (tmpLeadingJet)
795       {
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)
798         {
799           jetValid = kFALSE;
800           break;
801         }
802       }
803     }
804
805     if(!jetValid)
806       continue;
807    
808     tmpRhos[jetCountAccepted] = jet->Pt() / jet->Area();
809     tmpAreas[jetCountAccepted] = jet->Area();
810     jetCountAccepted++;
811   }
812
813   if (jetCountAccepted > 0)
814   {
815     rhoMedian = TMath::Median(jetCountAccepted, tmpRhos);
816     areaMean   = TMath::Mean(jetCountAccepted, tmpAreas);
817   }
818   #ifdef DEBUGMODE
819     AliInfo("Got KT background 2 density.");
820   #endif
821 }
822
823
824 //________________________________________________________________________
825 Int_t AliAnalysisTaskChargedJetsPA::GetRCBackgroundDensity(Int_t numberExcludeLeadingJets, Double_t& rhoMean, Double_t& rhoMedian, Double_t etaMin, Double_t etaMax, Int_t numberRandCones)
826 {
827   #ifdef DEBUGMODE
828     AliInfo("Getting RC background density.");
829   #endif
830
831   if(numberRandCones == 0)
832     numberRandCones = fNumberRandCones;
833
834   std::vector<AliEmcalJet> tmpCones(numberRandCones);
835
836   // Setting invalid values
837   rhoMean = -1.0;
838   rhoMedian = -1.0;
839
840   // Exclude UP TO numberExcludeLeadingJets
841   if (fNumberSignalJets < 2)
842     numberExcludeLeadingJets = fNumberSignalJets;
843
844   // Search given amount of RCs
845   Int_t numAcceptedRCs = 0;
846   for(Int_t i=0;i<numberRandCones;i++)
847   {
848     Double_t tmpRandConeEta = 0.0;
849     Double_t tmpRandConePhi = 0.0;
850     Double_t excludedJetEta = 0.0;
851     Double_t excludedJetPhi = 0.0;
852
853     // Search random cone in acceptance with no overlap with already excluded jets (leading jets and random cones)
854     Bool_t coneValid = kTRUE;
855
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
859     else
860       tmpRandConeEta = etaMin + fRandom->Rndm()*(etaMax-etaMin);
861
862     tmpRandConePhi = fRandom->Rndm()*TMath::TwoPi();
863
864     // Go through all excluded leading jets and check if there's an overlap
865      
866     for(Int_t j=0;j<numberExcludeLeadingJets;j++)
867     {
868       AliEmcalJet* tmpJet = NULL;
869
870       if (j==0)
871         tmpJet = fFirstLeadingJet;
872       else if (j==1)
873         tmpJet = fSecondLeadingJet;
874       else
875         AliFatal("Trying to exclude more than 2 jets in RC background -- not implemented.");
876
877       excludedJetPhi = tmpJet->Phi();
878       excludedJetEta = tmpJet->Eta();
879       Double_t tmpDeltaPhi = GetDeltaPhi(tmpRandConePhi, excludedJetPhi);
880       
881       if ( tmpDeltaPhi*tmpDeltaPhi + TMath::Abs(tmpRandConeEta-excludedJetEta)*TMath::Abs(tmpRandConeEta-excludedJetEta) <= fRandConeRadius*fRandConeRadius)
882       {
883         coneValid = kFALSE;
884         break;
885       }
886     }
887
888     // RC is accepted, so save it
889     if(coneValid)
890     {
891       AliEmcalJet tmpJet(GetConePt(tmpRandConeEta, tmpRandConePhi, fRandConeRadius), tmpRandConeEta, tmpRandConePhi, 0.0);
892       tmpCones[numAcceptedRCs] = tmpJet;
893       numAcceptedRCs++;
894     }
895   }
896
897   // Calculate Rho and the mean from the RCs (no excluded jets are considered!)
898   if(numAcceptedRCs > 0)
899   {
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());
903
904     rhoMean = TMath::Mean(tmpRho.begin(), tmpRho.end());
905     rhoMedian = 0.0; // NOT IMPLEMENTED because TMath::Median is not working with iterators
906   }
907     
908   #ifdef DEBUGMODE
909     AliInfo("Got RC background density.");
910   #endif
911   return numAcceptedRCs;
912 }
913
914 //________________________________________________________________________
915 void AliAnalysisTaskChargedJetsPA::GetTrackBackgroundDensity(Int_t numberExcludeLeadingJets, Double_t& rhoMean, Double_t& area, Double_t etaMin, Double_t etaMax)
916 {
917   #ifdef DEBUGMODE
918     AliInfo("Getting track background density.");
919   #endif
920
921   Double_t summedTracksPt = 0.0;
922
923   if ((etaMin == 0) && (etaMax == 0))
924   {
925     etaMin = -fTrackEtaWindow;
926     etaMax = +fTrackEtaWindow;
927   }
928
929   // Setting invalid values
930   rhoMean = -1.0;
931   area = -1.0;
932   // Exclude UP TO numberExcludeLeadingJets
933   if (fNumberSignalJets < 2)
934     numberExcludeLeadingJets = fNumberSignalJets;
935
936
937   Int_t   trackCount = fTrackArray->GetEntries();
938   Int_t   trackCountAccepted = 0;
939   for (Int_t i = 0; i < trackCount; i++)
940   {
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))
945       {
946         for (Int_t j = 0; j < numberExcludeLeadingJets; j++)
947         {
948           AliEmcalJet* tmpJet = NULL;
949           if (j==0)
950             tmpJet = fFirstLeadingJet;
951           else if (j==1)
952             tmpJet = fSecondLeadingJet;
953           else
954             AliFatal("Trying to exclude more than 2 jets in track background -- not implemented.");
955
956           if (IsTrackInCone(tmpTrack, tmpJet->Eta(), tmpJet->Phi(), fTrackBackgroundConeRadius))
957           {
958             trackValid = kFALSE;
959             break;
960           }
961         }
962         if (trackValid)
963         {
964           // Add track pt to array
965           summedTracksPt = summedTracksPt + tmpTrack->Pt();
966           trackCountAccepted++;
967         }
968       }
969   }
970
971   if (trackCountAccepted > 0)
972   {
973     Double_t tmpArea = 0.0;
974
975     tmpArea = (2.0*fTrackEtaWindow) * TMath::TwoPi() * (etaMax-etaMin)/(2.0*fTrackEtaWindow); // area of the used eta strip
976     
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()));
982    
983     rhoMean = summedTracksPt/tmpArea;
984     area  = tmpArea;
985   }
986
987   #ifdef DEBUGMODE
988     AliInfo("Got track background density.");
989   #endif
990 }
991
992 //________________________________________________________________________
993 void AliAnalysisTaskChargedJetsPA::GetTrackBackgroundDensity(Int_t numberExcludeLeadingJets, Double_t& rhoMean, Double_t& area, AliEmcalJet* excludeJet1, AliEmcalJet* excludeJet2, Bool_t doSearchPerpendicular)
994 {
995   #ifdef DEBUGMODE
996     AliInfo("Getting track background density.");
997   #endif
998
999   // Setting invalid values
1000   Double_t summedTracksPt = 0.0;
1001   rhoMean = -1.0;
1002   area = -1.0;
1003
1004   Double_t tmpRadius = 0.0;
1005   if (doSearchPerpendicular)
1006     tmpRadius = 0.5*TMath::Pi(); // exclude 90 degrees around jets
1007   else
1008     tmpRadius = fSignalJetRadius;
1009     
1010   numberExcludeLeadingJets = 2; // dijet is excluded here in any case
1011
1012
1013
1014   if (!fTrackArray || !fJetArray)
1015   {
1016     AliError("Could not get the track/jet array to calculate track rho!");
1017     return;
1018   }
1019
1020   Int_t   trackCount = fTrackArray->GetEntries();
1021   Int_t   trackCountAccepted = 0;
1022   for (Int_t i = 0; i < trackCount; i++)
1023   {
1024     AliVTrack* tmpTrack = static_cast<AliVTrack*>(fTrackArray->At(i));
1025     if (IsTrackInAcceptance(tmpTrack))
1026     {
1027       if (IsTrackInCone(tmpTrack, excludeJet1->Eta(), excludeJet1->Phi(), tmpRadius))
1028         continue;
1029
1030       if (numberExcludeLeadingJets > 1)
1031         if (IsTrackInCone(tmpTrack, excludeJet2->Eta(), excludeJet2->Phi(), tmpRadius))
1032           continue;
1033
1034         // Add track pt to array
1035         summedTracksPt = summedTracksPt + tmpTrack->Pt();
1036         trackCountAccepted++;
1037     }
1038   }
1039
1040   if (trackCountAccepted > 0)
1041   {
1042     Double_t tmpArea = 2.0*fTrackEtaWindow*TMath::TwoPi()  - 2*(tmpRadius*tmpRadius * TMath::Pi()); //TPC area - excluding jet area
1043     rhoMean = summedTracksPt/tmpArea;
1044     area = tmpArea;
1045   }
1046
1047   #ifdef DEBUGMODE
1048     AliInfo("Got track background density.");
1049   #endif
1050 }
1051
1052 //________________________________________________________________________
1053 void AliAnalysisTaskChargedJetsPA::Calculate(AliVEvent* event)
1054 {
1055   #ifdef DEBUGMODE
1056     AliInfo("Starting Calculate().");
1057   #endif
1058   ////////////////////// NOTE: initialization & casting
1059
1060   if (!event) {
1061     AliError("??? Event pointer == 0 ???");
1062     return;
1063   }
1064  
1065   if (!fInitialized)
1066     ExecOnce(); // Get tracks, jets, background from arrays if not already given + Init Histos
1067   
1068   // Additional cuts
1069   FillHistogram("hNumberEvents", 0.5); // number of events before manual cuts
1070   
1071   if ((event->GetPrimaryVertex()->GetNContributors() < 1) || (TMath::Abs(event->GetPrimaryVertex()->GetZ()) > fVertexWindow))
1072     return;
1073
1074   if (TMath::Sqrt(event->GetPrimaryVertex()->GetX()*event->GetPrimaryVertex()->GetX() + event->GetPrimaryVertex()->GetY()*event->GetPrimaryVertex()->GetY())  > fVertexMaxR)
1075     return;
1076
1077   FillHistogram("hNumberEvents", 1.5); // number of events after manual cuts
1078
1079   #ifdef DEBUGMODE
1080     AliInfo("Calculate()::Init done.");
1081   #endif
1082
1083   ////////////////////// NOTE: Get Centrality, (Leading)Signal jets and Background
1084
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");
1091
1092   // Get jets
1093   if (fAnalyzeBackground || fAnalyzeJets)
1094     GetSignalJets();
1095
1096   // Get background
1097
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)
1110   {
1111
1112     // Get backgrounds in bins of eta
1113     for(Int_t i = 0; i<fBackgroundEtaBins; i++)
1114     {
1115       // scheme: etaMin = RangeMin + l*binN; etaMax = RangeMin + l*(binN+1)
1116
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);
1123
1124     }
1125     Int_t tmpNRCs = 0;
1126
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]);
1133
1134   }
1135
1136   #ifdef DEBUGMODE
1137     AliInfo("Calculate()::Centrality&SignalJets&Background-Calculation done.");
1138   #endif
1139   ////////////////////// NOTE: Pythia histograms
1140   if(fAnalyzePythia)
1141   {
1142     FillHistogram("hPythiaPtHard", GetPtHard());
1143     FillHistogram("hPythiaNTrials", GetPtHardBin()-0.1, fTrials);
1144     FillHistogram("hPythiaXSection", GetPtHardBin()-0.1, fCrossSection);
1145
1146     #ifdef DEBUGMODE
1147       AliInfo("Calculate()::Pythia done.");
1148     #endif
1149   }
1150
1151   ////////////////////// NOTE: Track & QA histograms
1152
1153   if (fAnalyzeQA)
1154   {
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);
1158
1159     Int_t trackCountAcc = 0;
1160     Int_t nTracks = fTrackArray->GetEntries();
1161     for (Int_t i = 0; i < nTracks; i++)
1162     {
1163       AliVTrack* track = static_cast<AliVTrack*>(fTrackArray->At(i));
1164       if (IsTrackInAcceptance(track))
1165       {
1166         FillHistogram("hTrackPhiEta", track->Phi(),track->Eta(), 1);
1167         FillHistogram("hTrackPt", track->Pt());
1168         FillHistogram("hTrackEta", track->Eta());
1169         FillHistogram("hTrackCharge", track->Charge());
1170         trackCountAcc++;
1171       }
1172     }
1173     FillHistogram("hTrackCountAcc", trackCountAcc, centralityPercentile);
1174
1175     if (fHasClusters)
1176     {
1177       Int_t clusterCountAcc = 0;
1178       Int_t nClusters = fClusterArray->GetEntries();
1179       for (Int_t i = 0; i < nClusters; i++)
1180       {
1181         AliVCluster* cluster = static_cast<AliVCluster*>(fClusterArray->At(i));
1182         if (IsClusterInAcceptance(cluster))
1183         {
1184           FillHistogram("hClusterE", cluster->E());
1185           clusterCountAcc++;
1186         }
1187       }
1188       FillHistogram("hClusterCountAcc", clusterCountAcc, centralityPercentile);
1189     }
1190   }
1191   #ifdef DEBUGMODE
1192     AliInfo("Calculate()::QA done.");
1193   #endif
1194
1195   ////////////////////// NOTE: Jet analysis and calculations
1196
1197   if (fAnalyzeJets)
1198   {
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());
1205
1206     // ### Dijets I ###
1207     if(fNumberSignalJets >= 2)
1208     {
1209       FillHistogram("hLeadingJetDeltaPhi", GetDeltaPhi(fFirstLeadingJet->Phi(), fSecondLeadingJet->Phi()));
1210       FillHistogram("hLeadingJetDeltaPhiPt", GetDeltaPhi(fFirstLeadingJet->Phi(), fSecondLeadingJet->Phi()), fFirstLeadingJet->Pt());
1211
1212       if (IsDijet(fFirstLeadingJet, fSecondLeadingJet)) // Gettin' the money
1213       {
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);
1220       }
1221     }
1222
1223     // SIGNAL JET ANALYSIS
1224     for (Int_t i = 0; i<fNumberSignalJets; i++)
1225     {
1226       AliEmcalJet* tmpJet = fSignalJets[i];
1227
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());
1235       
1236       // Background subtracted spectra
1237       
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]));
1244
1245
1246       Double_t tmpCorrFactor = GetJetBackgroundCorrFactor(tmpJet->Eta(), ktBackgroundRhoMedian[fBackgroundEtaBins]);
1247       FillHistogram("hAppliedEtaCorrectionFactor", tmpCorrFactor);
1248       tmpCorrFactor = GetJetBackgroundCorrFactor(tmpJet->Eta(), ktBackground2RhoMedian[fBackgroundEtaBins]);
1249       FillHistogram("hAppliedEtaCorrectionFactor2", tmpCorrFactor);
1250
1251       // Signal jet vs. signal jet
1252       for (Int_t j = i+1; j<fNumberSignalJets; j++)
1253       {
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()));
1257
1258         // ### Dijets II ###
1259         if (IsDijet(tmpJet, tmpJet2)) // Gettin' the money
1260         {
1261           FillHistogram("hDijet2ConstituentsPt", tmpJet->Pt()); FillHistogram("hDijet2ConstituentsPt", tmpJet2->Pt());
1262           FillHistogram("hDijet2LeadingJetPt", fFirstLeadingJet->Pt());
1263           FillHistogram("hDijet2PtCorrelation", tmpJet->Pt(), tmpJet2->Pt());
1264         }
1265       }
1266     }
1267   } //endif AnalyzeJets
1268
1269   #ifdef DEBUGMODE
1270     AliInfo("Calculate()::Jets done.");
1271   #endif
1272   ////////////////////// NOTE: Background analysis
1273
1274   if (fAnalyzeBackground)
1275   {
1276
1277     Int_t leadingJetIds[] = {-1, -1};
1278     GetLeadingJets(fBackgroundJetArray, &leadingJetIds[0], kFALSE);
1279
1280     for (Int_t i = 0; i < fBackgroundJetArray->GetEntries(); i++)
1281     {
1282       AliEmcalJet* jet = static_cast<AliEmcalJet*>(fBackgroundJetArray->At(i));
1283       if (!jet)
1284       {
1285         AliError(Form("%s: Could not receive kt jet %d", GetName(), i));
1286         continue;
1287       }
1288       if (!IsBackgroundJetInAcceptance(jet))
1289         continue;
1290       if (!((jet->Eta() >= -fBackgroundJetEtaWindow) && (jet->Eta() < fBackgroundJetEtaWindow)))
1291         continue;
1292       
1293       FillHistogram("hKTJetPhiEta", jet->Phi(),jet->Eta());
1294       if(i==leadingJetIds[0])
1295         FillHistogram("hKTLeadingJetPhiEta", jet->Phi(),jet->Eta());
1296             
1297     }
1298
1299     // ############# RC, Track, and KT background calculations
1300     Double_t etaMin = 0;
1301     for (Int_t i=0;i<fBackgroundEtaBins;i++)
1302     {
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.)
1309       {
1310         FillHistogram("hRCBackgroundMostCentral", etaMin, rcBackgroundRhoMean[i]);
1311         FillHistogram("hTrackBackgroundMostCentral", etaMin, trackBackgroundRhoMean[i]);
1312         FillHistogram("hKTBackgroundMostCentral", etaMin, ktBackgroundRhoMedian[i]);
1313         FillHistogram("hKTBackground2MostCentral", etaMin, ktBackground2RhoMedian[i]);
1314       }
1315       else if(centralityPercentile >= 80.)
1316       {
1317         FillHistogram("hRCBackgroundMostPeripheral", etaMin, rcBackgroundRhoMean[i]);
1318         FillHistogram("hTrackBackgroundMostPeripheral", etaMin, trackBackgroundRhoMean[i]);
1319         FillHistogram("hKTBackgroundMostPeripheral", etaMin, ktBackgroundRhoMedian[i]);
1320         FillHistogram("hKTBackground2MostPeripheral", etaMin, ktBackground2RhoMedian[i]);
1321       }
1322     }
1323
1324     FillHistogram("hRCBackgroundVsCentrality", rcBackgroundRhoMean[fBackgroundEtaBins], centralityPercentile);
1325     FillHistogram("hTrackBackgroundVsCentrality", trackBackgroundRhoMean[fBackgroundEtaBins], centralityPercentile);
1326     FillHistogram("hKTBackgroundVsCentrality", ktBackgroundRhoMedian[fBackgroundEtaBins], centralityPercentile);
1327     FillHistogram("hKTBackground2VsCentrality", ktBackground2RhoMedian[fBackgroundEtaBins], centralityPercentile);
1328
1329     if (dijetBackground >= 0)
1330     {
1331       // Background in Dijet events
1332       FillHistogram("hDijetBackground", dijetBackground); 
1333       if(centralityPercentile <= 20.)
1334         FillHistogram("hDijetBackgroundMostCentral", dijetBackground); 
1335       FillHistogram("hDijetBackgroundVsCentrality", dijetBackground, centralityPercentile);
1336     }
1337     if (dijetBackgroundPerpendicular >= 0)
1338     {
1339       // Background in Dijet events
1340       FillHistogram("hDijetBackgroundPerpendicular", dijetBackgroundPerpendicular); 
1341       if(centralityPercentile <= 20.)
1342         FillHistogram("hDijetBackgroundPerpendicularMostCentral", dijetBackgroundPerpendicular); 
1343       FillHistogram("hDijetBackgroundPerpendicularVsCentrality", dijetBackgroundPerpendicular, centralityPercentile);
1344     }
1345
1346     // ########## Delta pT calculations (most central, kt is eta corrected)
1347     if (centralityPercentile <= 20.)
1348     {
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;
1353
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);
1358
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);
1363
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);
1368
1369       GetDeltaPt(tmpDeltaPtKT2Eta2Excl, ktBackground2RhoMedian[fKTDeltaPtEtaBin], 2, fKTDeltaPtEtaBin);
1370
1371       // kT Background
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);
1378
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);
1387
1388       // RC Background
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);
1395       // TR Background
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);
1402     }
1403   }
1404   
1405   #ifdef DEBUGMODE
1406     AliInfo("Calculate()::Background done.");
1407   #endif
1408
1409
1410 //________________________________________________________________________
1411 Bool_t AliAnalysisTaskChargedJetsPA::Notify()
1412 {
1413   // Implemented Notify() to read the cross sections
1414   // and number of trials from pyxsec.root
1415   // 
1416   #ifdef DEBUGMODE
1417     AliInfo("Notify started.");
1418   #endif
1419
1420   if(fAnalyzePythia)
1421   {
1422     TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
1423     TFile *currFile = tree->GetCurrentFile();
1424
1425     TString file(currFile->GetName());
1426
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,"");
1431     }
1432     else {
1433       // not an archive take the basename....
1434       file.ReplaceAll(gSystem->BaseName(file.Data()),"");
1435     }
1436    
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
1438     if(!fxsec){
1439       // next trial fetch the histgram file
1440       fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec_hists.root"));
1441       if(!fxsec){
1442           // not a severe condition but inciate that we have no information
1443         return kFALSE;
1444       }
1445       else{
1446         // find the tlist we want to be independtent of the name so use the Tkey
1447         TKey* key = (TKey*)fxsec->GetListOfKeys()->At(0); 
1448         if(!key){
1449           fxsec->Close();
1450           return kFALSE;
1451         }
1452         TList *list = dynamic_cast<TList*>(key->ReadObj());
1453         if(!list){
1454           fxsec->Close();
1455           return kFALSE;
1456         }
1457         fCrossSection = ((TProfile*)list->FindObject("h1Xsec"))->GetBinContent(1);
1458         fTrials  = ((TH1F*)list->FindObject("h1Trials"))->GetBinContent(1);
1459         fxsec->Close();
1460       }
1461     } // no tree pyxsec.root
1462     else {
1463       TTree *xtree = (TTree*)fxsec->Get("Xsection");
1464       if(!xtree){
1465         fxsec->Close();
1466         return kFALSE;
1467       }
1468       UInt_t   ntrials  = 0;
1469       Double_t  xsection  = 0;
1470       xtree->SetBranchAddress("xsection",&xsection);
1471       xtree->SetBranchAddress("ntrials",&ntrials);
1472       xtree->GetEntry(0);
1473       fTrials = ntrials;
1474       fCrossSection = xsection;
1475       fxsec->Close();
1476     }
1477     #ifdef DEBUGMODE
1478       AliInfo("Notify ended.");
1479     #endif
1480   }
1481   return kTRUE;
1482 }
1483
1484 //________________________________________________________________________
1485 inline void AliAnalysisTaskChargedJetsPA::FillHistogram(const char * key, Double_t x)
1486 {
1487   TH1* tmpHist = static_cast<TH1*>(fOutputList->FindObject(GetHistoName(key)));
1488   if(!tmpHist)
1489   {
1490     AliInfo(Form("Cannot find histogram <%s> ",key)) ;
1491     return;
1492   }
1493
1494   tmpHist->Fill(x);
1495 }
1496
1497 //________________________________________________________________________
1498 inline void AliAnalysisTaskChargedJetsPA::FillHistogram(const char * key, Double_t x, Double_t y)
1499 {
1500   TH1* tmpHist = static_cast<TH1*>(fOutputList->FindObject(GetHistoName(key)));
1501   if(!tmpHist)
1502   {
1503     AliInfo(Form("Cannot find histogram <%s> ",key));
1504     return;
1505   }
1506
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
1511 }
1512
1513 //________________________________________________________________________
1514 inline void AliAnalysisTaskChargedJetsPA::FillHistogram(const char * key, Double_t x, Double_t y, Double_t add)
1515 {
1516   TH2* tmpHist = static_cast<TH2*>(fOutputList->FindObject(GetHistoName(key)));
1517   if(!tmpHist)
1518   {
1519     AliInfo(Form("Cannot find histogram <%s> ",key));
1520     return;
1521   }
1522   
1523   tmpHist->Fill(x,y,add);
1524 }
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)
1527 {
1528   T* tmpHist = new T(GetHistoName(name), GetHistoName(title), xBins, xMin, xMax);
1529
1530   tmpHist->GetXaxis()->SetTitle(xTitle);
1531   tmpHist->GetYaxis()->SetTitle(yTitle);
1532   tmpHist->SetOption(options);
1533   tmpHist->SetMarkerStyle(kFullCircle);
1534   tmpHist->Sumw2();
1535
1536   fHistList->Add(tmpHist);
1537   fHistCount++;
1538   
1539   return tmpHist;
1540 }
1541
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)
1544 {
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);
1551   tmpHist->Sumw2();
1552
1553   fHistList->Add(tmpHist);
1554   fHistCount++;
1555
1556   return tmpHist;
1557 }
1558
1559 //________________________________________________________________________
1560 void AliAnalysisTaskChargedJetsPA::Terminate(Option_t *)
1561 {
1562   PostData(1, fOutputList);
1563
1564   // Mandatory
1565   fOutputList = dynamic_cast<TList*> (GetOutputData(1)); // '1' refers to the output slot
1566   if (!fOutputList) {
1567     printf("ERROR: Output list not available\n");
1568     return;
1569   }
1570 }
1571
1572 //________________________________________________________________________
1573 AliAnalysisTaskChargedJetsPA::~AliAnalysisTaskChargedJetsPA()
1574 {
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()) {
1578     delete fOutputList;
1579   }
1580 }
1581
1582 //________________________________________________________________________
1583 void AliAnalysisTaskChargedJetsPA::UserCreateOutputObjects()
1584 {
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.
1588
1589   fRandom = new TRandom3(0);
1590   
1591   fOutputList = new TList();
1592   fOutputList->SetOwner(); // otherwise it produces leaks in merging
1593
1594   PostData(1, fOutputList);
1595 }
1596
1597 //________________________________________________________________________
1598 void AliAnalysisTaskChargedJetsPA::UserExec(Option_t *) 
1599 {
1600   #ifdef DEBUGMODE
1601     AliInfo("UserExec() started.");
1602   #endif
1603
1604   Calculate(InputEvent());
1605         
1606   PostData(1, fOutputList);
1607 }