Charged jets (pPb): Added histograms
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / UserTasks / AliAnalysisTaskChargedJetsPA.cxx
1 #ifndef ALIANALYSISTASKSE_H
2 #include <Riostream.h>
3 #include <TROOT.h>
4 #include <TFile.h>
5 #include <TChain.h>
6 #include <TTree.h>
7 #include <TKey.h>
8 #include <TProfile.h>
9 #include <TProfile2D.h>
10 #include <TH1F.h>
11 #include <TH2F.h>
12 #include <TCanvas.h>
13 #include <TList.h>
14 #include <TClonesArray.h>
15 #include <TObject.h>
16 #include <TMath.h>
17 #include <TSystem.h>
18 #include <TInterpreter.h>
19 #include <TH1.h>
20 #include "AliAnalysisTask.h"
21 #include "AliCentrality.h"
22 #include "AliStack.h"
23 #include "AliESDEvent.h"
24 #include "AliESDInputHandler.h"
25 #include "AliAODEvent.h"
26 #include "AliAODHandler.h"
27 #include "AliAnalysisManager.h"
28 #include "AliAnalysisTaskSE.h"
29 #endif
30
31 #include <time.h>
32 #include <TRandom3.h>
33 #include "AliGenPythiaEventHeader.h"
34 #include "AliAODMCHeader.h"
35 #include "AliMCEvent.h"
36 #include "AliLog.h"
37 #include <AliEmcalJet.h>
38 #include <AliPicoTrack.h>
39 #include "AliVEventHandler.h"
40 #include "AliVParticle.h"
41 #include "AliAODMCParticle.h"
42 #include "AliAnalysisUtils.h"
43 #include "AliRhoParameter.h"
44 #include "TVector3.h"
45
46 #include "AliAnalysisTaskChargedJetsPA.h"
47 using std::min;
48
49 //TODO: Not accessing the particles when using MC
50 //TODO: FillHistogram can be done better with virtual TH1(?)
51 ClassImp(AliAnalysisTaskChargedJetsPA)
52
53 // ######################################################################################## DEFINE HISTOGRAMS
54 void AliAnalysisTaskChargedJetsPA::Init()
55 {
56   #ifdef DEBUGMODE
57     AliInfo("Creating histograms.");
58   #endif
59
60   TH1D* tmpHisto = AddHistogram1D<TH1D>("hNumberEvents", "Number of events (0 = before cuts, 1 = after cuts)", "", 2, 0, 2, "stage","N^{Events}/cut");
61   tmpHisto->GetXaxis()->SetBinLabel(1, "Before cuts");
62   tmpHisto->GetXaxis()->SetBinLabel(2, "After cuts");
63
64   tmpHisto = AddHistogram1D<TH1D>("hEventAcceptance", "Accepted events (0 = before cuts, 1 = after pile up, 2 = after vertex)", "", 3, 0, 3, "stage","N^{Events}/cut");
65   tmpHisto->GetXaxis()->SetBinLabel(1, "Before cuts");
66   tmpHisto->GetXaxis()->SetBinLabel(2, "After pile up");
67   tmpHisto->GetXaxis()->SetBinLabel(3, "After vertex");
68
69   tmpHisto = AddHistogram1D<TH1D>("hTrackAcceptance", "Accepted tracks (0 = before cuts, 1 = after eta, 2 = after pT)", "", 3, 0, 3, "stage","N^{Tracks}/cut");
70   tmpHisto->GetXaxis()->SetBinLabel(1, "Before cuts");
71   tmpHisto->GetXaxis()->SetBinLabel(2, "After eta");
72   tmpHisto->GetXaxis()->SetBinLabel(3, "After p_{T}");
73
74   tmpHisto = AddHistogram1D<TH1D>("hJetAcceptance", "Accepted jets (0 = before cuts, 1 = after eta, 2 = after pT, 3 = after area)", "", 4, 0, 4, "stage","N^{Jets}/cut");
75   tmpHisto->GetXaxis()->SetBinLabel(1, "Before cuts");
76   tmpHisto->GetXaxis()->SetBinLabel(2, "After eta");
77   tmpHisto->GetXaxis()->SetBinLabel(3, "After p_{T}");
78   tmpHisto->GetXaxis()->SetBinLabel(4, "After area");
79
80   // NOTE: Jet histograms
81   if (fAnalyzeJets)
82   {
83     // ######## Jet spectra
84     TH2* tmpHisto2D = AddHistogram2D<TH2D>("hJetPtCutStages", "Jets p_{T} distribution", "", 500, -50., 200., 4, 0, 4, "p_{T} (GeV/c)","Cut stage","dN^{Jets}/dp_{T}");
85     tmpHisto2D->GetYaxis()->SetBinLabel(1, "Before cuts");
86     tmpHisto2D->GetYaxis()->SetBinLabel(2, "After eta");
87     tmpHisto2D->GetYaxis()->SetBinLabel(3, "After p_{T}");
88     tmpHisto2D->GetYaxis()->SetBinLabel(4, "After area");
89     AddHistogram2D<TH2D>("hJetPtBgrdSubtractedKTImprovedCMS", "Jets p_{T} distribution, KT background (Improved CMS) subtracted", "", 500, -50., 200., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Jets}/dp_{T}");
90     AddHistogram2D<TH2D>("hJetPtBgrdSubtractedExternal_Phi1", "Jets p_{T} distribution, external background (Improved CMS) subtracted (1st part of azimuth)", "", 500, -50., 200., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Jets}/dp_{T}");    
91     AddHistogram2D<TH2D>("hJetPtBgrdSubtractedExternal_Phi2", "Jets p_{T} distribution, external background (Improved CMS) subtracted (2nd part of azimuth)", "", 500, -50., 200., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Jets}/dp_{T}");    
92
93     AddHistogram2D<TH2D>("hJetPtBgrdSubtractedPP", "Jets p_{T} distribution, pp background subtracted", "", 500, -50., 200., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Jets}/dp_{T}");
94     AddHistogram2D<TH2D>("hJetPtBgrdSubtractedExternal", "Jets p_{T} distribution, external bgrd. subtracted", "", 500, -50., 200., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Jets}/dp_{T}");    
95
96
97
98     AddHistogram2D<TProfile2D>("hJetPtSubtractedRhoKTImprovedCMS", "Mean subtracted KT (CMS w/o signal) background from jets", "COLZ", 600, 0, 150, fNumberOfCentralityBins, 0, 100, "Jet p_{T}", "Centrality", "#rho mean");
99     AddHistogram2D<TProfile2D>("hJetPtSubtractedRhoExternal", "Mean subtracted KT (External) background from jets", "COLZ", 600, 0, 150, fNumberOfCentralityBins, 0, 100, "Jet p_{T}", "Centrality", "#rho mean");
100     AddHistogram2D<TProfile2D>("hJetPtSubtractedRhoPP", "Mean subtracted KT (pp from Michal) background from jets", "COLZ", 600, 0, 150, fNumberOfCentralityBins, 0, 100, "Jet p_{T}", "Centrality", "#rho mean");
101
102     if(fAnalyzeDeprecatedBackgrounds)
103     {
104       AddHistogram2D<TH2D>("hJetPtBgrdSubtractedTR", "Jets p_{T} distribution, TR background (Cone R=0.6 around jets excluded) subtracted", "", 500, -50., 200., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Jets}/dp_{T}");
105       AddHistogram2D<TH2D>("hJetPtBgrdSubtractedKTPbPb", "Jets p_{T} distribution, KT background (PbPb w/o ghosts) subtracted", "", 500, -50., 200., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Jets}/dp_{T}");
106       AddHistogram2D<TH2D>("hJetPtBgrdSubtractedKTPbPbWithGhosts", "Jets p_{T} distribution, KT background (PbPb w/ ghosts) subtracted", "", 500, -50., 200., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Jets}/dp_{T}");
107       AddHistogram2D<TH2D>("hJetPtBgrdSubtractedKTCMS", "Jets p_{T} distribution, KT background (CMS) subtracted", "", 500, -50., 200., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Jets}/dp_{T}");    
108       AddHistogram2D<TH2D>("hJetPtBgrdSubtractedKTMean", "Jets p_{T} distribution, KT background (Mean) subtracted", "", 500, -50., 200., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Jets}/dp_{T}");    
109       AddHistogram2D<TH2D>("hJetPtBgrdSubtractedKTTrackLike", "Jets p_{T} distribution, KT background (track-like) subtracted", "", 500, -50., 200., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Jets}/dp_{T}");
110     }
111
112     // ######## Jet profiles
113     if(fAnalyzeJetProfile)
114     {
115       AddHistogram2D<TH2D>("hJetProfile10GeV", "Jet profile, cone p_{T}/jet p_{T} vs. jet radius, jet p_{T} > 10 GeV", "", 12, 0, 0.6,200, 0., 2., "Cone radius","dN^{Jets}/dR", "Ratio");
116       AddHistogram2D<TH2D>("hJetProfile20GeV", "Jet profile, cone p_{T}/jet p_{T} vs. jet radius, jet p_{T} > 20 GeV", "", 12, 0, 0.6,200, 0., 2., "Cone radius","dN^{Jets}/dR", "Ratio");
117       AddHistogram2D<TH2D>("hJetProfile30GeV", "Jet profile, cone p_{T}/jet p_{T} vs. jet radius, jet p_{T} > 30 GeV", "", 12, 0, 0.6,200, 0., 2., "Cone radius","dN^{Jets}/dR", "Ratio");
118       AddHistogram2D<TH2D>("hJetProfile40GeV", "Jet profile, cone p_{T}/jet p_{T} vs. jet radius, jet p_{T} > 40 GeV", "", 12, 0, 0.6,200, 0., 2., "Cone radius","dN^{Jets}/dR", "Ratio");
119       AddHistogram2D<TH2D>("hJetProfile50GeV", "Jet profile, cone p_{T}/jet p_{T} vs. jet radius, jet p_{T} > 50 GeV", "", 12, 0, 0.6,200, 0., 2., "Cone radius","dN^{Jets}/dR", "Ratio");
120       AddHistogram2D<TH2D>("hJetProfile60GeV", "Jet profile, cone p_{T}/jet p_{T} vs. jet radius, jet p_{T} > 60 GeV", "", 12, 0, 0.6,200, 0., 2., "Cone radius","dN^{Jets}/dR", "Ratio");
121       AddHistogram2D<TH2D>("hJetProfile70GeV", "Jet profile, cone p_{T}/jet p_{T} vs. jet radius, jet p_{T} > 70 GeV", "", 12, 0, 0.6,200, 0., 2., "Cone radius","dN^{Jets}/dR", "Ratio");
122     }
123     // ######## Jet stuff
124     AddHistogram2D<TH2D>("hJetConstituentPt", "Jet constituents p_{T} distribution", "", 500, -50., 200., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Tracks}/dp_{T}");
125     AddHistogram1D<TH1D>("hJetCountAll", "Number of Jets", "", 200, 0., 200., "N jets","dN^{Events}/dN^{Jets}");
126     AddHistogram1D<TH1D>("hJetCountAccepted", "Number of accepted Jets", "", 200, 0., 200., "N jets","dN^{Events}/dN^{Jets}");
127     AddHistogram2D<TH2D>("hJetCount", "Correlation jets/accepted jets", "", 200, 0., 200., 200, 0., 200., "N jets","N jets accepted", "d^{2}N^{Events}/dN^{Jets dN^{Jets, acc}}");
128     AddHistogram1D<TH1D>("hLeadingJetPt", "Leading jet p_{T}", "", 500, -50., 200., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
129     AddHistogram1D<TH1D>("hSecondLeadingJetPt", "Second leading jet p_{T}", "", 500, -50., 200., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
130     AddHistogram1D<TH1D>("hCorrectedLeadingJetPt", "Corrected leading jet p_{T}", "", 500, -50., 200., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
131     AddHistogram1D<TH1D>("hCorrectedSecondLeadingJetPt", "Corrected second leading jet p_{T}", "", 500, -50., 200., "p_{T} (GeV/c)","dN^{Jets}/dp_{T}");
132     AddHistogram1D<TH1D>("hJetDeltaPhi", "Jets combinatorial #Delta #phi", "", 250, 0., TMath::Pi(), "#Delta #phi","dN^{Jets}/d(#Delta #phi)");
133     AddHistogram1D<TH1D>("hLeadingJetDeltaPhi", "1st and 2nd leading jet #Delta #phi", "", 250, 0., TMath::Pi(), "#Delta #phi","dN^{Jets}/d(#Delta #phi)");
134   }
135
136   // NOTE: Jet background histograms
137   if (fAnalyzeBackground)
138   {
139     // ########## Default background estimates
140     AddHistogram2D<TH2D>("hKTBackgroundImprovedCMS", "KT background density (Improved CMS approach)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
141     AddHistogram2D<TH2D>("hKTBackgroundImprovedCMSExternal", "KT background density (Improved CMS approach from external task)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
142     AddHistogram2D<TH2D>("hKTBackgroundImprovedCMSExternal20GeV", "KT background density (Improved CMS approach from external task, jet p_{T} > 20 GeV)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
143     AddHistogram2D<TH2D>("hPPBackground", "PP background density (Michals approach)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
144
145     AddHistogram2D<TH2D>("hDeltaPtExternalBgrd", "Background fluctuations #delta p_{T} (KT, External)", "", 1801, -40.0, 80.0, fNumberOfCentralityBins, 0, 100, "#delta p_{T} (GeV/c)","Centrality","dN^{Jets}/d#delta p_{T}");
146     AddHistogram2D<TH2D>("hDeltaPtExternalBgrdVsPt", "Background fluctuations #delta p_{T} (KT, External, in p_{T} bins)", "", 1801, -40.0, 80.0, 200, 0, 200, "#delta p_{T} (GeV/c)","Raw jet p_{T}","dN^{Jets}/d#delta p_{T}");
147     AddHistogram2D<TH2D>("hDeltaPtKTImprovedCMS", "Background fluctuations #delta p_{T} (KT, Improved CMS-like)", "", 1801, -40.0, 80.0, fNumberOfCentralityBins, 0, 100, "#delta p_{T} (GeV/c)","Centrality","dN^{Jets}/d#delta p_{T}");
148     AddHistogram2D<TH2D>("hDeltaPtKTImprovedCMSPartialExclusion", "Background fluctuations #delta p_{T} (KT, Improved CMS-like, partial jet exclusion)", "", 1801, -40.0, 80.0, fNumberOfCentralityBins, 0, 100, "#delta p_{T} (GeV/c)","Centrality","dN^{Jets}/d#delta p_{T}");
149     AddHistogram2D<TH2D>("hDeltaPtKTImprovedCMSPartialExclusion_Signal", "Background fluctuations #delta p_{T} (KT, Improved CMS-like, partial jet exclusion w/ 1/N_{sig} probability)", "", 1801, -40.0, 80.0, fNumberOfCentralityBins, 0, 100, "#delta p_{T} (GeV/c)","Centrality","dN^{Jets}/d#delta p_{T}");
150     AddHistogram2D<TH2D>("hDeltaPtKTImprovedCMSFullExclusion", "Background fluctuations #delta p_{T} (KT, Improved CMS-like, full leading jet exclusion)", "", 1801, -40.0, 80.0, fNumberOfCentralityBins, 0, 100, "#delta p_{T} (GeV/c)","Centrality","dN^{Jets}/d#delta p_{T}");
151     AddHistogram2D<TH2D>("hDeltaPtNoBackground", "Background fluctuations #delta p_{T} (No background)", "", 1801, -40.0, 80.0, fNumberOfCentralityBins, 0, 100, "#delta p_{T} (GeV/c)","Centrality","dN^{Jets}/d#delta p_{T}");
152
153     AddHistogram1D<TProfile>("hKTMeanBackgroundImprovedCMS", "KT background mean (Improved CMS approach)", "", 100, 0, 100, "Centrality", "#rho mean");
154
155     if(fAnalyzeDeprecatedBackgrounds)
156     {
157       // ########## Different background estimates
158       AddHistogram2D<TH2D>("hKTBackgroundPbPb", "KT background density (PbPb approach, no ghosts)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
159       AddHistogram2D<TH2D>("hKTBackgroundPbPbWithGhosts", "KT background density (PbPb approach w/ ghosts)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
160       AddHistogram2D<TH2D>("hKTBackgroundCMS", "KT background density (CMS approach)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
161       AddHistogram2D<TH2D>("hKTBackgroundMean", "KT background density (Mean approach)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
162       AddHistogram2D<TH2D>("hKTBackgroundTrackLike", "KT background density (Track-like approach)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
163
164       AddHistogram2D<TH2D>("hTRBackgroundNoExcl", "TR background density (No signal excluded)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
165       AddHistogram2D<TH2D>("hTRBackgroundCone02", "TR background density (Cones R=0.2 around signal jets excluded)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
166       AddHistogram2D<TH2D>("hTRBackgroundCone04", "TR background density (Cones R=0.4 around signal jets excluded)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
167       AddHistogram2D<TH2D>("hTRBackgroundCone06", "TR background density (Cones R=0.6 around signal jets excluded)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
168       AddHistogram2D<TH2D>("hTRBackgroundCone08", "TR background density (Cones R=0.8 around signal jets excluded)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
169       AddHistogram2D<TH2D>("hTRBackgroundExact",  "TR background density (signal jets exactly excluded)", "LEGO2", 400, 0., 40., fNumberOfCentralityBins, 0, 100, "#rho (GeV/c)","Centrality", "dN^{Events}/d#rho");
170
171       // ########## Delta Pt
172       AddHistogram2D<TH2D>("hDeltaPtKTPbPb", "Background fluctuations #delta p_{T} (KT, PbPb w/o ghosts)", "", 1801, -40.0, 80.0, fNumberOfCentralityBins, 0, 100, "#delta p_{T} (GeV/c)","Centrality","dN^{Jets}/d#delta p_{T}");
173       AddHistogram2D<TH2D>("hDeltaPtKTPbPbWithGhosts", "Background fluctuations #delta p_{T} (KT, PbPb w/ ghosts)", "", 1801, -40.0, 80.0, fNumberOfCentralityBins, 0, 100, "#delta p_{T} (GeV/c)","Centrality","dN^{Jets}/d#delta p_{T}");
174       AddHistogram2D<TH2D>("hDeltaPtKTCMS", "Background fluctuations #delta p_{T} (KT, CMS-like)", "", 1801, -40.0, 80.0, fNumberOfCentralityBins, 0, 100, "#delta p_{T} (GeV/c)","Centrality","dN^{Jets}/d#delta p_{T}");
175       AddHistogram2D<TH2D>("hDeltaPtKTMean", "Background fluctuations #delta p_{T} (KT, Mean)", "", 1801, -40.0, 80.0, fNumberOfCentralityBins, 0, 100, "#delta p_{T} (GeV/c)","Centrality","dN^{Jets}/d#delta p_{T}");
176       AddHistogram2D<TH2D>("hDeltaPtKTTrackLike", "Background fluctuations #delta p_{T} (KT, track-like)", "", 1801, -40.0, 80.0, fNumberOfCentralityBins, 0, 100, "#delta p_{T} (GeV/c)","Centrality","dN^{Jets}/d#delta p_{T}");
177       AddHistogram2D<TH2D>("hDeltaPtTR", "Background fluctuations #delta p_{T} (TR, cone R=0.6)", "", 1801, -40.0, 80.0, fNumberOfCentralityBins, 0, 100,  "#delta p_{T} (GeV/c)","Centrality","dN^{Jets}/d#delta p_{T}");
178
179       // ########## Profiles for background means vs. centrality
180       AddHistogram1D<TProfile>("hKTMeanBackgroundPbPb", "KT background mean (PbPb approach w/o ghosts)", "", fNumberOfCentralityBins, 0, 100, "Centrality", "#rho mean");
181       AddHistogram1D<TProfile>("hKTMeanBackgroundPbPbWithGhosts", "KT background mean (PbPb approach)", "", fNumberOfCentralityBins, 0, 100, "Centrality", "#rho mean");
182       AddHistogram1D<TProfile>("hKTMeanBackgroundCMS", "KT background mean (CMS approach)", "", fNumberOfCentralityBins, 0, 100, "Centrality", "#rho mean");
183       AddHistogram1D<TProfile>("hKTMeanBackgroundMean", "KT background mean (Mean approach)", "",  fNumberOfCentralityBins, 0, 100, "Centrality", "#rho mean");
184       AddHistogram1D<TProfile>("hKTMeanBackgroundTPC", "KT background mean (Track-like approach)", "", fNumberOfCentralityBins, 0, 100, "Centrality", "#rho mean");
185       AddHistogram1D<TProfile>("hTRMeanBackground", "TR background mean", "", fNumberOfCentralityBins, 0, 100, "Centrality", "#rho mean");
186     }
187   }
188   // NOTE: Jet constituent correlations
189   if(fAnalyzeMassCorrelation)
190   {
191     AddHistogram1D<TH1D>("hJetMassFromConstituents", "Jet mass by mass of charged constituents", "", 200, 0., 200., "N jets","dN^{Jets}/dN^{mass}");
192     AddHistogram1D<TH1D>("hJetMass", "Jet mass from fastjet", "", 200, 0., 200., "N jets","dN^{Jets}/dN^{mass}");
193
194     AddHistogram2D<TH2D>("hJetPtVsMass", "Correlation jet pt/summed constituent jet", "", 400, 0., 100., 400, 0., 100., "p_{T}","Mass", "d^{2}N}/dN^{p_{T}dM");
195     AddHistogram2D<TH2D>("hJetPtVsJetMass", "Correlation jet pt/summed jet", "", 400, 0., 100., 400, 0., 100., "p_{T}","Mass", "d^{2}N}/dN^{p_{T}dM");
196
197     AddHistogram2D<TH2D>("hJetPtVsProtonCount", "Correlation jet pt/amount of proton in jet", "", 400, 0., 100., 20, 0., 1., "p_{T}","Proton percentage", "d^{2}N}/dN^{p_{T}dPercentage");
198     AddHistogram2D<TH2D>("hJetPtVsPionCount", "Correlation jet pt/amount of pion in jet", "", 400, 0., 100., 20, 0., 1., "p_{T}","Pion percentage", "d^{2}N}/dN^{p_{T}dPercentage");
199     AddHistogram2D<TH2D>("hJetPtVsKaonCount", "Correlation jet pt/amount of kaon in jet", "", 400, 0., 100., 20, 0., 1., "p_{T}","Kaon percentage", "d^{2}N}/dN^{p_{T}dPercentage");
200     AddHistogram2D<TH2D>("hJetPtVsElectronCount", "Correlation jet pt/amount of proton in jet", "", 400, 0., 100., 20, 0., 1., "p_{T}","Electron percentage", "d^{2}N}/dN^{p_{T}dPercentage");
201     AddHistogram2D<TH2D>("hJetPtVsOthersCount", "Correlation jet pt/amount of other particles in jet", "", 400, 0., 100., 20, 0., 1., "p_{T}","Others percentage", "d^{2}N}/dN^{p_{T}dPercentage");
202     AddHistogram1D<TH1D>("hJetConstituentProtonCount", "Proton count in jets", "", 100, 0., 100., "N protons","dN^{Jets}/dN^{protons}");
203     AddHistogram1D<TH1D>("hJetConstituentPionCount", "Pion count in jets", "", 100, 0., 100., "N pions","dN^{Jets}/dN^{pions}");
204     AddHistogram1D<TH1D>("hJetConstituentKaonCount", "Kaon count in jets", "", 100, 0., 100., "N kaons","dN^{Jets}/dN^{kaons}");
205     AddHistogram1D<TH1D>("hJetConstituentElectronCount", "Electron count in jets", "", 100, 0., 100., "N electrons","dN^{Jets}/dN^{electrons}");
206     AddHistogram1D<TH1D>("hJetConstituentOthersCount", "Others count in jets", "", 100, 0., 100., "N others","dN^{Jets}/dN^{others}");
207
208     AddHistogram2D<TH2D>("hJetPtVsMass_6_14", "Correlation jet pt/summed constituent jet", "", 400, 0., 100., 400, 0., 100., "p_{T}","Mass", "d^{2}N}/dN^{p_{T}dM");
209     AddHistogram2D<TH2D>("hJetPtVsJetMass_6_14", "Correlation jet pt/summed jet", "", 400, 0., 100., 400, 0., 100., "p_{T}","Mass", "d^{2}N}/dN^{p_{T}dM");
210     AddHistogram2D<TH2D>("hJetPtVsProtonCount_6_14", "Correlation jet pt/amount of proton in jet", "", 400, 0., 100., 8, 6., 14., "p_{T}","Proton count", "d^{2}N}/dN^{p_{T}dN");
211     AddHistogram2D<TH2D>("hJetPtVsPionCount_6_14", "Correlation jet pt/amount of pion in jet", "", 400, 0., 100., 8, 6., 14., "p_{T}","Pion count", "d^{2}N}/dN^{p_{T}dN");
212     AddHistogram2D<TH2D>("hJetPtVsKaonCount_6_14", "Correlation jet pt/amount of kaon in jet", "", 400, 0., 100., 8, 6., 14., "p_{T}","Kaon count", "d^{2}N}/dN^{p_{T}dN");
213     AddHistogram2D<TH2D>("hJetPtVsElectronCount_6_14", "Correlation jet pt/amount of proton in jet", "", 400, 0., 100., 8, 6., 14., "p_{T}","Electron count", "d^{2}N}/dN^{p_{T}dN");
214     AddHistogram2D<TH2D>("hJetPtVsOthersCount_6_14", "Correlation jet pt/amount of other particles in jet", "", 400, 0., 100., 8, 6., 14., "p_{T}","Others count", "d^{2}N}/dN^{p_{T}dN");
215     AddHistogram1D<TH1D>("hJetConstituentProtonCount_6_14", "Proton count in jets", "", 100, 0., 100., "N protons","dN^{Jets}/dN^{protons}");
216     AddHistogram1D<TH1D>("hJetConstituentPionCount_6_14", "Pion count in jets", "", 100, 0., 100., "N pions","dN^{Jets}/dN^{pions}");
217     AddHistogram1D<TH1D>("hJetConstituentKaonCount_6_14", "Kaon count in jets", "", 100, 0., 100., "N kaons","dN^{Jets}/dN^{kaons}");
218     AddHistogram1D<TH1D>("hJetConstituentElectronCount_6_14", "Electron count in jets", "", 100, 0., 100., "N electrons","dN^{Jets}/dN^{electrons}");
219     AddHistogram1D<TH1D>("hJetConstituentOthersCount_6_14", "Others count in jets", "", 100, 0., 100., "N others","dN^{Jets}/dN^{others}");
220
221     AddHistogram2D<TH2D>("hJetPtVsMass_2_X", "Correlation jet pt/summed constituent jet", "", 400, 0., 100., 400, 0., 100., "p_{T}","Mass", "d^{2}N}/dN^{p_{T}dM");
222     AddHistogram2D<TH2D>("hJetPtVsJetMass_2_X", "Correlation jet pt/summed jet", "", 400, 0., 100., 400, 0., 100., "p_{T}","Mass", "d^{2}N}/dN^{p_{T}dM");
223     AddHistogram2D<TH2D>("hJetPtVsProtonCount_2_X", "Correlation jet pt/amount of proton in jet", "", 400, 0., 100., 68, 2., 70., "p_{T}","Proton count", "d^{2}N}/dN^{p_{T}dN");
224     AddHistogram2D<TH2D>("hJetPtVsPionCount_2_X", "Correlation jet pt/amount of pion in jet", "", 400, 0., 100., 68, 2., 70., "p_{T}","Pion count", "d^{2}N}/dN^{p_{T}dN");
225     AddHistogram2D<TH2D>("hJetPtVsKaonCount_2_X", "Correlation jet pt/amount of kaon in jet", "", 400, 0., 100., 68, 2., 70., "p_{T}","Kaon count", "d^{2}N}/dN^{p_{T}dN");
226     AddHistogram2D<TH2D>("hJetPtVsElectronCount_2_X", "Correlation jet pt/amount of proton in jet", "", 400, 0., 100., 68, 2., 70., "p_{T}","Electron count", "d^{2}N}/dN^{p_{T}dN");
227     AddHistogram2D<TH2D>("hJetPtVsOthersCount_2_X", "Correlation jet pt/amount of other particles in jet", "", 400, 0., 100., 68, 2., 70., "p_{T}","Others count", "d^{2}N}/dN^{p_{T}dN");
228     AddHistogram1D<TH1D>("hJetConstituentProtonCount_2_X", "Proton count in jets", "", 100, 0., 100., "N protons","dN^{Jets}/dN^{protons}");
229     AddHistogram1D<TH1D>("hJetConstituentPionCount_2_X", "Pion count in jets", "", 100, 0., 100., "N pions","dN^{Jets}/dN^{pions}");
230     AddHistogram1D<TH1D>("hJetConstituentKaonCount_2_X", "Kaon count in jets", "", 100, 0., 100., "N kaons","dN^{Jets}/dN^{kaons}");
231     AddHistogram1D<TH1D>("hJetConstituentElectronCount_2_X", "Electron count in jets", "", 100, 0., 100., "N electrons","dN^{Jets}/dN^{electrons}");
232     AddHistogram1D<TH1D>("hJetConstituentOthersCount_2_X", "Others count in jets", "", 100, 0., 100., "N others","dN^{Jets}/dN^{others}");
233
234     AddHistogram2D<TH2D>("hJetPtVsMass_2_7", "Correlation jet pt/summed constituent jet", "", 400, 0., 100., 400, 0., 100., "p_{T}","Mass", "d^{2}N}/dN^{p_{T}dM");
235     AddHistogram2D<TH2D>("hJetPtVsJetMass_2_7", "Correlation jet pt/summed jet", "", 400, 0., 100., 400, 0., 100., "p_{T}","Mass", "d^{2}N}/dN^{p_{T}dM");
236     AddHistogram2D<TH2D>("hJetPtVsProtonCount_2_7", "Correlation jet pt/amount of proton in jet", "", 400, 0., 100., 6, 2., 8., "p_{T}","Proton count", "d^{2}N}/dN^{p_{T}dN");
237     AddHistogram2D<TH2D>("hJetPtVsPionCount_2_7", "Correlation jet pt/amount of pion in jet", "", 400, 0., 100., 6, 2., 8., "p_{T}","Pion count", "d^{2}N}/dN^{p_{T}dN");
238     AddHistogram2D<TH2D>("hJetPtVsKaonCount_2_7", "Correlation jet pt/amount of kaon in jet", "", 400, 0., 100., 6, 2., 8., "p_{T}","Kaon count", "d^{2}N}/dN^{p_{T}dN");
239     AddHistogram2D<TH2D>("hJetPtVsElectronCount_2_7", "Correlation jet pt/amount of proton in jet", "", 400, 0., 100., 6, 2., 8., "p_{T}","Electron count", "d^{2}N}/dN^{p_{T}dN");
240     AddHistogram2D<TH2D>("hJetPtVsOthersCount_2_7", "Correlation jet pt/amount of other particles in jet", "", 400, 0., 100., 6, 2., 8., "p_{T}","Others count", "d^{2}N}/dN^{p_{T}dN");
241     AddHistogram1D<TH1D>("hJetConstituentProtonCount_2_7", "Proton count in jets", "", 100, 0., 100., "N protons","dN^{Jets}/dN^{protons}");
242     AddHistogram1D<TH1D>("hJetConstituentPionCount_2_7", "Pion count in jets", "", 100, 0., 100., "N pions","dN^{Jets}/dN^{pions}");
243     AddHistogram1D<TH1D>("hJetConstituentKaonCount_2_7", "Kaon count in jets", "", 100, 0., 100., "N kaons","dN^{Jets}/dN^{kaons}");
244     AddHistogram1D<TH1D>("hJetConstituentElectronCount_2_7", "Electron count in jets", "", 100, 0., 100., "N electrons","dN^{Jets}/dN^{electrons}");
245     AddHistogram1D<TH1D>("hJetConstituentOthersCount_2_7", "Others count in jets", "", 100, 0., 100., "N others","dN^{Jets}/dN^{others}");
246
247   }
248
249   // NOTE: Track & Cluster & QA histograms
250   if (fAnalyzeQA)
251   {
252     AddHistogram1D<TH1D>("hVertexX", "X distribution of the vertex", "", 2000, -1., 1., "#Delta x(cm)","dN^{Events}/dx");
253     AddHistogram1D<TH1D>("hVertexY", "Y distribution of the vertex", "", 2000, -1., 1., "#Delta y(cm)","dN^{Events}/dy");
254     AddHistogram2D<TH2D>("hVertexXY", "XY distribution of the vertex", "COLZ", 500, -1., 1., 500, -1., 1.,"#Delta x(cm)", "#Delta y(cm)","dN^{Events}/dxdy");
255     AddHistogram1D<TH1D>("hVertexZBeforeVertexCut", "Z distribution of the vertex (before std. vertex cut)", "", 200, -20., 20., "#Delta z(cm)","dN^{Events}/dz");
256     AddHistogram1D<TH1D>("hVertexZAfterVertexCut", "Z distribution of the vertex (after std. vertex cut)", "", 200, -20., 20., "#Delta z(cm)","dN^{Events}/dz");
257     AddHistogram1D<TH1D>("hVertexR", "R distribution of the vertex", "", 100, 0., 1., "#Delta r(cm)","dN^{Events}/dr");
258     AddHistogram1D<TH1D>("hCentralityV0M", "Centrality distribution V0M", "", fNumberOfCentralityBins, 0., 100., "Centrality","dN^{Events}");
259     AddHistogram1D<TH1D>("hCentralityV0A", "Centrality distribution V0A", "", fNumberOfCentralityBins, 0., 100., "Centrality","dN^{Events}");
260     AddHistogram1D<TH1D>("hCentralityV0C", "Centrality distribution V0C", "", fNumberOfCentralityBins, 0., 100., "Centrality","dN^{Events}");
261     AddHistogram1D<TH1D>("hCentralityZNA", "Centrality distribution ZNA", "", fNumberOfCentralityBins, 0., 100., "Centrality","dN^{Events}");
262     AddHistogram1D<TH1D>("hCentrality", Form("Centrality distribution %s", fCentralityType.Data()), "", fNumberOfCentralityBins, 0., 100., "Centrality","dN^{Events}");
263
264
265     AddHistogram2D<TH2D>("hTrackCountAcc", "Number of tracks in acceptance vs. centrality", "LEGO2", 750, 0., 750., fNumberOfCentralityBins, 0, 100, "N tracks","Centrality", "dN^{Events}/dN^{Tracks}");
266     AddHistogram2D<TH2D>("hTrackPt", "Tracks p_{T} distribution", "", 1000, 0., 250., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)", "Centrality", "dN^{Tracks}/dp_{T}");
267     AddHistogram2D<TH2D>("hTrackPtNegEta", "Tracks p_{T} distribution (negative #eta)", "", 1000, 0., 250., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Tracks}/dp_{T}");
268     AddHistogram2D<TH2D>("hTrackPtPosEta", "Tracks p_{T} distribution (positive #eta)", "", 1000, 0., 250., fNumberOfCentralityBins, 0, 100, "p_{T} (GeV/c)","Centrality","dN^{Tracks}/dp_{T}");
269     AddHistogram1D<TH1D>("hTrackCharge", "Charge", "", 11, -5, 5, "Charge (e)","dN^{Tracks}/dq");
270     AddHistogram1D<TH1D>("hTrackPhi", "Track #phi distribution", "", 360, 0, TMath::TwoPi(), "#phi","dN^{Tracks}/d#phi");
271     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)");
272     AddHistogram2D<TH2D>("hTrackPtPhiEta", "Track p_{T} angular distribution", "LEGO2", 100, 0., 2*TMath::Pi(),100, -2.5, 2.5, "#phi","#eta","dp_{T}^{Tracks}/(d#phi d#eta)");
273
274     AddHistogram2D<TH2D>("hTrackPhiPtCut", "Track #phi distribution for different pT cuts", "LEGO2", 360, 0, TMath::TwoPi(), 20, 0, 20, "#phi", "p_{T} lower cut", "dN^{Tracks}/d#phi dp_{T}");
275     AddHistogram2D<TH2D>("hTrackPhiTrackType", "Track #phi distribution for different track types", "LEGO2", 360, 0, TMath::TwoPi(), 3, 0, 3, "#phi", "Label", "dN^{Tracks}/d#phi");
276     AddHistogram2D<TH2D>("hTrackEta", "Track #eta distribution", "COLZ", 180, fMinEta, fMaxEta, fNumberOfCentralityBins, 0., 100., "#eta", "Centrality", "dN^{Tracks}/d#eta");
277     if (fAnalyzeJets)
278     {
279       // ######## Jet QA
280       AddHistogram1D<TH1D>("hRawJetArea", "Jets area distribution w/o area cut", "", 200, 0., 2., "Area","dN^{Jets}/dA");
281       AddHistogram2D<TH2D>("hJetArea", "Jets area distribution", "COLZ", 200, 0., 2.,  500, -50., 200, "Area","Jet p_{T}","dN^{Jets}/dA");
282       AddHistogram2D<TH2D>("hRawJetPhiEta", "Raw Jets angular distribution w/o #eta cut", "LEGO2", 360, 0., 2*TMath::Pi(),100, -1.0, 1.0, "#phi","#eta","dN^{Jets}/(d#phi d#eta)");
283       AddHistogram2D<TH2D>("hJetEta", "Jets #eta distribution", "COLZ", 180, fMinEta, fMaxEta, fNumberOfCentralityBins, 0., 100., "#eta", "Centrality", "dN^{Jets}/d#eta");
284       AddHistogram2D<TH2D>("hJetEta2GeVTracks", "Jets #eta distribution, track p_{T} > 2 GeV", "COLZ", 180, fMinEta, fMaxEta, fNumberOfCentralityBins, 0., 100., "#eta", "Centrality", "dN^{Jets}/d#eta");
285       AddHistogram2D<TH2D>("hJetEta4GeVTracks", "Jets #eta distribution, track p_{T} > 4 GeV", "COLZ", 180, fMinEta, fMaxEta, fNumberOfCentralityBins, 0., 100., "#eta", "Centrality", "dN^{Jets}/d#eta");
286       AddHistogram2D<TH2D>("hJetPhiEta", "Jets angular distribution", "LEGO2", 360, 0., 2*TMath::Pi(),100, -1.0, 1.0, "#phi","#eta","dN^{Jets}/(d#phi d#eta)");
287       AddHistogram2D<TH2D>("hJetPtPhiEta", "Jets p_{T} angular distribution", "LEGO2", 360, 0., 2*TMath::Pi(),100, -1.0, 1.0, "#phi","#eta","dp_{T}^{Jets}/(d#phi d#eta)");
288       AddHistogram2D<TH2D>("hJetPtVsConstituentCount", "Jets number of constituents vs. jet p_{T}", "COLZ", 400, 0., 200., 100, 0., 100., "p_{T}","N^{Tracks}","dN^{Jets}/(dp_{T} dN^{tracks})");
289     }
290   }
291
292   // register Histograms
293   for (Int_t i = 0; i < fHistCount; i++)
294   {
295     fOutputList->Add(fHistList->At(i));
296   }
297   
298   PostData(1,fOutputList); // important for merging
299
300 }
301
302 //________________________________________________________________________
303 AliAnalysisTaskChargedJetsPA::AliAnalysisTaskChargedJetsPA(const char *name, const char* trackArrayName, const char* jetArrayName, const char* backgroundJetArrayName) : AliAnalysisTaskSE(name), fOutputList(0), fAnalyzeJets(1), fAnalyzeJetProfile(1), fAnalyzeQA(1), fAnalyzeBackground(1), fAnalyzeDeprecatedBackgrounds(1), fAnalyzePythia(0), fAnalyzeMassCorrelation(0), fHasTracks(0), fHasJets(0), fHasBackgroundJets(0), fIsKinematics(0), fUseDefaultVertexCut(1), fUsePileUpCut(1), fSetCentralityToOne(0), fNoExternalBackground(0), fBackgroundForJetProfile(0), fPartialAnalysisNParts(1), fPartialAnalysisIndex(0), fJetArray(0), fTrackArray(0), fBackgroundJetArray(0), fJetArrayName(0), fTrackArrayName(0), fBackgroundJetArrayName(0), fNumPtHardBins(11), fUsePtHardBin(-1), fRhoTaskName(), fNcoll(6.88348), fRandConeRadius(0.4), fSignalJetRadius(0.4), fBackgroundJetRadius(0.4), fTRBackgroundConeRadius(0.6), fNumberRandCones(8), fNumberExcludedJets(-1), fDijetMaxAngleDeviation(10.0), fBackgroundJetEtaWindow(0.5), fMinEta(-0.9), fMaxEta(0.9), fMinJetEta(-0.5), fMaxJetEta(0.5), fMinTrackPt(0.150), fMinJetPt(0.15), fMinJetArea(0.5), fMinBackgroundJetPt(0.0), fMinDijetLeadingPt(10.0), fNumberOfCentralityBins(20), fCentralityType("V0A"), fFirstLeadingJet(0), fSecondLeadingJet(0), fNumberSignalJets(0), fCrossSection(0.0), fTrials(0.0), fRandom(0), fHelperClass(0), fInitialized(0), fTaskInstanceCounter(0), fHistList(0), fHistCount(0), fIsDEBUG(0), fEventCounter(0)
304 {
305   #ifdef DEBUGMODE
306     AliInfo("Calling constructor.");
307   #endif
308
309   // Every instance of this task gets his own number
310   static Int_t instance = 0;
311   fTaskInstanceCounter = instance;
312   instance++;
313
314   fTrackArrayName = new TString(trackArrayName);
315   if (fTrackArrayName->Contains("MCParticles") || fTrackArrayName->Contains("mcparticles"))
316     fIsKinematics = kTRUE;
317
318   fJetArrayName = new TString(jetArrayName);
319   if (strcmp(fJetArrayName->Data(),"") == 0)
320   {
321     fAnalyzeJets = kFALSE;
322     fAnalyzeJetProfile = kFALSE;
323   }
324   else
325     fAnalyzeJets = kTRUE;
326     
327   fBackgroundJetArrayName = new TString(backgroundJetArrayName);
328   if (strcmp(fBackgroundJetArrayName->Data(),"") == 0)
329     fAnalyzeBackground = kFALSE;
330   else
331     fAnalyzeBackground = kTRUE;
332
333   DefineOutput(1, TList::Class());
334  
335   fHistList = new TList();
336
337   for(Int_t i=0;i<1024;i++)
338     fSignalJets[i] = NULL;
339
340
341   #ifdef DEBUGMODE
342     AliInfo("Constructor done.");
343   #endif
344   
345 }
346
347 //________________________________________________________________________
348 inline Double_t AliAnalysisTaskChargedJetsPA::GetConePt(Double_t eta, Double_t phi, Double_t radius)
349 {
350   Double_t tmpConePt = 0.0;
351
352   for (Int_t i = 0; i < fTrackArray->GetEntries(); i++)
353   {
354     AliVTrack* tmpTrack = static_cast<AliVTrack*>(fTrackArray->At(i));
355     if (IsTrackInAcceptance(tmpTrack))
356       if(IsTrackInCone(tmpTrack, eta, phi, radius))
357         tmpConePt = tmpConePt + tmpTrack->Pt();
358   }
359   return tmpConePt;
360 }
361
362 //________________________________________________________________________
363 inline Double_t AliAnalysisTaskChargedJetsPA::GetCorrectedConePt(Double_t eta, Double_t phi, Double_t radius, Double_t background)
364 {
365   Double_t tmpConePt = 0.0;
366
367   for (Int_t i = 0; i < fTrackArray->GetEntries(); i++)
368   {
369     AliVTrack* tmpTrack = static_cast<AliVTrack*>(fTrackArray->At(i));
370     if (IsTrackInAcceptance(tmpTrack))
371       if(IsTrackInCone(tmpTrack, eta, phi, radius))
372         tmpConePt = tmpConePt + tmpTrack->Pt();
373   }
374   Double_t realConeArea = (2.0*(fMaxEta-fMinEta)) * TMath::TwoPi() * MCGetOverlapCircleRectancle(eta, phi, radius, fMinEta, fMaxEta, 0., TMath::TwoPi());
375   tmpConePt -= background * realConeArea; // subtract background
376
377   return tmpConePt;
378 }
379
380
381 //________________________________________________________________________
382 inline Double_t AliAnalysisTaskChargedJetsPA::GetPtHard()
383 {
384   #ifdef DEBUGMODE
385     AliInfo("Starting GetPtHard.");
386   #endif
387   AliGenPythiaEventHeader* pythiaHeader = dynamic_cast<AliGenPythiaEventHeader*>(MCEvent()->GenEventHeader());
388   if (MCEvent()) 
389     if (!pythiaHeader)
390     {
391       // Check if AOD
392       AliAODMCHeader* aodMCH = dynamic_cast<AliAODMCHeader*>(InputEvent()->FindListObject(AliAODMCHeader::StdBranchName()));
393
394       if (aodMCH)
395       {
396         for(UInt_t i = 0;i<aodMCH->GetNCocktailHeaders();i++)
397         {
398           pythiaHeader = dynamic_cast<AliGenPythiaEventHeader*>(aodMCH->GetCocktailHeader(i));
399           if (pythiaHeader) break;
400         }
401       }
402     }
403
404   #ifdef DEBUGMODE
405     AliInfo("Ending GetPtHard.");
406   #endif
407   if (pythiaHeader)
408     return pythiaHeader->GetPtHard();
409
410   AliWarning(Form("In task %s: GetPtHard() failed!", GetName()));
411   return -1.0;
412 }
413
414
415 //________________________________________________________________________
416 inline Double_t AliAnalysisTaskChargedJetsPA::GetPythiaTrials()
417 {
418   #ifdef DEBUGMODE
419     AliInfo("Starting GetPythiaTrials.");
420   #endif
421   AliGenPythiaEventHeader* pythiaHeader = dynamic_cast<AliGenPythiaEventHeader*>(MCEvent()->GenEventHeader());
422   if (MCEvent()) 
423     if (!pythiaHeader)
424     {
425       // Check if AOD
426       AliAODMCHeader* aodMCH = dynamic_cast<AliAODMCHeader*>(InputEvent()->FindListObject(AliAODMCHeader::StdBranchName()));
427
428       if (aodMCH)
429       {
430         for(UInt_t i = 0;i<aodMCH->GetNCocktailHeaders();i++)
431         {
432           pythiaHeader = dynamic_cast<AliGenPythiaEventHeader*>(aodMCH->GetCocktailHeader(i));
433           if (pythiaHeader) break;
434         }
435       }
436     }
437
438   #ifdef DEBUGMODE
439     AliInfo("Ending GetPythiaTrials.");
440   #endif
441   if (pythiaHeader)
442     return pythiaHeader->Trials();
443
444   AliWarning(Form("In task %s: GetPythiaTrials() failed!", GetName()));
445   return -1.0;
446 }
447
448
449
450 //________________________________________________________________________
451 inline Int_t AliAnalysisTaskChargedJetsPA::GetPtHardBin()
452 {
453   #ifdef DEBUGMODE
454     AliInfo("Starting GetPtHardBin.");
455   #endif
456   // ########## PT HARD BIN EDGES
457   const Int_t kPtHardLowerEdges[] =  { 0, 5,11,21,36,57, 84,117,152,191,234};
458   const Int_t kPtHardHigherEdges[] = { 5,11,21,36,57,84,117,152,191,234,1000000};
459
460   Int_t tmpPtHardBin = 0;
461   Double_t tmpPtHard = GetPtHard();
462  
463   for (tmpPtHardBin = 0; tmpPtHardBin <= fNumPtHardBins; tmpPtHardBin++)
464     if (tmpPtHard >= kPtHardLowerEdges[tmpPtHardBin] && tmpPtHard < kPtHardHigherEdges[tmpPtHardBin])
465       break;
466
467   #ifdef DEBUGMODE
468     AliInfo("Ending GetPtHardBin.");
469   #endif
470   return tmpPtHardBin;
471 }
472
473 //________________________________________________________________________
474 Double_t AliAnalysisTaskChargedJetsPA::GetExternalRho()
475 {
476   // Get rho from event.
477   AliRhoParameter *rho = 0;
478   if (!fRhoTaskName.IsNull()) {
479     rho = dynamic_cast<AliRhoParameter*>(InputEvent()->FindListObject(fRhoTaskName.Data()));
480     if (!rho) {
481       AliWarning(Form("%s: Could not retrieve rho with name %s!", GetName(), fRhoTaskName.Data())); 
482       return 0;
483     }
484   }
485   else
486     return 0;
487
488   return (rho->GetVal());
489 }
490
491
492 //________________________________________________________________________
493 inline Bool_t AliAnalysisTaskChargedJetsPA::IsTrackInCone(AliVTrack* track, Double_t eta, Double_t phi, Double_t radius)
494 {
495   // This is to use a full cone in phi even at the edges of phi (2pi -> 0) (0 -> 2pi)
496   Double_t trackPhi = 0.0;
497   if (track->Phi() > (TMath::TwoPi() - (radius-phi)))
498     trackPhi = track->Phi() - TMath::TwoPi();
499   else if (track->Phi() < (phi+radius - TMath::TwoPi()))
500     trackPhi = track->Phi() + TMath::TwoPi();
501   else
502     trackPhi = track->Phi();
503   
504   if ( TMath::Abs(trackPhi-phi)*TMath::Abs(trackPhi-phi) + TMath::Abs(track->Eta()-eta)*TMath::Abs(track->Eta()-eta) <= radius*radius)
505     return kTRUE;
506   
507   return kFALSE;
508 }
509
510 //________________________________________________________________________
511 inline Bool_t AliAnalysisTaskChargedJetsPA::IsTrackInJet(AliEmcalJet* jet, Int_t trackIndex)
512 {
513   for (Int_t i = 0; i < jet->GetNumberOfTracks(); ++i)
514   {
515     Int_t jetTrack = jet->TrackAt(i);
516     if (jetTrack == trackIndex)
517       return kTRUE;
518   }
519   return kFALSE;
520 }
521
522 //________________________________________________________________________
523 inline Bool_t AliAnalysisTaskChargedJetsPA::IsJetOverlapping(AliEmcalJet* jet1, AliEmcalJet* jet2)
524 {
525   for (Int_t i = 0; i < jet1->GetNumberOfTracks(); ++i)
526   {
527     Int_t jet1Track = jet1->TrackAt(i);
528     for (Int_t j = 0; j < jet2->GetNumberOfTracks(); ++j)
529     {
530       Int_t jet2Track = jet2->TrackAt(j);
531       if (jet1Track == jet2Track)
532         return kTRUE;
533     }
534   }
535   return kFALSE;
536 }
537
538 //________________________________________________________________________
539 inline Bool_t AliAnalysisTaskChargedJetsPA::IsEventInAcceptance(AliVEvent* event)
540 {
541   if (!event)
542     return kFALSE;
543
544   FillHistogram("hEventAcceptance", 0.5); // number of events before manual cuts
545   if(fUsePileUpCut)
546     if(fHelperClass->IsPileUpEvent(event))
547       return kFALSE;
548
549   FillHistogram("hEventAcceptance", 1.5); // number of events after pileup cuts
550
551   if(fAnalyzeQA)
552     FillHistogram("hVertexZBeforeVertexCut",event->GetPrimaryVertex()->GetZ());
553
554   if(fUseDefaultVertexCut)
555   {
556     if(!fHelperClass->IsVertexSelected2013pA(event))
557       return kFALSE;
558   }
559   else // Failsafe vertex cut
560   {
561     if(!event->GetPrimaryVertex() || (TMath::Abs(event->GetPrimaryVertex()->GetZ()) > 10.0) || (event->GetPrimaryVertex()->GetNContributors()<1)) 
562       return kFALSE;
563   }
564
565   if(fAnalyzeQA)
566     FillHistogram("hVertexZAfterVertexCut",event->GetPrimaryVertex()->GetZ());
567
568   FillHistogram("hEventAcceptance", 2.5); // number of events after vertex cut
569
570   return kTRUE;
571 }
572
573 //________________________________________________________________________
574 inline Bool_t AliAnalysisTaskChargedJetsPA::IsTrackInAcceptance(AliVParticle* track)
575 {
576   FillHistogram("hTrackAcceptance", 0.5);
577   if (track != 0)
578   {
579     if(fIsKinematics)
580     {
581       // TODO: Only working for AOD MC
582       if((!track->Charge()) || (!(static_cast<AliAODMCParticle*>(track))->IsPhysicalPrimary()) )
583         return kFALSE;
584     }
585     if ((track->Eta() < fMaxEta) && (track->Eta() >= fMinEta))
586     {
587       FillHistogram("hTrackAcceptance", 1.5);
588       if (track->Pt() >= fMinTrackPt)
589       {
590         FillHistogram("hTrackAcceptance", 2.5);
591         return kTRUE;
592       }
593     }
594   }
595   return kFALSE;
596 }
597
598 //________________________________________________________________________
599 inline Bool_t AliAnalysisTaskChargedJetsPA::IsBackgroundJetInAcceptance(AliEmcalJet *jet)
600 {   
601   if (jet != 0)
602     if (TMath::Abs(jet->Eta()) <= fBackgroundJetEtaWindow)
603       if (jet->Pt() >= fMinBackgroundJetPt)
604         return kTRUE;
605
606   return kFALSE;
607 }
608
609 //________________________________________________________________________
610 inline Bool_t AliAnalysisTaskChargedJetsPA::IsSignalJetInAcceptance(AliEmcalJet *jet)
611 {   
612   FillHistogram("hJetAcceptance", 0.5);
613   if (jet != 0)
614     if ((jet->Eta() >= fMinJetEta) && (jet->Eta() < fMaxJetEta))
615     {
616       FillHistogram("hJetAcceptance", 1.5);
617       if (jet->Pt() >= fMinJetPt)
618       {
619         FillHistogram("hJetAcceptance", 2.5);
620         if (jet->Area() >= fMinJetArea)
621         {
622           FillHistogram("hJetAcceptance", 3.5);
623           return kTRUE;
624         }
625       }
626     }
627   return kFALSE;
628 }
629
630 //________________________________________________________________________
631 inline Bool_t AliAnalysisTaskChargedJetsPA::IsDijet(AliEmcalJet *jet1, AliEmcalJet *jet2)
632 {   
633   // Output from GetDeltaPhi is < pi in any case
634   if ((jet1 != 0) && (jet2 != 0))
635     if((TMath::Pi() - GetDeltaPhi(jet1->Phi(),jet2->Phi())) < fDijetMaxAngleDeviation)
636       if ((jet1->Pt() > fMinDijetLeadingPt) && (jet2->Pt() > fMinDijetLeadingPt)) //TODO: Introduce recoil jet?
637         return kTRUE;
638
639   return kFALSE;
640 }
641
642 //________________________________________________________________________
643 void AliAnalysisTaskChargedJetsPA::ExecOnce()
644 {
645   #ifdef DEBUGMODE
646     AliInfo("Starting ExecOnce.");
647   #endif
648   fInitialized = kTRUE;
649
650   // Check for track array
651   if (strcmp(fTrackArrayName->Data(), "") != 0)
652   {
653     fTrackArray = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fTrackArrayName->Data()));
654     fHasTracks = kTRUE;
655     if (!fTrackArray) 
656     {
657       AliWarning(Form("%s: Could not retrieve tracks %s! This is OK, if tracks are not demanded.", GetName(), fTrackArrayName->Data())); 
658       fHasTracks = kFALSE;
659     } 
660     else
661     {
662       TClass *cl = fTrackArray->GetClass();
663       if (!cl->GetBaseClass("AliVParticle"))
664       {
665         AliError(Form("%s: Collection %s does not contain AliVParticle objects!", GetName(), fTrackArrayName->Data())); 
666         fTrackArray = 0;
667         fHasTracks = kFALSE;
668       }
669     }
670   }
671
672   // Check for jet array
673   if (strcmp(fJetArrayName->Data(), "") != 0)
674   {
675     fJetArray = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fJetArrayName->Data()));
676     fHasJets = kTRUE;
677
678     if (!fJetArray) 
679     {
680       AliWarning(Form("%s: Could not retrieve jets %s! This is OK, if jets are not demanded.", GetName(), fJetArrayName->Data())); 
681       fHasJets = kFALSE;
682     } 
683     else
684     {
685       if (!fJetArray->GetClass()->GetBaseClass("AliEmcalJet")) 
686       {
687         AliError(Form("%s: Collection %s does not contain AliEmcalJet objects!", GetName(), fJetArrayName->Data())); 
688         fJetArray = 0;
689         fHasJets = kFALSE;
690       }
691     }
692   }
693
694   // Check for background object
695   if (strcmp(fBackgroundJetArrayName->Data(), "") != 0)
696   {
697     fBackgroundJetArray = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fBackgroundJetArrayName->Data()));
698     fHasBackgroundJets = kTRUE;
699     if (!fBackgroundJetArray)
700     {
701       AliInfo(Form("%s: Could not retrieve background jets %s! This is OK, if background is not demanded.", GetName(), fBackgroundJetArrayName->Data())); 
702       fHasBackgroundJets = kFALSE;
703     }
704   }
705
706   // Look, if initialization is OK
707   if (!fHasTracks && fAnalyzeBackground)
708   {
709     AliError(Form("%s: Tracks NOT successfully casted although demanded! Deactivating background analysis",GetName()));
710     fAnalyzeBackground = kFALSE;
711   }
712   if ((!fHasJets && fAnalyzeJets) || (!fHasJets && fAnalyzeBackground))
713   {
714     AliError(Form("%s: Jets NOT successfully casted although demanded!  Deactivating jet- and background analysis",GetName()));
715     fAnalyzeJets = kFALSE;
716     fAnalyzeBackground = kFALSE;
717   }
718   if (!fHasBackgroundJets && fAnalyzeBackground)
719   {
720     AliError(Form("%s: Background NOT successfully casted although demanded!  Deactivating background analysis",GetName()));
721     fAnalyzeBackground = kFALSE;
722   }
723
724   // Initialize helper class (for vertex selection & pile up correction)
725   fHelperClass = new AliAnalysisUtils();
726   fHelperClass->SetCutOnZVertexSPD(kFALSE);
727   // Histogram init
728   Init();
729
730   #ifdef DEBUGMODE
731     AliInfo("ExecOnce done.");
732   #endif
733
734 }
735
736 //________________________________________________________________________
737 void AliAnalysisTaskChargedJetsPA::GetSignalJets()
738 {
739   // Reset vars
740   fFirstLeadingJet = NULL;
741   fSecondLeadingJet = NULL;
742   fNumberSignalJets = 0;
743
744   TList tmpJets;
745   for (Int_t i = 0; i < fJetArray->GetEntries(); i++)
746   {
747     AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJetArray->At(i));
748     if (!jet)
749     {
750       AliError(Form("%s: Could not receive jet %d", GetName(), i));
751       continue;
752     }
753     if (!IsSignalJetInAcceptance(jet))
754       continue;
755
756     for (Int_t j = 0; j <= tmpJets.GetEntries(); j++)
757     {
758       if (j>tmpJets.GetEntries()-1) // When passed last item add the jet at the end
759       {
760         tmpJets.Add(jet);
761         break;
762       }
763
764       AliEmcalJet* listJet = static_cast<AliEmcalJet*>(tmpJets.At(j));
765      
766       if(jet->Pt() < listJet->Pt()) // Insert jet before that one in list if pt smaller
767       {
768         tmpJets.AddAt(jet, j);
769         break;
770       }
771     }
772   }
773
774   for (Int_t i = 0; i < tmpJets.GetEntries(); i++)
775   {
776     AliEmcalJet* jet = static_cast<AliEmcalJet*>(tmpJets.At(i));
777     fSignalJets[fNumberSignalJets] = jet;
778     fNumberSignalJets++;
779   }
780
781   Int_t leadingJets[]   = {-1, -1};
782   GetLeadingJets(fJetArray, &leadingJets[0], kTRUE);
783   
784   if (leadingJets[0] > -1)
785     fFirstLeadingJet  = static_cast<AliEmcalJet*>(fJetArray->At(leadingJets[0]));
786   if (leadingJets[1] > -1)
787     fSecondLeadingJet = static_cast<AliEmcalJet*>(fJetArray->At(leadingJets[1]));
788 }
789
790 //________________________________________________________________________
791 Int_t AliAnalysisTaskChargedJetsPA::GetLeadingJets(TClonesArray* jetArray, Int_t* jetIDArray, Bool_t isSignalJets)
792 {
793 // Writes first two leading jets into already registered array jetIDArray
794
795   if (!jetArray)
796   {
797     AliError("Could not get the jet array to get leading jets from it!");
798     return 0;
799   }
800
801   Float_t maxJetPts[] = {0, 0};
802   jetIDArray[0] = -1;
803   jetIDArray[1] = -1;
804
805   Int_t jetCount = jetArray->GetEntries();
806   Int_t jetCountAccepted = 0;
807
808   for (Int_t i = 0; i < jetCount; i++)
809   {
810     AliEmcalJet* jet = static_cast<AliEmcalJet*>(jetArray->At(i));
811     if (!jet) 
812     {
813       AliError(Form("%s: Could not receive jet %d", GetName(), i));
814       continue;
815     }
816
817     if(isSignalJets)
818     {
819       if (!IsSignalJetInAcceptance(jet)) continue;
820     }
821     else
822     {
823       if (!IsBackgroundJetInAcceptance(jet)) continue;
824     }    
825
826     if (jet->Pt() > maxJetPts[0]) 
827     {
828       maxJetPts[1] = maxJetPts[0];
829       jetIDArray[1] = jetIDArray[0];
830       maxJetPts[0] = jet->Pt();
831       jetIDArray[0] = i;
832     }
833     else if (jet->Pt() > maxJetPts[1]) 
834     {
835       maxJetPts[1] = jet->Pt();
836       jetIDArray[1] = i;
837     }
838     jetCountAccepted++;
839   }
840   return jetCountAccepted;
841 }
842
843
844 //________________________________________________________________________
845 Double_t AliAnalysisTaskChargedJetsPA::GetCorrectedJetPt(AliEmcalJet* jet, Double_t background)
846 {
847   #ifdef DEBUGMODE
848     AliInfo("Getting corrected jet spectra.");
849   #endif
850
851   if(!jet)
852   {
853     AliError("Jet pointer passed to GetCorrectedJetPt() not valid!");
854     return -1.0;
855   }
856
857   Double_t correctedPt = -1.0;
858
859   // if the passed background is not valid, do not subtract it
860   if(background < 0)
861     background = 0;
862
863   // Subtract background
864   correctedPt = jet->Pt() - background * jet->Area();
865
866   #ifdef DEBUGMODE
867     AliInfo("Got corrected jet spectra.");
868   #endif 
869
870   return correctedPt;
871 }
872
873
874
875 //________________________________________________________________________
876 Double_t AliAnalysisTaskChargedJetsPA::GetDeltaPt(Double_t rho, Double_t leadingJetExclusionProbability)
877 {
878   #ifdef DEBUGMODE
879     AliInfo("Getting Delta Pt.");
880   #endif
881
882   // Define an invalid delta pt
883   Double_t deltaPt = -10000.0;
884
885   // Define eta range
886   Double_t etaMin, etaMax;
887   etaMin = fMinEta+fRandConeRadius;
888   etaMax = fMaxEta-fRandConeRadius;
889
890   // Define random cone
891   Bool_t coneValid = kTRUE;
892   Double_t tmpRandConeEta = etaMin + fRandom->Rndm()*(etaMax-etaMin);
893   Double_t tmpRandConePhi = fRandom->Rndm()*TMath::TwoPi();
894
895   // if there is a jet, check for overlap if demanded
896   if(leadingJetExclusionProbability)
897   {
898     AliEmcalJet* tmpLeading = dynamic_cast<AliEmcalJet*>(fJetArray->At(0));
899     // Get leading jet (regardless of pT)
900     for (Int_t i = 1; i<fJetArray->GetEntries(); i++)
901     {
902       AliEmcalJet* tmpJet = static_cast<AliEmcalJet*>(fJetArray->At(i));
903       // if jet is in acceptance and higher, take as new leading
904       if (tmpJet)
905         if ( ((tmpJet->Eta() >= fMinJetEta) && (tmpJet->Eta() < fMaxJetEta)) && (tmpJet->Area() >= fMinJetArea))
906           if((!tmpLeading) || (tmpJet->Pt() > tmpLeading->Pt()))
907             tmpLeading = tmpJet;
908     }
909     if(tmpLeading)
910     {
911       Double_t excludedJetPhi = tmpLeading->Phi();
912       Double_t excludedJetEta = tmpLeading->Eta();
913       Double_t tmpDeltaPhi = GetDeltaPhi(tmpRandConePhi, excludedJetPhi);
914
915       // Check, if cone has overlap with jet
916       if ( tmpDeltaPhi*tmpDeltaPhi + TMath::Abs(tmpRandConeEta-excludedJetEta)*TMath::Abs(tmpRandConeEta-excludedJetEta) <= fRandConeRadius*fRandConeRadius)
917       {
918         // Define probability to exclude the RC
919         Double_t probability = leadingJetExclusionProbability;
920
921         // Only exclude cone with a given probability
922         if (fRandom->Rndm()<=probability)
923           coneValid = kFALSE;
924       }
925     }
926   }
927
928
929   // Get the cones' pt and calculate delta pt
930   if (coneValid)
931     deltaPt = GetConePt(tmpRandConeEta,tmpRandConePhi,fRandConeRadius) - (rho*fRandConeRadius*fRandConeRadius*TMath::Pi());
932
933   return deltaPt;
934   #ifdef DEBUGMODE
935     AliInfo("Got Delta Pt.");
936   #endif
937 }
938
939 //________________________________________________________________________
940 void AliAnalysisTaskChargedJetsPA::GetKTBackgroundDensityAll(Int_t numberExcludeLeadingJets, Double_t& rhoPbPb, Double_t& rhoPbPbWithGhosts, Double_t& rhoCMS, Double_t& rhoImprovedCMS, Double_t& rhoMean, Double_t& rhoTrackLike)
941 {
942   #ifdef DEBUGMODE
943     AliInfo("Getting ALL KT background density.");
944   #endif
945
946   static Double_t tmpRhoPbPb[1024];
947   static Double_t tmpRhoPbPbWithGhosts[1024];
948   static Double_t tmpRhoMean[1024];
949   static Double_t tmpRhoCMS[1024];
950   static Double_t tmpRhoImprovedCMS[1024];
951   Double_t tmpCoveredArea = 0.0;
952   Double_t tmpSummedArea = 0.0;
953   Double_t tmpPtTrackLike = 0.0;
954   Double_t tmpAreaTrackLike = 0.0;
955
956   // Setting invalid values
957   rhoPbPb = 0.0;
958   rhoPbPbWithGhosts = 0.0;
959   rhoCMS = 0.0;
960   rhoImprovedCMS = 0.0;
961   rhoMean = 0.0;
962   rhoTrackLike = 0.0;
963
964   Int_t rhoPbPbJetCount = 0;
965   Int_t rhoPbPbWithGhostsJetCount = 0;
966   Int_t rhoCMSJetCount = 0;
967   Int_t rhoImprovedCMSJetCount = 0;
968   Int_t rhoMeanJetCount = 0;
969
970
971   // Find 2 leading KT jets for the original PbPb approach
972   Int_t leadingKTJets[]   = {-1, -1};
973   GetLeadingJets(fBackgroundJetArray, &leadingKTJets[0], kFALSE);
974
975   // Exclude UP TO numberExcludeLeadingJets
976   if(numberExcludeLeadingJets==-1)
977     numberExcludeLeadingJets = fNumberSignalJets;
978   if (fNumberSignalJets < numberExcludeLeadingJets)
979     numberExcludeLeadingJets = fNumberSignalJets;
980
981   for (Int_t i = 0; i < fBackgroundJetArray->GetEntries(); i++)
982   {
983     AliEmcalJet* backgroundJet = static_cast<AliEmcalJet*>(fBackgroundJetArray->At(i));
984
985     if (!backgroundJet)
986     {
987       AliError(Form("%s: Could not receive jet %d", GetName(), i));
988       continue;
989     } 
990
991     tmpSummedArea += backgroundJet->Area();
992     if(backgroundJet->Pt() > 0.150)
993       tmpCoveredArea += backgroundJet->Area();
994
995     if (!IsBackgroundJetInAcceptance(backgroundJet))
996       continue;
997
998     // Search for overlap with signal jets
999     Bool_t isOverlapping = kFALSE;
1000     for(Int_t j=0;j<numberExcludeLeadingJets;j++)
1001     {
1002       AliEmcalJet* signalJet = fSignalJets[j];
1003      
1004       if(signalJet->Pt() >= 5.0)
1005         if(IsJetOverlapping(signalJet, backgroundJet))
1006         {
1007           isOverlapping = kTRUE;
1008           break;
1009         }
1010     }
1011
1012     Double_t tmpRho = 0.0;
1013     if(backgroundJet->Area())
1014       tmpRho = backgroundJet->Pt() / backgroundJet->Area();
1015
1016     // PbPb approach (take ghosts into account)
1017     if ((i != leadingKTJets[0]) && (i != leadingKTJets[1]))
1018     {
1019       tmpRhoPbPbWithGhosts[rhoPbPbWithGhostsJetCount] = tmpRho;
1020       rhoPbPbWithGhostsJetCount++;
1021     }
1022
1023     if(backgroundJet->Pt() > 0.150)
1024     {
1025       // CMS approach: don't take ghosts into acount
1026       tmpRhoCMS[rhoCMSJetCount] = tmpRho;
1027       rhoCMSJetCount++;
1028
1029       // Improved CMS approach: like CMS but excluding signal
1030       if(!isOverlapping)
1031       {
1032         tmpRhoImprovedCMS[rhoImprovedCMSJetCount] = tmpRho;
1033         rhoImprovedCMSJetCount++;
1034       }
1035
1036       // PbPb w/o ghosts approach (just neglect ghosts)
1037       if ((i != leadingKTJets[0]) && (i != leadingKTJets[1]))
1038       {  
1039         tmpRhoPbPb[rhoPbPbJetCount] = tmpRho;
1040         rhoPbPbJetCount++;
1041       }
1042     }
1043
1044     // (no overlap with signal jets)
1045     if(!isOverlapping)
1046     {
1047       // Mean approach
1048       tmpRhoMean[rhoMeanJetCount] = tmpRho;
1049       rhoMeanJetCount++;
1050       
1051       // Track like approach approach
1052       tmpPtTrackLike += backgroundJet->Pt();
1053       tmpAreaTrackLike += backgroundJet->Area();
1054     }
1055
1056   }
1057
1058   if (tmpAreaTrackLike > 0)
1059     rhoTrackLike = tmpPtTrackLike/tmpAreaTrackLike;
1060   if (rhoPbPbJetCount > 0)
1061     rhoPbPb = TMath::Median(rhoPbPbJetCount, tmpRhoPbPb);
1062   if (rhoPbPbWithGhostsJetCount > 0)
1063     rhoPbPbWithGhosts = TMath::Median(rhoPbPbWithGhostsJetCount, tmpRhoPbPbWithGhosts);
1064   if (rhoCMSJetCount > 0)
1065     rhoCMS = TMath::Median(rhoCMSJetCount, tmpRhoCMS) * tmpCoveredArea/tmpSummedArea;
1066   if (rhoImprovedCMSJetCount > 0)
1067   {
1068     rhoImprovedCMS = TMath::Median(rhoImprovedCMSJetCount, tmpRhoImprovedCMS) * tmpCoveredArea/tmpSummedArea;
1069   }
1070   if (rhoMeanJetCount > 0)
1071     rhoMean = TMath::Mean(rhoMeanJetCount, tmpRhoMean);
1072
1073   #ifdef DEBUGMODE
1074     AliInfo("Got ALL KT background density.");
1075   #endif
1076 }
1077
1078 //________________________________________________________________________
1079 void AliAnalysisTaskChargedJetsPA::GetKTBackgroundDensity(Int_t numberExcludeLeadingJets, Double_t& rhoImprovedCMS)
1080 {
1081   #ifdef DEBUGMODE
1082     AliInfo("Getting KT background density.");
1083   #endif
1084
1085   static Double_t tmpRhoImprovedCMS[1024];
1086   Double_t tmpCoveredArea = 0.0;
1087   Double_t tmpSummedArea = 0.0;
1088
1089   // Setting invalid values
1090   rhoImprovedCMS = 0.0;
1091
1092   Int_t rhoImprovedCMSJetCount = 0;
1093
1094   // Exclude UP TO numberExcludeLeadingJets
1095   if(numberExcludeLeadingJets==-1)
1096     numberExcludeLeadingJets = fNumberSignalJets;
1097   if (fNumberSignalJets < numberExcludeLeadingJets)
1098     numberExcludeLeadingJets = fNumberSignalJets;
1099
1100   for (Int_t i = 0; i < fBackgroundJetArray->GetEntries(); i++)
1101   {
1102     AliEmcalJet* backgroundJet = static_cast<AliEmcalJet*>(fBackgroundJetArray->At(i));
1103
1104     if (!backgroundJet)
1105     {
1106       AliError(Form("%s: Could not receive jet %d", GetName(), i));
1107       continue;
1108     } 
1109
1110     // Search for overlap with signal jets
1111     Bool_t isOverlapping = kFALSE;
1112     for(Int_t j=0;j<numberExcludeLeadingJets;j++)
1113     {
1114       AliEmcalJet* signalJet = fSignalJets[j];
1115       if(signalJet->Pt() >= 5.0)     
1116         if(IsJetOverlapping(signalJet, backgroundJet))
1117         {
1118           isOverlapping = kTRUE;
1119           break;
1120         }
1121     }
1122
1123     tmpSummedArea += backgroundJet->Area();
1124     if(backgroundJet->Pt() > 0.150)
1125       tmpCoveredArea += backgroundJet->Area();
1126
1127     if (!IsBackgroundJetInAcceptance(backgroundJet))
1128       continue;
1129
1130     Double_t tmpRho = backgroundJet->Pt() / backgroundJet->Area();
1131
1132     if(backgroundJet->Pt() > 0.150)
1133       if(!isOverlapping)
1134       {
1135         tmpRhoImprovedCMS[rhoImprovedCMSJetCount] = tmpRho;
1136         rhoImprovedCMSJetCount++;
1137       }
1138   }
1139
1140   if (rhoImprovedCMSJetCount > 0)
1141   {
1142     rhoImprovedCMS = TMath::Median(rhoImprovedCMSJetCount, tmpRhoImprovedCMS) * tmpCoveredArea/tmpSummedArea;
1143   }
1144   #ifdef DEBUGMODE
1145     AliInfo("Got KT background density.");
1146   #endif
1147 }
1148
1149
1150 //________________________________________________________________________
1151 void AliAnalysisTaskChargedJetsPA::GetTRBackgroundDensity(Int_t numberExcludeLeadingJets, Double_t& rhoNoExclusion, Double_t& rhoConeExclusion02, Double_t& rhoConeExclusion04, Double_t& rhoConeExclusion06, Double_t& rhoConeExclusion08, Double_t& rhoExactExclusion)
1152 {
1153   #ifdef DEBUGMODE
1154     AliInfo("Getting TR background density.");
1155   #endif
1156
1157   Double_t summedTracksPtCone04 = 0.0;
1158   Double_t summedTracksPtCone02 = 0.0;
1159   Double_t summedTracksPtCone06 = 0.0;
1160   Double_t summedTracksPtCone08 = 0.0;
1161   Double_t summedTracksPtWithinJets = 0.0;
1162   Double_t summedTracksPt = 0.0;
1163   
1164   // Setting invalid values
1165   rhoNoExclusion = 0.0;
1166   rhoConeExclusion02 = 0.0;
1167   rhoConeExclusion04 = 0.0;
1168   rhoConeExclusion06 = 0.0;
1169   rhoConeExclusion08 = 0.0; 
1170   rhoExactExclusion  = 0.0;
1171
1172   // Exclude UP TO numberExcludeLeadingJets
1173   if(numberExcludeLeadingJets==-1)
1174     numberExcludeLeadingJets = fNumberSignalJets;
1175   if (fNumberSignalJets < numberExcludeLeadingJets)
1176     numberExcludeLeadingJets = fNumberSignalJets;
1177
1178   Int_t fSignalJetCount5GeV = 0;
1179   for(Int_t j=0;j<numberExcludeLeadingJets;j++)
1180   {
1181     AliEmcalJet* signalJet = fSignalJets[j];
1182     if(signalJet->Pt() < 5.0)
1183       continue;
1184     fSignalJetCount5GeV++;
1185   }
1186
1187   for (Int_t i = 0; i < fTrackArray->GetEntries(); i++)
1188   {
1189     AliVTrack* tmpTrack = static_cast<AliVTrack*>(fTrackArray->At(i));
1190     Bool_t trackWithinJet = kFALSE; Bool_t trackWithin02Cone = kFALSE; Bool_t trackWithin04Cone = kFALSE; Bool_t trackWithin06Cone = kFALSE; Bool_t trackWithin08Cone = kFALSE;
1191
1192     if (IsTrackInAcceptance(tmpTrack))
1193     {
1194       // Check if tracks overlaps with jet
1195       for(Int_t j=0;j<numberExcludeLeadingJets;j++)
1196       {
1197         AliEmcalJet* signalJet = fSignalJets[j];
1198
1199         if(signalJet->Pt() < 5.0)
1200           continue;
1201
1202         // Exact jet exclusion
1203         if (IsTrackInJet(signalJet, i))
1204           trackWithinJet = kTRUE;
1205
1206         // Cone exclusions
1207         if (IsTrackInCone(tmpTrack, signalJet->Eta(), signalJet->Phi(), 0.2))
1208         {
1209           trackWithin02Cone = kTRUE;
1210           trackWithin04Cone = kTRUE;
1211           trackWithin06Cone = kTRUE;
1212           trackWithin08Cone = kTRUE;
1213           break;
1214         }
1215         else if (IsTrackInCone(tmpTrack, signalJet->Eta(), signalJet->Phi(), 0.4))
1216         {
1217           trackWithin04Cone = kTRUE;
1218           trackWithin06Cone = kTRUE;
1219           trackWithin08Cone = kTRUE;
1220         }
1221         else if (IsTrackInCone(tmpTrack, signalJet->Eta(), signalJet->Phi(), 0.6))
1222         {
1223           trackWithin06Cone = kTRUE;
1224           trackWithin08Cone = kTRUE;
1225         }
1226         else if (IsTrackInCone(tmpTrack, signalJet->Eta(), signalJet->Phi(), 0.8))
1227         {
1228           trackWithin08Cone = kTRUE;
1229         }
1230       }
1231
1232       if(!trackWithin08Cone)
1233       {
1234         summedTracksPtCone08 += tmpTrack->Pt();
1235       }
1236       if(!trackWithin06Cone)
1237       {
1238         summedTracksPtCone06 += tmpTrack->Pt();
1239       }
1240       if(!trackWithin04Cone)
1241       {
1242         summedTracksPtCone04 += tmpTrack->Pt();
1243       }
1244       if(!trackWithin02Cone)
1245       {
1246         summedTracksPtCone02 += tmpTrack->Pt();
1247       }
1248       if(!trackWithinJet)
1249       {
1250         summedTracksPtWithinJets += tmpTrack->Pt();
1251       }
1252       summedTracksPt += tmpTrack->Pt();
1253
1254     }
1255   }
1256
1257   // Calculate the correct area where the tracks were taking from
1258
1259   Double_t tmpFullTPCArea = (2.0*(fMaxEta-fMinEta)) * TMath::TwoPi();
1260   Double_t tmpAreaCone02     = tmpFullTPCArea;
1261   Double_t tmpAreaCone04     = tmpFullTPCArea;
1262   Double_t tmpAreaCone06     = tmpFullTPCArea;
1263   Double_t tmpAreaCone08     = tmpFullTPCArea;
1264   Double_t tmpAreaWithinJets = tmpFullTPCArea;
1265   std::vector<Double_t> tmpEtas(fSignalJetCount5GeV);
1266   std::vector<Double_t> tmpPhis(fSignalJetCount5GeV);
1267
1268   Int_t iSignal = 0;
1269   for(Int_t i=0;i<numberExcludeLeadingJets;i++)
1270   {
1271     AliEmcalJet* tmpJet = fSignalJets[i];
1272
1273     if(tmpJet->Pt() < 5.0)
1274       continue;
1275
1276     tmpEtas[iSignal] = tmpJet->Eta();
1277     tmpPhis[iSignal] = tmpJet->Phi();
1278     tmpAreaWithinJets -= tmpJet->Area();
1279
1280     iSignal++;
1281   }
1282
1283   tmpAreaCone02 -= tmpFullTPCArea * MCGetOverlapMultipleCirclesRectancle(fSignalJetCount5GeV, tmpEtas, tmpPhis, 0.2, fMinEta, fMaxEta, 0., TMath::TwoPi());
1284   tmpAreaCone04 -= tmpFullTPCArea * MCGetOverlapMultipleCirclesRectancle(fSignalJetCount5GeV, tmpEtas, tmpPhis, 0.4, fMinEta, fMaxEta, 0., TMath::TwoPi());
1285   tmpAreaCone06 -= tmpFullTPCArea * MCGetOverlapMultipleCirclesRectancle(fSignalJetCount5GeV, tmpEtas, tmpPhis, 0.6, fMinEta, fMaxEta, 0., TMath::TwoPi());
1286   tmpAreaCone08 -= tmpFullTPCArea * MCGetOverlapMultipleCirclesRectancle(fSignalJetCount5GeV, tmpEtas, tmpPhis, 0.8, fMinEta, fMaxEta, 0., TMath::TwoPi());
1287  
1288   rhoConeExclusion02 = summedTracksPtCone02/tmpAreaCone02;
1289   rhoConeExclusion04 = summedTracksPtCone04/tmpAreaCone04;
1290   rhoConeExclusion06 = summedTracksPtCone06/tmpAreaCone06;
1291   rhoConeExclusion08 = summedTracksPtCone08/tmpAreaCone08;
1292   rhoExactExclusion  = summedTracksPtWithinJets/tmpAreaWithinJets;
1293   rhoNoExclusion     = summedTracksPt/tmpFullTPCArea;
1294
1295
1296   #ifdef DEBUGMODE
1297     AliInfo("Got TR background density.");
1298   #endif
1299 }
1300
1301 //________________________________________________________________________
1302 void AliAnalysisTaskChargedJetsPA::GetTRBackgroundDensity(Int_t numberExcludeLeadingJets, Double_t& rhoMean, Double_t& area, AliEmcalJet* excludeJet1, AliEmcalJet* excludeJet2, Bool_t doSearchPerpendicular)
1303 {
1304   #ifdef DEBUGMODE
1305     AliInfo("Getting TR background density.");
1306   #endif
1307
1308   // Setting invalid values
1309   Double_t summedTracksPt = 0.0;
1310   rhoMean = 0.0;
1311   area = -1.0;
1312
1313   Double_t tmpRadius = 0.0;
1314   if (doSearchPerpendicular)
1315     tmpRadius = 0.4*TMath::Pi(); // exclude 90 degrees around jets
1316   else
1317     tmpRadius = 0.8;
1318     
1319   numberExcludeLeadingJets = 2; // dijet is excluded here in any case
1320
1321
1322
1323   if (!fTrackArray || !fJetArray)
1324   {
1325     AliError("Could not get the track/jet array to calculate track rho!");
1326     return;
1327   }
1328
1329   Int_t   trackCount = fTrackArray->GetEntries();
1330   Int_t   trackCountAccepted = 0;
1331   for (Int_t i = 0; i < trackCount; i++)
1332   {
1333     AliVTrack* tmpTrack = static_cast<AliVTrack*>(fTrackArray->At(i));
1334     if (IsTrackInAcceptance(tmpTrack))
1335     {
1336       if (IsTrackInCone(tmpTrack, excludeJet1->Eta(), excludeJet1->Phi(), tmpRadius))
1337         continue;
1338
1339       if (numberExcludeLeadingJets > 1)
1340         if (IsTrackInCone(tmpTrack, excludeJet2->Eta(), excludeJet2->Phi(), tmpRadius))
1341           continue;
1342
1343         // Add track pt to array
1344         summedTracksPt = summedTracksPt + tmpTrack->Pt();
1345         trackCountAccepted++;
1346     }
1347   }
1348
1349   if (trackCountAccepted > 0)
1350   {
1351     Double_t tmpArea = 2.0*(fMaxEta-fMinEta)*TMath::TwoPi()  - 2*(tmpRadius*tmpRadius * TMath::Pi()); //TPC area - excluding jet area
1352     rhoMean = summedTracksPt/tmpArea;
1353     area = tmpArea;
1354   }
1355
1356   #ifdef DEBUGMODE
1357     AliInfo("Got TR background density.");
1358   #endif
1359 }
1360
1361 //________________________________________________________________________
1362 void AliAnalysisTaskChargedJetsPA::GetPPBackgroundDensity(Double_t& background)
1363 {
1364   // This is the background that was used for the pp 7 TeV ALICE paper
1365   // The background is estimated using the leading jet
1366
1367   background = 0;
1368
1369   // Get leading jet
1370   Int_t leadingJets[]   = {-1, -1};
1371   GetLeadingJets(fJetArray, &leadingJets[0], kTRUE);
1372   AliEmcalJet* jet = NULL;
1373   if(leadingJets[0] != -1)
1374     jet = static_cast<AliEmcalJet*>(fJetArray->At(leadingJets[0]));
1375   else
1376     return;
1377
1378   Double_t jetMom[3] = { jet->Px(), jet->Py(), jet->Pz() };
1379   TVector3 jet3mom1(jetMom);
1380   TVector3 jet3mom2(jetMom);
1381
1382   jet3mom1.RotateZ(TMath::Pi());
1383   jet3mom2.RotateZ(-TMath::Pi());
1384
1385   for (int i = 0; i < fTrackArray->GetEntries(); i++)
1386   {
1387     AliVTrack* track = static_cast<AliVTrack*>(fTrackArray->At(i));
1388     if (!IsTrackInAcceptance(track))
1389       continue;
1390
1391     Double_t trackMom[3] = { track->Px(), track->Py(), track->Pz() };
1392     TVector3 track3mom(trackMom);
1393
1394     Double_t dR1 = jet3mom1.DeltaR(track3mom);
1395     Double_t dR2 = jet3mom2.DeltaR(track3mom);
1396
1397     if (dR1 <= fSignalJetRadius || dR2 <= fSignalJetRadius)
1398       background += track3mom.Pt();
1399   }
1400
1401   background /= (2 * TMath::Pi() * fSignalJetRadius * fSignalJetRadius);
1402 }
1403
1404 //________________________________________________________________________
1405 void AliAnalysisTaskChargedJetsPA::Calculate(AliVEvent* event)
1406 {
1407   #ifdef DEBUGMODE
1408     AliInfo("Starting Calculate().");
1409   #endif
1410   ////////////////////// NOTE: initialization & casting
1411
1412   fEventCounter++;
1413
1414   // Check, if analysis should be done in pt hard bins
1415   if(fUsePtHardBin != -1)
1416     if(GetPtHardBin() != fUsePtHardBin)
1417       return;
1418
1419   // This is to take only every Nth event
1420   if((fEventCounter+fPartialAnalysisIndex) % fPartialAnalysisNParts != 0)
1421     return;
1422
1423   FillHistogram("hNumberEvents",0.5);
1424
1425   if(!IsEventInAcceptance(event))
1426     return;
1427
1428   FillHistogram("hNumberEvents",1.5);
1429
1430   #ifdef DEBUGMODE
1431     AliInfo("Calculate()::Init done.");
1432   #endif
1433
1434   ////////////////////// NOTE: Get Centrality, (Leading)Signal jets and Background
1435
1436   // Get centrality
1437   AliCentrality* tmpCentrality = NULL;
1438   tmpCentrality = event->GetCentrality();
1439   Double_t centralityPercentile = -1.0;
1440   Double_t centralityPercentileV0A = 0.0;
1441   Double_t centralityPercentileV0C = 0.0;
1442   Double_t centralityPercentileV0M = 0.0;
1443   Double_t centralityPercentileZNA = 0.0;
1444   if (tmpCentrality != NULL)
1445   {
1446     centralityPercentile = tmpCentrality->GetCentralityPercentile(fCentralityType.Data());
1447     centralityPercentileV0A = tmpCentrality->GetCentralityPercentile("V0A");
1448     centralityPercentileV0C = tmpCentrality->GetCentralityPercentile("V0C");
1449     centralityPercentileV0M = tmpCentrality->GetCentralityPercentile("V0M");
1450     centralityPercentileZNA = tmpCentrality->GetCentralityPercentile("ZNA");
1451   }
1452
1453   if((centralityPercentile < 0.0) || (centralityPercentile > 100.0))
1454     AliWarning(Form("Centrality value not valid (c=%E)",centralityPercentile)); 
1455
1456   if(fSetCentralityToOne)
1457     centralityPercentile = 1.0;
1458
1459   // Get jets
1460   if (fAnalyzeBackground || fAnalyzeJets)
1461     GetSignalJets();
1462
1463   // Get background estimates
1464   Double_t              backgroundKTImprovedCMS = -1.0;
1465   Double_t              backgroundKTImprovedCMSExternal = -1.0;
1466   Double_t              backgroundKTPbPb = -1.0;
1467   Double_t              backgroundKTPbPbWithGhosts = -1.0;
1468   Double_t              backgroundKTCMS = -1.0;
1469   Double_t              backgroundKTMean = -1.0;
1470   Double_t              backgroundKTTrackLike = -1.0;
1471   Double_t              backgroundTRNoExcl = -1.0;
1472   Double_t              backgroundTRCone02 = -1.0;
1473   Double_t              backgroundTRCone04 = -1.0;
1474   Double_t              backgroundTRCone06 = -1.0;
1475   Double_t              backgroundTRCone08 = -1.0;
1476   Double_t              backgroundTRExact  = -1.0;
1477   Double_t              backgroundPP       = -1.0;
1478   Double_t              backgroundJetProfile = -1.0;
1479   // Calculate background for different jet exclusions
1480
1481   if (fAnalyzeBackground)
1482   {
1483
1484     if(fAnalyzeDeprecatedBackgrounds)
1485       GetKTBackgroundDensityAll    (fNumberExcludedJets, backgroundKTPbPb, backgroundKTPbPbWithGhosts, backgroundKTCMS, backgroundKTImprovedCMS, backgroundKTMean, backgroundKTTrackLike);
1486     else
1487       GetKTBackgroundDensity       (fNumberExcludedJets, backgroundKTImprovedCMS);
1488
1489     if(fAnalyzeDeprecatedBackgrounds)
1490       GetTRBackgroundDensity    (fNumberExcludedJets, backgroundTRNoExcl, backgroundTRCone02, backgroundTRCone04, backgroundTRCone06, backgroundTRCone08, backgroundTRExact);
1491
1492     // pp background
1493     GetPPBackgroundDensity(backgroundPP);
1494
1495     backgroundKTImprovedCMSExternal = GetExternalRho();
1496     if(fNoExternalBackground)
1497       backgroundKTImprovedCMSExternal = 0;
1498
1499     if(fBackgroundForJetProfile==0)
1500       backgroundJetProfile = backgroundKTImprovedCMSExternal;
1501     else if(fBackgroundForJetProfile==1)
1502       backgroundJetProfile = backgroundKTImprovedCMS;
1503     else if(fBackgroundForJetProfile==2)
1504       backgroundJetProfile = backgroundKTCMS;
1505     else if(fBackgroundForJetProfile==3)
1506       backgroundJetProfile = backgroundPP;
1507     else if(fBackgroundForJetProfile==4)
1508       backgroundJetProfile = backgroundTRCone06;
1509     else if(fBackgroundForJetProfile==5)
1510       backgroundJetProfile = 0;
1511   }
1512
1513   #ifdef DEBUGMODE
1514     AliInfo("Calculate()::Centrality&SignalJets&Background-Calculation done.");
1515   #endif
1516
1517   if (fAnalyzeQA)
1518   {
1519     FillHistogram("hVertexX",event->GetPrimaryVertex()->GetX());
1520     FillHistogram("hVertexY",event->GetPrimaryVertex()->GetY());
1521     FillHistogram("hVertexXY",event->GetPrimaryVertex()->GetX(), event->GetPrimaryVertex()->GetY());
1522     FillHistogram("hVertexR",TMath::Sqrt(event->GetPrimaryVertex()->GetX()*event->GetPrimaryVertex()->GetX() + event->GetPrimaryVertex()->GetY()*event->GetPrimaryVertex()->GetY()));
1523     FillHistogram("hCentralityV0M",centralityPercentileV0M);
1524     FillHistogram("hCentralityV0A",centralityPercentileV0A);
1525     FillHistogram("hCentralityV0C",centralityPercentileV0C);
1526     FillHistogram("hCentralityZNA",centralityPercentileZNA);
1527     FillHistogram("hCentrality",centralityPercentile);
1528
1529     Int_t trackCountAcc = 0;
1530     Int_t nTracks = fTrackArray->GetEntries();
1531     for (Int_t i = 0; i < nTracks; i++)
1532     {
1533       AliVTrack* track = static_cast<AliVTrack*>(fTrackArray->At(i));
1534
1535       if (track != 0)
1536         if (track->Pt() >= fMinTrackPt)
1537         {
1538           FillHistogram("hTrackPhiEta", track->Phi(),track->Eta(), 1);
1539           FillHistogram("hTrackPtPhiEta", track->Phi(),track->Eta(), track->Pt());
1540         }
1541
1542       if (IsTrackInAcceptance(track))
1543       {
1544         FillHistogram("hTrackPt", track->Pt(), centralityPercentile);
1545         if(track->Eta() >= 0)
1546           FillHistogram("hTrackPtPosEta", track->Pt(), centralityPercentile);
1547         else
1548           FillHistogram("hTrackPtNegEta", track->Pt(), centralityPercentile);
1549                 
1550         FillHistogram("hTrackEta", track->Eta(), centralityPercentile);
1551         FillHistogram("hTrackPhi", track->Phi());
1552         
1553         if(static_cast<AliPicoTrack*>(track))
1554           FillHistogram("hTrackPhiTrackType", track->Phi(), (static_cast<AliPicoTrack*>(track))->GetTrackType());
1555
1556         for(Int_t j=0;j<20;j++)
1557           if(track->Pt() > j)
1558             FillHistogram("hTrackPhiPtCut", track->Phi(), track->Pt());
1559
1560         FillHistogram("hTrackCharge", track->Charge());
1561         trackCountAcc++;
1562       }
1563     }
1564     FillHistogram("hTrackCountAcc", trackCountAcc, centralityPercentile);
1565
1566
1567 /*
1568     // This is code that is run locally to get some special events
1569     TFile* fileOutput = new TFile("SpecialEvents.root", "UPDATE");
1570     if(fSecondLeadingJet&&(fSecondLeadingJet->Pt()>10.))
1571     {
1572       cout << "Event found\n";
1573       TH2D* tmpEvent = new TH2D(Form("Event%lu", fEventCounter), "", 40, -0.9, 0.9, 40, 0., TMath::TwoPi());
1574       tmpEvent->GetXaxis()->SetTitle("#eta");
1575       tmpEvent->GetYaxis()->SetTitle("#phi");
1576       tmpEvent->SetOption("LEGO2");
1577       tmpEvent->Sumw2();
1578
1579       for (Int_t i = 0; i < nTracks; i++)
1580       {
1581         AliVTrack* track = static_cast<AliVTrack*>(fTrackArray->At(i));
1582
1583         if (IsTrackInAcceptance(track))
1584         {
1585           tmpEvent->Fill(track->Eta(), track->Phi(), track->Pt());
1586         }
1587       }
1588       tmpEvent->Write(0, kOverwrite);
1589     }
1590     fileOutput->Close();
1591 */
1592   }
1593   #ifdef DEBUGMODE
1594     AliInfo("Calculate()::QA done.");
1595   #endif
1596
1597   ////////////////////// NOTE: Jet analysis and calculations
1598
1599   if (fAnalyzeJets)
1600   {
1601     for (Int_t i = 0; i<fJetArray->GetEntries(); i++)
1602     {
1603       AliEmcalJet* tmpJet = static_cast<AliEmcalJet*>(fJetArray->At(i));
1604       if (!tmpJet)
1605         continue;
1606
1607       // ### RAW JET ANALYSIS
1608       if (tmpJet->Area() >= fMinJetArea)
1609         FillHistogram("hRawJetPhiEta", tmpJet->Phi(), tmpJet->Eta());
1610       if ((tmpJet->Eta() >= fMinJetEta) && (tmpJet->Eta() < fMaxJetEta))
1611         FillHistogram("hRawJetArea", tmpJet->Area());
1612
1613       // Jet pt for different area cut
1614       FillHistogram("hJetPtCutStages", tmpJet->Pt(), 0.5);
1615       if ((tmpJet->Eta() >= fMinJetEta) && (tmpJet->Eta() < fMaxJetEta))
1616       {
1617         FillHistogram("hJetPtCutStages", tmpJet->Pt(), 1.5);
1618         if (tmpJet->Pt() >= fMinJetPt)
1619         {
1620           FillHistogram("hJetPtCutStages", tmpJet->Pt(), 2.5);
1621           if (tmpJet->Area() >= fMinJetArea)
1622           {
1623             FillHistogram("hJetPtCutStages", tmpJet->Pt(), 3.5);
1624           }
1625         }
1626       }
1627
1628
1629       if(IsSignalJetInAcceptance(tmpJet))
1630       {
1631       // ### SIGNAL JET ANALYSIS
1632         // Jet spectra
1633         FillHistogram("hJetPtBgrdSubtractedKTImprovedCMS", GetCorrectedJetPt(tmpJet, backgroundKTImprovedCMS), centralityPercentile);
1634         FillHistogram("hJetPtBgrdSubtractedPP", GetCorrectedJetPt(tmpJet, backgroundPP), centralityPercentile);
1635         FillHistogram("hJetPtBgrdSubtractedExternal", GetCorrectedJetPt(tmpJet, backgroundKTImprovedCMSExternal), centralityPercentile);
1636         if(tmpJet->Phi() >= TMath::Pi())
1637           FillHistogram("hJetPtBgrdSubtractedExternal_Phi2", GetCorrectedJetPt(tmpJet, backgroundKTImprovedCMSExternal), centralityPercentile);
1638         else          
1639           FillHistogram("hJetPtBgrdSubtractedExternal_Phi1", GetCorrectedJetPt(tmpJet, backgroundKTImprovedCMSExternal), centralityPercentile);
1640
1641         FillHistogram("hJetPtSubtractedRhoKTImprovedCMS", tmpJet->Pt(), centralityPercentile, tmpJet->Pt() - GetCorrectedJetPt(tmpJet, backgroundKTImprovedCMS));
1642         FillHistogram("hJetPtSubtractedRhoExternal", tmpJet->Pt(), centralityPercentile, tmpJet->Pt() - GetCorrectedJetPt(tmpJet, backgroundKTImprovedCMSExternal));
1643         FillHistogram("hJetPtSubtractedRhoPP", tmpJet->Pt(), centralityPercentile, tmpJet->Pt() - GetCorrectedJetPt(tmpJet, backgroundPP));
1644         FillHistogram("hDeltaPtExternalBgrdVsPt", GetDeltaPt(backgroundKTImprovedCMSExternal), GetCorrectedJetPt(tmpJet, backgroundKTImprovedCMSExternal));
1645
1646
1647         if(fAnalyzeDeprecatedBackgrounds)
1648         {
1649           FillHistogram("hJetPtBgrdSubtractedTR", GetCorrectedJetPt(tmpJet, backgroundTRCone06), centralityPercentile);
1650           FillHistogram("hJetPtBgrdSubtractedKTPbPb", GetCorrectedJetPt(tmpJet, backgroundKTPbPb), centralityPercentile);
1651           FillHistogram("hJetPtBgrdSubtractedKTPbPbWithGhosts", GetCorrectedJetPt(tmpJet, backgroundKTPbPbWithGhosts), centralityPercentile);
1652           FillHistogram("hJetPtBgrdSubtractedKTCMS", GetCorrectedJetPt(tmpJet, backgroundKTCMS), centralityPercentile);
1653           FillHistogram("hJetPtBgrdSubtractedKTMean", GetCorrectedJetPt(tmpJet, backgroundKTMean), centralityPercentile);
1654           FillHistogram("hJetPtBgrdSubtractedKTTrackLike", GetCorrectedJetPt(tmpJet, backgroundKTTrackLike), centralityPercentile);
1655         }
1656
1657         // Jet profile analysis
1658         if(TMath::Abs(tmpJet->Eta()) <= 0.3)
1659         {
1660           if(tmpJet->Pt()>=70.0)
1661           {
1662             FillHistogram("hJetProfile70GeV", 0.05-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.05, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1663             FillHistogram("hJetProfile70GeV", 0.10-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.10, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1664             FillHistogram("hJetProfile70GeV", 0.15-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.15, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1665             FillHistogram("hJetProfile70GeV", 0.20-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.20, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1666             FillHistogram("hJetProfile70GeV", 0.25-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.25, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1667             FillHistogram("hJetProfile70GeV", 0.30-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.30, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1668             FillHistogram("hJetProfile70GeV", 0.35-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.35, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1669             FillHistogram("hJetProfile70GeV", 0.40-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.40, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1670             FillHistogram("hJetProfile70GeV", 0.45-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.45, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1671             FillHistogram("hJetProfile70GeV", 0.50-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.50, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1672             FillHistogram("hJetProfile70GeV", 0.55-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.55, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1673             FillHistogram("hJetProfile70GeV", 0.60-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.60, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1674           }
1675           else if(GetCorrectedJetPt(tmpJet, backgroundJetProfile)>=60.0)
1676           {
1677             FillHistogram("hJetProfile60GeV", 0.05-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.05, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1678             FillHistogram("hJetProfile60GeV", 0.10-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.10, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1679             FillHistogram("hJetProfile60GeV", 0.15-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.15, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1680             FillHistogram("hJetProfile60GeV", 0.20-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.20, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1681             FillHistogram("hJetProfile60GeV", 0.25-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.25, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1682             FillHistogram("hJetProfile60GeV", 0.30-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.30, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1683             FillHistogram("hJetProfile60GeV", 0.35-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.35, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1684             FillHistogram("hJetProfile60GeV", 0.40-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.40, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1685             FillHistogram("hJetProfile60GeV", 0.45-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.45, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1686             FillHistogram("hJetProfile60GeV", 0.50-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.50, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1687             FillHistogram("hJetProfile60GeV", 0.55-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.55, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1688             FillHistogram("hJetProfile60GeV", 0.60-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.60, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1689           }
1690           else if(GetCorrectedJetPt(tmpJet, backgroundJetProfile)>=50.0)
1691           {
1692             FillHistogram("hJetProfile50GeV", 0.05-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.05, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1693             FillHistogram("hJetProfile50GeV", 0.10-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.10, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1694             FillHistogram("hJetProfile50GeV", 0.15-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.15, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1695             FillHistogram("hJetProfile50GeV", 0.20-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.20, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1696             FillHistogram("hJetProfile50GeV", 0.25-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.25, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1697             FillHistogram("hJetProfile50GeV", 0.30-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.30, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1698             FillHistogram("hJetProfile50GeV", 0.35-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.35, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1699             FillHistogram("hJetProfile50GeV", 0.40-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.40, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1700             FillHistogram("hJetProfile50GeV", 0.45-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.45, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1701             FillHistogram("hJetProfile50GeV", 0.50-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.50, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1702             FillHistogram("hJetProfile50GeV", 0.55-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.55, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1703             FillHistogram("hJetProfile50GeV", 0.60-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.60, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1704           }
1705           else if(GetCorrectedJetPt(tmpJet, backgroundJetProfile)>=40.0)
1706           {
1707             FillHistogram("hJetProfile40GeV", 0.05-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.05, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1708             FillHistogram("hJetProfile40GeV", 0.10-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.10, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1709             FillHistogram("hJetProfile40GeV", 0.15-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.15, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1710             FillHistogram("hJetProfile40GeV", 0.20-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.20, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1711             FillHistogram("hJetProfile40GeV", 0.25-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.25, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1712             FillHistogram("hJetProfile40GeV", 0.30-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.30, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1713             FillHistogram("hJetProfile40GeV", 0.35-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.35, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1714             FillHistogram("hJetProfile40GeV", 0.40-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.40, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1715             FillHistogram("hJetProfile40GeV", 0.45-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.45, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1716             FillHistogram("hJetProfile40GeV", 0.50-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.50, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1717             FillHistogram("hJetProfile40GeV", 0.55-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.55, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1718             FillHistogram("hJetProfile40GeV", 0.60-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.60, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1719           }
1720           else if(GetCorrectedJetPt(tmpJet, backgroundJetProfile)>=30.0)
1721           {
1722             FillHistogram("hJetProfile30GeV", 0.05-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.05, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1723             FillHistogram("hJetProfile30GeV", 0.10-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.10, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1724             FillHistogram("hJetProfile30GeV", 0.15-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.15, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1725             FillHistogram("hJetProfile30GeV", 0.20-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.20, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1726             FillHistogram("hJetProfile30GeV", 0.25-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.25, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1727             FillHistogram("hJetProfile30GeV", 0.30-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.30, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1728             FillHistogram("hJetProfile30GeV", 0.35-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.35, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1729             FillHistogram("hJetProfile30GeV", 0.40-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.40, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1730             FillHistogram("hJetProfile30GeV", 0.45-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.45, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1731             FillHistogram("hJetProfile30GeV", 0.50-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.50, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1732             FillHistogram("hJetProfile30GeV", 0.55-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.55, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1733             FillHistogram("hJetProfile30GeV", 0.60-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.60, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1734           }
1735           else if(GetCorrectedJetPt(tmpJet, backgroundJetProfile)>=20.0)
1736           {
1737             FillHistogram("hJetProfile20GeV", 0.05-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.05, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1738             FillHistogram("hJetProfile20GeV", 0.10-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.10, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1739             FillHistogram("hJetProfile20GeV", 0.15-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.15, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1740             FillHistogram("hJetProfile20GeV", 0.20-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.20, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1741             FillHistogram("hJetProfile20GeV", 0.25-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.25, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1742             FillHistogram("hJetProfile20GeV", 0.30-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.30, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1743             FillHistogram("hJetProfile20GeV", 0.35-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.35, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1744             FillHistogram("hJetProfile20GeV", 0.40-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.40, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1745             FillHistogram("hJetProfile20GeV", 0.45-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.45, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1746             FillHistogram("hJetProfile20GeV", 0.50-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.50, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1747             FillHistogram("hJetProfile20GeV", 0.55-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.55, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1748             FillHistogram("hJetProfile20GeV", 0.60-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.60, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1749           }
1750           else if(GetCorrectedJetPt(tmpJet, backgroundJetProfile)>=10.0)
1751           {
1752             FillHistogram("hJetProfile10GeV", 0.05-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.05, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1753             FillHistogram("hJetProfile10GeV", 0.10-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.10, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1754             FillHistogram("hJetProfile10GeV", 0.15-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.15, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1755             FillHistogram("hJetProfile10GeV", 0.20-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.20, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1756             FillHistogram("hJetProfile10GeV", 0.25-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.25, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1757             FillHistogram("hJetProfile10GeV", 0.30-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.30, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1758             FillHistogram("hJetProfile10GeV", 0.35-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.35, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1759             FillHistogram("hJetProfile10GeV", 0.40-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.40, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1760             FillHistogram("hJetProfile10GeV", 0.45-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.45, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1761             FillHistogram("hJetProfile10GeV", 0.50-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.50, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1762             FillHistogram("hJetProfile10GeV", 0.55-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.55, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1763             FillHistogram("hJetProfile10GeV", 0.60-0.05/2, (GetCorrectedConePt(tmpJet->Eta(), tmpJet->Phi(), 0.60, backgroundJetProfile))/GetCorrectedJetPt(tmpJet, backgroundJetProfile));
1764           }
1765         }
1766         FillHistogram("hJetPtVsConstituentCount", tmpJet->Pt(),tmpJet->GetNumberOfTracks());
1767
1768         if((fAnalyzeQA) && (tmpJet->Pt() >= 5.0))
1769         {
1770           Double_t lowestTrackPt = 1e99;
1771           Double_t highestTrackPt = 0.0;
1772           for(Int_t j=0; j<tmpJet->GetNumberOfTracks(); j++)
1773           {
1774             FillHistogram("hJetConstituentPt", tmpJet->TrackAt(j, fTrackArray)->Pt(), centralityPercentile);
1775             // Find the lowest pT of a track in the jet
1776             if (tmpJet->TrackAt(j, fTrackArray)->Pt() < lowestTrackPt)
1777               lowestTrackPt = tmpJet->TrackAt(j, fTrackArray)->Pt();
1778             if (tmpJet->TrackAt(j, fTrackArray)->Pt() > highestTrackPt)
1779               highestTrackPt = tmpJet->TrackAt(j, fTrackArray)->Pt();
1780           }
1781           FillHistogram("hJetArea", tmpJet->Area(), tmpJet->Pt());
1782           // Signal jet vs. signal jet - "Combinatorial"
1783           for (Int_t j = i+1; j<fNumberSignalJets; j++)
1784             FillHistogram("hJetDeltaPhi", GetDeltaPhi(tmpJet->Phi(), fSignalJets[j]->Phi()));
1785
1786           FillHistogram("hJetPhiEta", tmpJet->Phi(),tmpJet->Eta());
1787           FillHistogram("hJetPtPhiEta", tmpJet->Phi(),tmpJet->Eta(),tmpJet->Pt());
1788           FillHistogram("hJetEta", tmpJet->Eta(), centralityPercentile);
1789
1790           if(lowestTrackPt>=2.0)
1791             FillHistogram("hJetEta2GeVTracks", tmpJet->Eta(), centralityPercentile);
1792           if(lowestTrackPt>=4.0)
1793             FillHistogram("hJetEta4GeVTracks", tmpJet->Eta(), centralityPercentile);
1794         }
1795
1796
1797         // Jet constituent mass analyses (for PYTHIA)
1798         if (fAnalyzeMassCorrelation)
1799         {
1800           Double_t jetMass = 0.0;
1801           Int_t constituentCount = tmpJet->GetNumberOfTracks();
1802           Int_t electronCount = 0;
1803           Int_t pionCount = 0;
1804           Int_t protonCount = 0;
1805           Int_t kaonCount = 0;
1806           Int_t othersCount = 0;
1807
1808           for(Int_t j=0; j<tmpJet->GetNumberOfTracks(); j++)
1809           {
1810             AliVParticle* tmpParticle = tmpJet->TrackAt(j, fTrackArray);
1811
1812             jetMass += tmpParticle->M();
1813             if(TMath::Abs(tmpParticle->PdgCode()) == 11) // electron, positron
1814               electronCount++;
1815             else if(TMath::Abs(tmpParticle->PdgCode()) == 211) // p-, p+
1816               pionCount++;
1817             else if(TMath::Abs(tmpParticle->PdgCode()) == 2212) // p, pbar
1818               protonCount++;
1819             else if(TMath::Abs(tmpParticle->PdgCode()) == 321) // kaon+,kaon-
1820               kaonCount++;
1821             else
1822               othersCount++;
1823           }
1824
1825           FillHistogram("hJetMassFromConstituents", jetMass);
1826           FillHistogram("hJetMass", tmpJet->M());
1827
1828           FillHistogram("hJetPtVsMass", tmpJet->Pt(), jetMass);
1829           FillHistogram("hJetPtVsJetMass", tmpJet->Pt(), tmpJet->M());
1830
1831           FillHistogram("hJetPtVsProtonCount", tmpJet->Pt(), protonCount/(static_cast<Double_t>(constituentCount)));
1832           FillHistogram("hJetPtVsPionCount", tmpJet->Pt(), pionCount/(static_cast<Double_t>(constituentCount)));
1833           FillHistogram("hJetPtVsKaonCount", tmpJet->Pt(), kaonCount/(static_cast<Double_t>(constituentCount)));
1834           FillHistogram("hJetPtVsElectronCount", tmpJet->Pt(), electronCount/(static_cast<Double_t>(constituentCount)));
1835           FillHistogram("hJetPtVsOthersCount", tmpJet->Pt(), othersCount/(static_cast<Double_t>(constituentCount)));
1836           FillHistogram("hJetConstituentProtonCount", protonCount);
1837           FillHistogram("hJetConstituentPionCount", pionCount);
1838           FillHistogram("hJetConstituentKaonCount", kaonCount);
1839           FillHistogram("hJetConstituentElectronCount", electronCount);
1840           FillHistogram("hJetConstituentOthersCount", othersCount);
1841
1842           // Results for the "normal" jet (avoiding single particle jets)
1843           if((constituentCount>=6) && (constituentCount<=14))
1844           {
1845             FillHistogram("hJetPtVsMass_6_14", tmpJet->Pt(), jetMass);
1846             FillHistogram("hJetPtVsJetMass_6_14", tmpJet->Pt(), tmpJet->M());
1847
1848             FillHistogram("hJetPtVsProtonCount_6_14", tmpJet->Pt(), protonCount);
1849             FillHistogram("hJetPtVsPionCount_6_14", tmpJet->Pt(), pionCount);
1850             FillHistogram("hJetPtVsKaonCount_6_14", tmpJet->Pt(), kaonCount);
1851             FillHistogram("hJetPtVsElectronCount_6_14", tmpJet->Pt(), electronCount);
1852             FillHistogram("hJetPtVsOthersCount_6_14", tmpJet->Pt(), othersCount);
1853             FillHistogram("hJetConstituentProtonCount_6_14", protonCount);
1854             FillHistogram("hJetConstituentPionCount_6_14", pionCount);
1855             FillHistogram("hJetConstituentKaonCount_6_14", kaonCount);
1856             FillHistogram("hJetConstituentElectronCount_6_14", electronCount);
1857             FillHistogram("hJetConstituentOthersCount_6_14", othersCount);
1858           }
1859           if((constituentCount>=2))
1860           {
1861             FillHistogram("hJetPtVsMass_2_X", tmpJet->Pt(), jetMass);
1862             FillHistogram("hJetPtVsJetMass_2_X", tmpJet->Pt(), tmpJet->M());
1863
1864             FillHistogram("hJetPtVsProtonCount_2_X", tmpJet->Pt(), protonCount);
1865             FillHistogram("hJetPtVsPionCount_2_X", tmpJet->Pt(), pionCount);
1866             FillHistogram("hJetPtVsKaonCount_2_X", tmpJet->Pt(), kaonCount);
1867             FillHistogram("hJetPtVsElectronCount_2_X", tmpJet->Pt(), electronCount);
1868             FillHistogram("hJetPtVsOthersCount_2_X", tmpJet->Pt(), othersCount);
1869             FillHistogram("hJetConstituentProtonCount_2_X", protonCount);
1870             FillHistogram("hJetConstituentPionCount_2_X", pionCount);
1871             FillHistogram("hJetConstituentKaonCount_2_X", kaonCount);
1872             FillHistogram("hJetConstituentElectronCount_2_X", electronCount);
1873             FillHistogram("hJetConstituentOthersCount_2_X", othersCount);
1874           }
1875           if((constituentCount>=2) && (constituentCount<=7))
1876           {
1877             FillHistogram("hJetPtVsMass_2_7", tmpJet->Pt(), jetMass);
1878             FillHistogram("hJetPtVsJetMass_2_7", tmpJet->Pt(), tmpJet->M());
1879
1880             FillHistogram("hJetPtVsProtonCount_2_7", tmpJet->Pt(), protonCount);
1881             FillHistogram("hJetPtVsPionCount_2_7", tmpJet->Pt(), pionCount);
1882             FillHistogram("hJetPtVsKaonCount_2_7", tmpJet->Pt(), kaonCount);
1883             FillHistogram("hJetPtVsElectronCount_2_7", tmpJet->Pt(), electronCount);
1884             FillHistogram("hJetPtVsOthersCount_2_7", tmpJet->Pt(), othersCount);
1885             FillHistogram("hJetConstituentProtonCount_2_7", protonCount);
1886             FillHistogram("hJetConstituentPionCount_2_7", pionCount);
1887             FillHistogram("hJetConstituentKaonCount_2_7", kaonCount);
1888             FillHistogram("hJetConstituentElectronCount_2_7", electronCount);
1889             FillHistogram("hJetConstituentOthersCount_2_7", othersCount);
1890           }
1891         }
1892       }
1893     }
1894
1895     // ### SOME JET PLOTS
1896     FillHistogram("hJetCountAll", fJetArray->GetEntries());
1897     FillHistogram("hJetCountAccepted", fNumberSignalJets);
1898     FillHistogram("hJetCount", fJetArray->GetEntries(), fNumberSignalJets);
1899     if (fFirstLeadingJet)
1900     {
1901       FillHistogram("hLeadingJetPt", fFirstLeadingJet->Pt());
1902       FillHistogram("hCorrectedLeadingJetPt", GetCorrectedJetPt(fFirstLeadingJet,backgroundKTImprovedCMSExternal));
1903     }
1904     if (fSecondLeadingJet)
1905     {
1906       FillHistogram("hSecondLeadingJetPt", fSecondLeadingJet->Pt());
1907       FillHistogram("hCorrectedSecondLeadingJetPt", GetCorrectedJetPt(fSecondLeadingJet,backgroundKTImprovedCMSExternal));
1908     }
1909   } //endif AnalyzeJets
1910
1911   #ifdef DEBUGMODE
1912     AliInfo("Calculate()::Jets done.");
1913   #endif
1914   ////////////////////// NOTE: Background analysis
1915
1916   if (fAnalyzeBackground)
1917   {
1918     // Calculate background in centrality classes
1919     FillHistogram("hKTBackgroundImprovedCMS", backgroundKTImprovedCMS, centralityPercentile);
1920     FillHistogram("hKTBackgroundImprovedCMSExternal", backgroundKTImprovedCMSExternal, centralityPercentile);
1921     if(fFirstLeadingJet && (fFirstLeadingJet->Pt()>=20.))
1922       FillHistogram("hKTBackgroundImprovedCMSExternal20GeV", backgroundKTImprovedCMSExternal, centralityPercentile);
1923
1924
1925     FillHistogram("hPPBackground", backgroundPP, centralityPercentile);
1926     FillHistogram("hKTMeanBackgroundImprovedCMS", centralityPercentile, backgroundKTImprovedCMS);
1927
1928     if(fAnalyzeDeprecatedBackgrounds)
1929     {
1930       FillHistogram("hKTBackgroundPbPb", backgroundKTPbPb, centralityPercentile);
1931       FillHistogram("hKTBackgroundPbPbWithGhosts", backgroundKTPbPbWithGhosts, centralityPercentile);
1932       FillHistogram("hKTBackgroundCMS", backgroundKTCMS, centralityPercentile);
1933       FillHistogram("hKTBackgroundMean", backgroundKTMean, centralityPercentile);
1934       FillHistogram("hKTBackgroundTrackLike", backgroundKTTrackLike, centralityPercentile);
1935
1936       FillHistogram("hTRBackgroundNoExcl", backgroundTRNoExcl, centralityPercentile);
1937       FillHistogram("hTRBackgroundCone02", backgroundTRCone02, centralityPercentile);
1938       FillHistogram("hTRBackgroundCone04", backgroundTRCone04, centralityPercentile);
1939       FillHistogram("hTRBackgroundCone06", backgroundTRCone06, centralityPercentile);
1940       FillHistogram("hTRBackgroundCone08", backgroundTRCone08, centralityPercentile);
1941       FillHistogram("hTRBackgroundExact", backgroundTRExact, centralityPercentile);
1942
1943       // Calculate background profiles in terms of centrality
1944       FillHistogram("hKTMeanBackgroundPbPb", centralityPercentile,  backgroundKTPbPb);
1945       FillHistogram("hKTMeanBackgroundPbPbWithGhosts", centralityPercentile,  backgroundKTPbPbWithGhosts);
1946       FillHistogram("hKTMeanBackgroundCMS", centralityPercentile, backgroundKTCMS);
1947       FillHistogram("hKTMeanBackgroundMean", centralityPercentile, backgroundKTMean);
1948       FillHistogram("hKTMeanBackgroundTPC", centralityPercentile, backgroundKTTrackLike);
1949       FillHistogram("hTRMeanBackground", centralityPercentile,  backgroundTRCone06);
1950     }
1951
1952
1953     // Calculate the delta pt
1954     Double_t tmpDeltaPtNoBackground = GetDeltaPt(0.0);
1955     Double_t tmpDeltaPtKTImprovedCMS = GetDeltaPt(backgroundKTImprovedCMS);
1956     Double_t tmpDeltaPtExternalBgrd = GetDeltaPt(backgroundKTImprovedCMSExternal);
1957
1958
1959     Double_t tmpDeltaPtKTImprovedCMSPartialExclusion = 0.0;
1960     if(fNcoll)
1961       tmpDeltaPtKTImprovedCMSPartialExclusion = GetDeltaPt(backgroundKTImprovedCMS, 1.0/fNcoll);
1962     else
1963       tmpDeltaPtKTImprovedCMSPartialExclusion = GetDeltaPt(backgroundKTImprovedCMS, 1.0);
1964
1965     Double_t tmpDeltaPtKTImprovedCMSPartialExclusion_Signal = 0.0;
1966     if(fNumberSignalJets)
1967       tmpDeltaPtKTImprovedCMSPartialExclusion_Signal = GetDeltaPt(backgroundKTImprovedCMS, 1.0/fNumberSignalJets);
1968     else
1969       tmpDeltaPtKTImprovedCMSPartialExclusion_Signal = GetDeltaPt(backgroundKTImprovedCMS, 1.0);
1970  
1971     Double_t tmpDeltaPtKTImprovedCMSFullExclusion = GetDeltaPt(backgroundKTImprovedCMS, 1.0);
1972
1973     Double_t tmpDeltaPtKTPbPb = 0;
1974     Double_t tmpDeltaPtKTPbPbWithGhosts = 0;
1975     Double_t tmpDeltaPtKTCMS = 0;
1976     Double_t tmpDeltaPtKTMean = 0;
1977     Double_t tmpDeltaPtKTTrackLike = 0;
1978     Double_t tmpDeltaPtTR = 0;
1979
1980     if(fAnalyzeDeprecatedBackgrounds)
1981     {
1982       tmpDeltaPtKTPbPb = GetDeltaPt(backgroundKTPbPb);
1983       tmpDeltaPtKTPbPbWithGhosts = GetDeltaPt(backgroundKTPbPbWithGhosts);
1984       tmpDeltaPtKTCMS = GetDeltaPt(backgroundKTCMS);
1985       tmpDeltaPtKTMean = GetDeltaPt(backgroundKTMean);
1986       tmpDeltaPtKTTrackLike = GetDeltaPt(backgroundKTTrackLike);
1987       tmpDeltaPtTR = GetDeltaPt(backgroundTRCone06);
1988     }
1989
1990     // If valid, fill the delta pt histograms
1991
1992     if(tmpDeltaPtExternalBgrd > -10000.0)
1993       FillHistogram("hDeltaPtExternalBgrd", tmpDeltaPtExternalBgrd, centralityPercentile);
1994     if(tmpDeltaPtKTImprovedCMS > -10000.0)
1995       FillHistogram("hDeltaPtKTImprovedCMS", tmpDeltaPtKTImprovedCMS, centralityPercentile);
1996     if(tmpDeltaPtKTImprovedCMSPartialExclusion > -10000.0)
1997       FillHistogram("hDeltaPtKTImprovedCMSPartialExclusion", tmpDeltaPtKTImprovedCMSPartialExclusion, centralityPercentile);
1998     if(tmpDeltaPtKTImprovedCMSPartialExclusion_Signal > -10000.0)
1999       FillHistogram("hDeltaPtKTImprovedCMSPartialExclusion_Signal", tmpDeltaPtKTImprovedCMSPartialExclusion_Signal, centralityPercentile);
2000     if(tmpDeltaPtKTImprovedCMSFullExclusion > -10000.0)
2001       FillHistogram("hDeltaPtKTImprovedCMSFullExclusion", tmpDeltaPtKTImprovedCMSFullExclusion, centralityPercentile);
2002
2003     if(tmpDeltaPtNoBackground > -10000.0)
2004       FillHistogram("hDeltaPtNoBackground", tmpDeltaPtNoBackground, centralityPercentile);
2005
2006
2007     if(fAnalyzeDeprecatedBackgrounds)
2008     {
2009       if(tmpDeltaPtKTPbPb > -10000.0)
2010         FillHistogram("hDeltaPtKTPbPb", tmpDeltaPtKTPbPb, centralityPercentile);
2011       if(tmpDeltaPtKTPbPbWithGhosts > -10000.0)
2012         FillHistogram("hDeltaPtKTPbPbWithGhosts", tmpDeltaPtKTPbPbWithGhosts, centralityPercentile);
2013       if(tmpDeltaPtKTCMS > -10000.0)
2014         FillHistogram("hDeltaPtKTCMS", tmpDeltaPtKTCMS, centralityPercentile);
2015       if(tmpDeltaPtKTMean > -10000.0)
2016         FillHistogram("hDeltaPtKTMean", tmpDeltaPtKTMean, centralityPercentile);
2017       if(tmpDeltaPtKTTrackLike > -10000.0)
2018         FillHistogram("hDeltaPtKTTrackLike", tmpDeltaPtKTTrackLike, centralityPercentile);
2019       if(tmpDeltaPtTR > -10000.0)
2020         FillHistogram("hDeltaPtTR", tmpDeltaPtTR, centralityPercentile);
2021     }
2022   }
2023   
2024   #ifdef DEBUGMODE
2025     AliInfo("Calculate()::Background done.");
2026   #endif
2027   
2028   #ifdef DEBUGMODE
2029     AliInfo("Calculate() done.");
2030   #endif
2031 }
2032
2033 //________________________________________________________________________
2034 Bool_t AliAnalysisTaskChargedJetsPA::UserNotify()
2035 {
2036   // Implemented Notify() to read the cross sections
2037   // and number of trials from pyxsec.root
2038   // 
2039   #ifdef DEBUGMODE
2040     AliInfo("UserNotify started.");
2041   #endif
2042
2043   if(fAnalyzePythia)
2044   {
2045     TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
2046     TFile *currFile = tree->GetCurrentFile();
2047
2048     TString file(currFile->GetName());
2049
2050     if(file.Contains("root_archive.zip#")){
2051       Ssiz_t pos1 = file.Index("root_archive",12,TString::kExact);
2052       Ssiz_t pos = file.Index("#",1,pos1,TString::kExact);
2053       file.Replace(pos+1,20,"");
2054     }
2055     else {
2056       // not an archive take the basename....
2057       file.ReplaceAll(gSystem->BaseName(file.Data()),"");
2058     }
2059    
2060     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
2061     if(!fxsec){
2062       // next trial fetch the histgram file
2063       fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec_hists.root"));
2064       if(!fxsec){
2065           // not a severe condition but inciate that we have no information
2066         return kFALSE;
2067       }
2068       else{
2069         // find the tlist we want to be independtent of the name so use the Tkey
2070         TKey* key = (TKey*)fxsec->GetListOfKeys()->At(0); 
2071         if(!key){
2072           fxsec->Close();
2073           return kFALSE;
2074         }
2075         TList *list = dynamic_cast<TList*>(key->ReadObj());
2076         if(!list){
2077           fxsec->Close();
2078           return kFALSE;
2079         }
2080         fCrossSection = ((TProfile*)list->FindObject("h1Xsec"))->GetBinContent(1);
2081         fTrials  = ((TH1F*)list->FindObject("h1Trials"))->GetBinContent(1);
2082         fxsec->Close();
2083       }
2084     } // no tree pyxsec.root
2085     else {
2086       TTree *xtree = (TTree*)fxsec->Get("Xsection");
2087       if(!xtree){
2088         fxsec->Close();
2089         return kFALSE;
2090       }
2091       UInt_t   ntrials  = 0;
2092       Double_t  xsection  = 0;
2093       xtree->SetBranchAddress("xsection",&xsection);
2094       xtree->SetBranchAddress("ntrials",&ntrials);
2095       xtree->GetEntry(0);
2096       fTrials = ntrials;
2097       fCrossSection = xsection;
2098       fxsec->Close();
2099     }
2100   }
2101   #ifdef DEBUGMODE
2102     AliInfo("UserNotify ended.");
2103   #endif
2104   return kTRUE;
2105 }
2106
2107
2108 //________________________________________________________________________
2109 inline Double_t AliAnalysisTaskChargedJetsPA::EtaToTheta(Double_t arg)
2110   {return 2.*atan(exp(-arg));} 
2111 //________________________________________________________________________
2112 inline Double_t AliAnalysisTaskChargedJetsPA::ThetaToEta(Double_t arg)
2113 {
2114   if ((arg > TMath::Pi()) || (arg < 0.0))
2115   {
2116     AliError(Form("ThetaToEta got wrong input! (%f)", arg));
2117     return 0.0;
2118   }
2119   return -log(tan(arg/2.));
2120 }
2121 //________________________________________________________________________
2122 inline Double_t AliAnalysisTaskChargedJetsPA::GetDeltaPhi(Double_t phi1, Double_t phi2)
2123   {return min(TMath::Abs(phi1-phi2),TMath::TwoPi() - TMath::Abs(phi1-phi2));}
2124
2125 //________________________________________________________________________
2126 Double_t AliAnalysisTaskChargedJetsPA::MCGetOverlapCircleRectancle(Double_t cPosX, Double_t cPosY, Double_t cRadius, Double_t rPosXmin, Double_t rPosXmax, Double_t rPosYmin, Double_t rPosYmax)
2127 {
2128   const Int_t kTests = 1000;
2129   Int_t hits = 0;
2130   TRandom3 randomGen(0);
2131  
2132   // Loop over kTests-many tests
2133   for (Int_t i=0; i<kTests; i++)
2134   {
2135     //Choose random position in rectangle for the tester
2136     Double_t tmpTestX = randomGen.Uniform(rPosXmin, rPosXmax);
2137     Double_t tmpTestY = randomGen.Uniform(rPosYmin, rPosYmax);
2138
2139     //Check, if tester is in circle. If yes, increment circle counter.
2140     Double_t tmpDistance = TMath::Sqrt( (tmpTestX - cPosX)*(tmpTestX - cPosX) + (tmpTestY - cPosY)*(tmpTestY - cPosY) );
2141     if(tmpDistance < cRadius)
2142       hits++;
2143   }
2144
2145   // return ratio
2146   return (static_cast<Double_t>(hits)/static_cast<Double_t>(kTests));
2147 }
2148
2149 //________________________________________________________________________
2150 Double_t AliAnalysisTaskChargedJetsPA::MCGetOverlapMultipleCirclesRectancle(Int_t numCircles, std::vector<Double_t> cPosX, std::vector<Double_t> cPosY, Double_t cRadius, Double_t rPosXmin, Double_t rPosXmax, Double_t rPosYmin, Double_t rPosYmax)
2151 {
2152
2153   const Int_t kTests = 1000;
2154   Int_t hits = 0;
2155   TRandom3 randomGen(0);
2156  
2157   // Loop over kTests-many tests
2158   for (Int_t i=0; i<kTests; i++)
2159   {
2160     //Choose random position in rectangle for the tester
2161     Double_t tmpTestX = randomGen.Uniform(rPosXmin, rPosXmax);
2162     Double_t tmpTestY = randomGen.Uniform(rPosYmin, rPosYmax);
2163
2164     //Check, if tester is in one of the circles. If yes, increment circle counter.
2165     for(Int_t j=0; j<numCircles; j++)
2166     {
2167       Double_t tmpDistance = TMath::Sqrt( (tmpTestX - cPosX[j])*(tmpTestX - cPosX[j]) + (tmpTestY - cPosY[j])*(tmpTestY - cPosY[j]) );
2168       if(tmpDistance < cRadius)
2169       {
2170         hits++;
2171         break;
2172       }
2173     }
2174   }
2175
2176   // return ratio
2177   return (static_cast<Double_t>(hits)/static_cast<Double_t>(kTests));
2178
2179 }
2180
2181 //________________________________________________________________________
2182 inline void AliAnalysisTaskChargedJetsPA::FillHistogram(const char * key, Double_t x)
2183 {
2184   TH1* tmpHist = static_cast<TH1*>(fOutputList->FindObject(GetHistoName(key)));
2185   if(!tmpHist)
2186   {
2187     AliError(Form("Cannot find histogram <%s> ",key)) ;
2188     return;
2189   }
2190
2191   tmpHist->Fill(x);
2192 }
2193
2194 //________________________________________________________________________
2195 inline void AliAnalysisTaskChargedJetsPA::FillHistogram(const char * key, Double_t x, Double_t y)
2196 {
2197   TH1* tmpHist = static_cast<TH1*>(fOutputList->FindObject(GetHistoName(key)));
2198   if(!tmpHist)
2199   {
2200     AliError(Form("Cannot find histogram <%s> ",key));
2201     return;
2202   }
2203
2204   if (tmpHist->IsA()->GetBaseClass("TH1"))
2205     static_cast<TH1*>(tmpHist)->Fill(x,y); // Fill x with y
2206   else if (tmpHist->IsA()->GetBaseClass("TH2"))
2207     static_cast<TH2*>(tmpHist)->Fill(x,y); // Fill x,y with 1
2208 }
2209
2210 //________________________________________________________________________
2211 inline void AliAnalysisTaskChargedJetsPA::FillHistogram(const char * key, Double_t x, Double_t y, Double_t add)
2212 {
2213   TH2* tmpHist = static_cast<TH2*>(fOutputList->FindObject(GetHistoName(key)));
2214   if(!tmpHist)
2215   {
2216     AliError(Form("Cannot find histogram <%s> ",key));
2217     return;
2218   }
2219   
2220   tmpHist->Fill(x,y,add);
2221 }
2222 //________________________________________________________________________
2223 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)
2224 {
2225   T* tmpHist = new T(GetHistoName(name), GetHistoName(title), xBins, xMin, xMax);
2226
2227   tmpHist->GetXaxis()->SetTitle(xTitle);
2228   tmpHist->GetYaxis()->SetTitle(yTitle);
2229   tmpHist->SetOption(options);
2230   tmpHist->SetMarkerStyle(kFullCircle);
2231   tmpHist->Sumw2();
2232
2233   fHistList->Add(tmpHist);
2234   fHistCount++;
2235   
2236   return tmpHist;
2237 }
2238
2239 //________________________________________________________________________
2240 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)
2241 {
2242   T* tmpHist = new T(GetHistoName(name), GetHistoName(title), xBins, xMin, xMax, yBins, yMin, yMax);
2243   tmpHist->GetXaxis()->SetTitle(xTitle);
2244   tmpHist->GetYaxis()->SetTitle(yTitle);
2245   tmpHist->GetZaxis()->SetTitle(zTitle);
2246   tmpHist->SetOption(options);
2247   tmpHist->SetMarkerStyle(kFullCircle);
2248   tmpHist->Sumw2();
2249
2250   fHistList->Add(tmpHist);
2251   fHistCount++;
2252
2253   return tmpHist;
2254 }
2255
2256 //________________________________________________________________________
2257 void AliAnalysisTaskChargedJetsPA::Terminate(Option_t *)
2258 {
2259   PostData(1, fOutputList);
2260
2261   // Mandatory
2262   fOutputList = dynamic_cast<TList*> (GetOutputData(1)); // '1' refers to the output slot
2263   if (!fOutputList) {
2264     printf("ERROR: Output list not available\n");
2265     return;
2266   }
2267 }
2268
2269 //________________________________________________________________________
2270 AliAnalysisTaskChargedJetsPA::~AliAnalysisTaskChargedJetsPA()
2271 {
2272   // Destructor. Clean-up the output list, but not the histograms that are put inside
2273   // (the list is owner and will clean-up these histograms). Protect in PROOF case.
2274   if (fOutputList && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
2275     delete fOutputList;
2276   }
2277 }
2278
2279 //________________________________________________________________________
2280 void AliAnalysisTaskChargedJetsPA::UserCreateOutputObjects()
2281 {
2282   // called once to create user defined output objects like histograms, plots etc. 
2283   // and to put it on the output list.
2284   // Note: Saving to file with e.g. OpenFile(0) is must be before creating other objects.
2285
2286   fRandom = new TRandom3(0);
2287
2288
2289   fOutputList = new TList();
2290   fOutputList->SetOwner(); // otherwise it produces leaks in merging
2291
2292   PostData(1, fOutputList);
2293 }
2294
2295 //________________________________________________________________________
2296 void AliAnalysisTaskChargedJetsPA::UserExec(Option_t *) 
2297 {
2298   #ifdef DEBUGMODE
2299     AliInfo("UserExec() started.");
2300   #endif
2301
2302   if (!InputEvent())
2303   {
2304     AliError("??? Event pointer == 0 ???");
2305     return;
2306   }
2307
2308   if (!fInitialized)
2309     ExecOnce(); // Get tracks, jets, background from arrays if not already given + Init Histos
2310   
2311   Calculate(InputEvent());
2312         
2313   PostData(1, fOutputList);
2314 }