]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskJetV2.h
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / UserTasks / AliAnalysisTaskJetV2.h
1 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. */
2 /* See cxx source for full Copyright notice */
3 /* $Id$ */
4
5 #ifndef AliAnalysisTaskJetV2_H
6 #define AliAnalysisTaskJetV2_H
7
8 #include <AliAnalysisTaskEmcalJet.h>
9 #include <AliEmcalJet.h>
10 #include <AliVEvent.h>
11 #include <AliVParticle.h>
12 #include <AliVCluster.h>
13 #include <TClonesArray.h>
14 #include <TMath.h>
15 #include <TArrayD.h>
16 #include <TRandom3.h>
17 #include <AliJetContainer.h>
18 #include <AliParticleContainer.h>
19
20 class TF1;
21 class THF1;
22 class THF2;
23 class TProfile;
24 class AliLocalRhoParameter;
25 class AliClusterContainer;
26 class AliVTrack;
27
28 class AliAnalysisTaskJetV2 : public AliAnalysisTaskEmcalJet {
29     public:
30          // enumerators
31         enum fitModulationType  { kNoFit, kV2, kV3, kCombined, kFourierSeries, kIntegratedFlow, kQC2, kQC4 }; // fit type
32         enum fitGoodnessTest    { kChi2ROOT, kChi2Poisson, kKolmogorov, kKolmogorovTOY, kLinearFit };
33         enum collisionType      { kPbPb, kPythia };                     // collision type
34         enum qcRecovery         { kFixedRho, kNegativeVn, kTryFit };    // how to deal with negative cn value for qcn value
35         enum runModeType        { kLocal, kGrid };                      // run mode type
36         enum dataType           { kESD, kAOD, kESDMC, kAODMC };         // data type
37         enum detectorType       { kTPC, kVZEROA, kVZEROC, kVZEROComb};  // detector that was used
38         enum analysisType       { kCharged, kFull };                    // analysis type
39         // constructors, destructor
40                                 AliAnalysisTaskJetV2();
41                                 AliAnalysisTaskJetV2(const char *name, runModeType type);
42         virtual                 ~AliAnalysisTaskJetV2();
43         // setting up the task and technical aspects
44         void                    ExecOnce();
45         Bool_t                  InitializeAnalysis();
46         virtual void            UserCreateOutputObjects();
47         virtual Bool_t          Run();
48         TH1F*                   BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Int_t c = -1, Bool_t append = kTRUE);
49         TH2F*                   BookTH2F(const char* name, const char* x, const char* y, Int_t binsx, Double_t minx, Double_t maxx, Int_t binsy, Double_t miny, Double_t maxy, Int_t c = -1, Bool_t append = kTRUE);
50         /* inline */    Double_t PhaseShift(Double_t x) const {  
51             while (x>=TMath::TwoPi())x-=TMath::TwoPi();
52             while (x<0.)x+=TMath::TwoPi();
53             return x; }
54         /* inline */    Double_t PhaseShift(Double_t x, Double_t n) const {
55             x = PhaseShift(x);
56             if(TMath::Nint(n)==2) while (x>TMath::Pi()) x-=TMath::Pi();
57             if(TMath::Nint(n)==3) {
58                 if(x>2.*TMath::TwoPi()/n) x = TMath::TwoPi() - x;
59                 if(x>TMath::TwoPi()/n) x = TMath::TwoPi()-(x+TMath::TwoPi()/n);
60             }
61             return x; }
62         /* inline */    Double_t ChiSquarePDF(Int_t ndf, Double_t x) const {
63             Double_t n(ndf/2.), denom(TMath::Power(2, n)*TMath::Gamma(n));
64             if (denom!=0)  return ((1./denom)*TMath::Power(x, n-1)*TMath::Exp(-x/2.)); 
65             return -999; }
66         // note that the cdf of the chisquare distribution is the normalized lower incomplete gamma function
67         /* inline */    Double_t ChiSquareCDF(Int_t ndf, Double_t x) const { return TMath::Gamma(ndf/2., x/2.); }
68         /* inline */    Double_t ChiSquare(TH1& histo, TF1* func) const {
69             // evaluate the chi2 using a poissonian error estimate on bins
70             Double_t chi2(0.);
71             for(Int_t i(0); i < histo.GetXaxis()->GetNbins(); i++) {
72                 if(histo.GetBinContent(i+1) <= 0.) continue;
73                 chi2 += TMath::Power((histo.GetBinContent(i+1)-func->Eval(histo.GetXaxis()->GetBinCenter(1+i))), 2)/histo.GetBinContent(i+1);
74             }
75            return chi2;
76         }
77         /* inline*/ Double_t KolmogorovTest(TH1F& histo, TF1* func) const {
78             // return the probability from a Kolmogorov test
79             TH1F test(histo);       // stack copy of test statistic
80             for(Int_t i(0); i < test.GetXaxis()->GetNbins(); i++) test.SetBinContent(i+1, func->Eval(test.GetXaxis()->GetBinCenter(1+i)));
81             if(fFitGoodnessTest == kKolmogorovTOY) return histo.TH1::KolmogorovTest((&test), "X");
82             return histo.TH1::KolmogorovTest((&test));
83         }
84  
85         // setters - analysis setup
86         void                    SetDebugMode(Int_t d)                           {fDebug = d;}
87         void                    SetRunToyMC(Bool_t t)                           {fRunToyMC = t; }
88         void                    SetAttachToEvent(Bool_t b)                      {fAttachToEvent = b;}
89         void                    SetFillHistograms(Bool_t b)                     {fFillHistograms = b;}
90         void                    SetFillQAHistograms(Bool_t qa)                  {fFillQAHistograms = qa;}
91         void                    SetReduceBinsXYByFactor(Float_t x, Float_t y)   {fReduceBinsXByFactor = x;
92                                                                                  fReduceBinsYByFactor = y;}
93         void                    SetNoEventWeightsForQC(Bool_t e)                {fNoEventWeightsForQC = e;}
94         void                    SetCentralityClasses(TArrayD* c)                {fCentralityClasses = c;}
95         void                    SetExpectedRuns(TArrayI* r)                     {fExpectedRuns = r;}
96         void                    SetExpectedSemiGoodRuns(TArrayI* r)             {fExpectedSemiGoodRuns = r;}
97         void                    SetIntegratedFlow(TH1F* i, TH1F* j)             {fUserSuppliedV2 = i;
98                                                                                  fUserSuppliedV3 = j; }
99         void                    SetOnTheFlyResCorrection(TH1F* r2, TH1F* r3)    {fUserSuppliedR2 = r2;
100                                                                                  fUserSuppliedR3 = r3; }
101         void                    SetNameRhoSmall(TString name)                   {fNameSmallRho = name; }
102         void                    SetRandomSeed(TRandom3* r)                      {if (fRandom) delete fRandom; fRandom = r; }
103         void                    SetModulationFit(TF1* fit);
104         void                    SetUseControlFit(Bool_t c);
105         void                    SetModulationFitMinMaxP(Float_t m, Float_t n)   {fMinPvalue = m; fMaxPvalue = n; }
106         void                    SetModulationFitType(fitModulationType type)    {fFitModulationType = type; }
107         void                    SetGoodnessTest(fitGoodnessTest test)           {fFitGoodnessTest = test; }
108         void                    SetQCnRecoveryType(qcRecovery type)             {fQCRecovery = type; }
109         void                    SetModulationFitOptions(TString opt)            {fFitModulationOptions = opt; }
110         void                    SetReferenceDetector(detectorType type)         {fDetectorType = type; }
111         void                    SetAnalysisType(analysisType type)              {fAnalysisType = type; }
112         void                    SetCollisionType(collisionType type)            {fCollisionType = type; }
113         void                    SetUsePtWeight(Bool_t w)                        {
114             fUsePtWeight = w; 
115             if(!fUsePtWeight) fUsePtWeightErrorPropagation = kFALSE; }
116         void                    SetUsePtWeightErrorPropagation(Bool_t w)        {fUsePtWeightErrorPropagation = w; }
117         void                    SetRunModeType(runModeType type)                {fRunModeType = type; }
118         void                    SetAbsVertexZ(Float_t v)                        {fAbsVertexZ = v; }
119         void                    SetMinDistanceRctoLJ(Float_t m)                 {fMinDisanceRCtoLJ = m; }
120         void                    SetMaxNoRandomCones(Int_t m)                    {fMaxCones = m; }
121         void                    SetExcludeLeadingJetsFromFit(Float_t n)         {fExcludeLeadingJetsFromFit = n; }
122         void                    SetRebinSwapHistoOnTheFly(Bool_t r)             {fRebinSwapHistoOnTheFly = r; }
123         void                    SetSaveThisPercentageOfFits(Float_t p)          {fPercentageOfFits = p; }
124         void                    SetUseV0EventPlaneFromHeader(Bool_t h)          {fUseV0EventPlaneFromHeader = h;}
125         void                    SetExplicitOutlierCutForYear(Int_t y)           {fExplicitOutlierCut = y;}
126         // getters - these are used as well by AliAnalyisTaskJetFlow, so be careful when changing them
127         TString                 GetJetsName() const                             {return GetJetContainer()->GetArrayName(); }
128         TString                 GetTracksName() const                           {return GetParticleContainer()->GetArrayName(); }
129         TString                 GetLocalRhoName() const                         {return fLocalRhoName; }
130         TArrayD*                GetCentralityClasses() const                    {return fCentralityClasses;}
131         TProfile*               GetResolutionParameters(Int_t h, Int_t c) const {return (h==2) ? fProfV2Resolution[c] : fProfV3Resolution[c];}
132         TList*                  GetOutputList() const                           {return fOutputList;}
133         AliLocalRhoParameter*   GetLocalRhoParameter() const                    {return fLocalRho;}
134         Double_t                GetJetRadius() const                            {return GetJetContainer()->GetJetRadius();}
135         /* inline */    AliEmcalJet* GetLeadingJet() {
136             // return pointer to the highest pt jet (before background subtraction) within acceptance
137             // only rudimentary cuts are applied on this level, hence the implementation outside of
138             // the framework
139             Int_t iJets(fJets->GetEntriesFast());
140             Double_t pt(0);
141             AliEmcalJet* leadingJet(0x0);
142             for(Int_t i(0); i < iJets; i++) {
143                 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
144                 if(!PassesSimpleCuts(jet)) continue;
145                 if(jet->Pt() > pt) {
146                    leadingJet = jet;
147                    pt = leadingJet->Pt();
148                 }
149             }
150             return leadingJet;
151         }
152         void                    ExecMe()                                {ExecOnce();}
153         AliAnalysisTaskJetV2* ReturnMe()                                {return this;}
154         // local cuts
155         void                    SetSoftTrackMinMaxPt(Float_t min, Float_t max)          {fSoftTrackMinPt = min; fSoftTrackMaxPt = max;}
156         void                    SetSemiGoodJetMinMaxPhi(Double_t a, Double_t b)         {fSemiGoodJetMinPhi = a; fSemiGoodJetMaxPhi = b;}
157         void                    SetSemiGoodTrackMinMaxPhi(Double_t a, Double_t b)       {fSemiGoodTrackMinPhi = a; fSemiGoodTrackMaxPhi = b;}
158         // numerical evaluations
159         void                    CalculateEventPlaneVZERO(Double_t vzero[2][2]) const;
160         void                    CalculateEventPlaneTPC(Double_t* tpc);
161         void                    CalculateEventPlaneCombinedVZERO(Double_t* comb) const;
162         void                    CalculateEventPlaneResolution(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc);
163         void                    CalculateRandomCone(
164                 Float_t &pt, 
165                 Float_t &eta, 
166                 Float_t &phi, 
167                 AliParticleContainer* tracksCont,
168                 AliClusterContainer* clusterCont = 0x0,
169                 AliEmcalJet* jet = 0x0
170                 ) const;
171         Double_t                CalculateQC2(Int_t harm);
172         Double_t                CalculateQC4(Int_t harm);
173         // helper calculations for the q-cumulant analysis, also used by AliAnalyisTaskJetFlow
174         void                    QCnQnk(Int_t n, Int_t k, Double_t &reQ, Double_t &imQ);
175         void                    QCnDiffentialFlowVectors(
176             TClonesArray* pois, TArrayD* ptBins, Bool_t vpart, Double_t* repn, Double_t* impn, 
177             Double_t *mp, Double_t *reqn, Double_t *imqn, Double_t* mq, Int_t n);
178         Double_t                QCnS(Int_t i, Int_t j);
179         Double_t                QCnM();
180         Double_t                QCnM11();
181         Double_t                QCnM1111();
182         Bool_t                  QCnRecovery(Double_t psi2, Double_t psi3);
183         // analysis details
184         Bool_t                  CorrectRho(Double_t psi2, Double_t psi3);
185         // event and track selection, also used by AliAnalyisTaskJetFlow
186         /* inline */    Bool_t PassesCuts(AliVParticle* track) const    { return AcceptTrack(track, 0); }
187         /* inline */    Bool_t PassesCuts(AliEmcalJet* jet)             { return AcceptJet(jet, 0); }
188         /* inline */    Bool_t PassesCuts(AliVCluster* clus) const      { return AcceptCluster(clus, 0); }
189         /* inline */    Bool_t PassesSimpleCuts(AliEmcalJet* jet)       {
190             Float_t minPhi(GetJetContainer()->GetJetPhiMin()), maxPhi(GetJetContainer()->GetJetPhiMax());
191             Float_t minEta(GetJetContainer()->GetJetEtaMin()), maxEta(GetJetContainer()->GetJetEtaMax());
192             return (jet && jet->Pt() > 1. && jet->Eta() > minEta && jet->Eta() < maxEta && jet->Phi() > minPhi && jet->Phi() < maxPhi && jet->Area() > .557*GetJetRadius()*GetJetRadius()*TMath::Pi());
193         }
194         Bool_t                  PassesCuts(AliVEvent* event);
195         Bool_t                  PassesCuts(Int_t year);
196         Bool_t                  PassesCuts(const AliVCluster* track) const;
197         // filling histograms
198         void                    FillHistogramsAfterSubtraction(Double_t psi2, Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc);
199         void                    FillTrackHistograms() const;
200         void                    FillClusterHistograms() const;
201         void                    FillEventPlaneHistograms(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc) const;
202         void                    FillRhoHistograms();
203         void                    FillDeltaPtHistograms(Double_t psi2) const; 
204         void                    FillJetHistograms(Double_t psi2);
205         void                    FillQAHistograms(AliVTrack* vtrack) const;
206         void                    FillQAHistograms(AliVEvent* vevent);
207         void                    FillAnalysisSummaryHistogram() const;
208         virtual void            Terminate(Option_t* option);
209         // interface methods for the output file
210         void                    SetOutputList(TList* l) {fOutputList = l;}
211         TH1F*                   GetResolutionFromOuptutFile(detectorType detector, Int_t h = 2, TArrayD* c = 0x0);
212         TH1F*                   CorrectForResolutionDiff(TH1F* v, detectorType detector, TArrayD* cen, Int_t c, Int_t h = 2);
213         TH1F*                   CorrectForResolutionInt(TH1F* v, detectorType detector, TArrayD* cen, Int_t h = 2);
214         TH1F*                   GetDifferentialQC(TProfile* refCumulants, TProfile* diffCumlants, TArrayD* ptBins, Int_t h);
215     private:
216         // analysis flags and settings
217         Int_t                   fDebug;                 // debug level (0 none, 1 fcn calls, 2 verbose)
218         Bool_t                  fRunToyMC;              // run toy mc for fit routine
219         Bool_t                  fLocalInit;             //! is the analysis initialized?
220         Bool_t                  fAttachToEvent;         // attach local rho to the event
221         Bool_t                  fFillHistograms;        // fill histograms
222         Bool_t                  fFillQAHistograms;      // fill qa histograms
223         Float_t                 fReduceBinsXByFactor;   // reduce the bins on x-axis of histo's by this much
224         Float_t                 fReduceBinsYByFactor;   // reduce the bins on y-axis of histo's by this much
225         Bool_t                  fNoEventWeightsForQC;   // don't store event weights for qc analysis
226         TArrayD*                fCentralityClasses;     //-> centrality classes (maximum 10)
227         TArrayI*                fExpectedRuns;          //-> array of expected run numbers, used for QA
228         TArrayI*                fExpectedSemiGoodRuns;  //-> array of expected semi-good runs, used for cuts and QA
229         TH1F*                   fUserSuppliedV2;        // histo with integrated v2
230         TH1F*                   fUserSuppliedV3;        // histo with integrated v3
231         TH1F*                   fUserSuppliedR2;        // correct the extracted v2 with this r
232         TH1F*                   fUserSuppliedR3;        // correct the extracted v3 with this r
233         AliParticleContainer*   fTracksCont;            //! tracks
234         AliClusterContainer*    fClusterCont;           //! cluster container
235         AliJetContainer*        fJetsCont;              //! jets
236         AliEmcalJet*            fLeadingJet;            //! leading jet
237         // members
238         Int_t                   fNAcceptedTracks;       //! number of accepted tracks
239         Int_t                   fNAcceptedTracksQCn;    //! accepted tracks for QCn
240         fitModulationType       fFitModulationType;     // fit modulation type
241         fitGoodnessTest         fFitGoodnessTest;       // fit goodness test type
242         qcRecovery              fQCRecovery;            // recovery type for e-by-e qc method
243         Bool_t                  fUsePtWeight;           // use dptdphi instead of dndphi
244         Bool_t                  fUsePtWeightErrorPropagation;   // recalculate the bin errors in case of pt weighting 
245         detectorType            fDetectorType;          // type of detector used for modulation fit
246         analysisType            fAnalysisType;          // analysis type (full or charged jets)
247         TString                 fFitModulationOptions;  // fit options for modulation fit
248         runModeType             fRunModeType;           // run mode type 
249         dataType                fDataType;              // datatype 
250         collisionType           fCollisionType;         // collision type
251         TRandom3*               fRandom;                //-> dont use gRandom to not interfere with other tasks
252         Int_t                   fRunNumber;             //! current runnumber (for QA and jet, track selection)
253         Int_t                   fMappedRunNumber;       //! mapped runnumer (for QA)
254         Int_t                   fInCentralitySelection; //! centrality bin
255         TF1*                    fFitModulation;         //-> modulation fit for rho
256         TF1*                    fFitControl;            //-> control fit
257         Float_t                 fMinPvalue;             // minimum value of p
258         Float_t                 fMaxPvalue;             // maximum value of p
259         TString                 fNameSmallRho;          // name of small rho
260         AliRhoParameter*        fCachedRho;             //! temp cache for rho pointer
261         // additional jet cuts (most are inherited)
262         Float_t                 fSoftTrackMinPt;        // min pt for soft tracks
263         Float_t                 fSoftTrackMaxPt;        // max pt for soft tracks
264         Double_t                fSemiGoodJetMinPhi;     // min phi for semi good tpc runs
265         Double_t                fSemiGoodJetMaxPhi;     // max phi for semi good tpc runs
266         Double_t                fSemiGoodTrackMinPhi;   // min phi for semi good tpc runs
267         Double_t                fSemiGoodTrackMaxPhi;   // max phi for semi good tpc runs
268         // event cuts
269         Float_t                 fAbsVertexZ;            // cut on zvertex
270         // general qa histograms
271         TH1F*                   fHistCentrality;        //! accepted centrality
272         TH1F*                   fHistVertexz;           //! accepted verte
273         TH2F*                   fHistRunnumbersPhi;     //! run numbers averaged phi
274         TH2F*                   fHistRunnumbersEta;     //! run numbers averaged eta
275         TH1F*                   fHistPvalueCDFROOT;     //! pdf value of chisquare p
276         TH2F*                   fHistPvalueCDFROOTCent; //! p value versus centrlaity from root
277         TH2F*                   fHistChi2ROOTCent;      //! reduced chi2 from ROOT, centrality correlation
278         TH2F*                   fHistPChi2Root;         //! correlation p value and reduced chi2
279         TH1F*                   fHistPvalueCDF;         //! cdf value of chisquare p
280         TH2F*                   fHistPvalueCDFCent;     //! p value vs centrality
281         TH2F*                   fHistChi2Cent;          //! reduced chi2, centrlaity correlation
282         TH2F*                   fHistPChi2;             //! correlation p value and reduced chi2
283         TH1F*                   fHistKolmogorovTest;    //! KolmogorovTest value
284         TH2F*                   fHistKolmogorovTestCent;//! KolmogorovTest value, centrality correlation
285         TH2F*                   fHistPKolmogorov;       //! p value vs kolmogorov value
286         TH2F*                   fHistRhoStatusCent;     //! status of rho as function of centrality
287         TH1F*                   fHistUndeterminedRunQA; //! undetermined run QA
288         // general settings
289         Float_t                 fMinDisanceRCtoLJ;      // min distance between rc and leading jet
290         Int_t                   fMaxCones;              // max number of random cones
291         Float_t                 fExcludeLeadingJetsFromFit;    // exclude n leading jets from fit
292         Bool_t                  fRebinSwapHistoOnTheFly;       // rebin swap histo on the fly
293         Float_t                 fPercentageOfFits;      // save this percentage of fits
294         Bool_t                  fUseV0EventPlaneFromHeader;    // use the vzero event plane from the header
295         Int_t                   fExplicitOutlierCut;    // cut on correlation of tpc and global multiplicity
296         // transient object pointers
297         TList*                  fOutputList;            //! output list
298         TList*                  fOutputListGood;        //! output list for local analysis
299         TList*                  fOutputListBad;         //! output list for local analysis
300         TH1F*                   fHistAnalysisSummary;   //! analysis summary
301         TH1F*                   fHistSwap;              //! swap histogram
302         TProfile*               fProfV2;                //! extracted v2
303         TProfile*               fProfV2Cumulant;        //! v2 cumulant
304         TProfile*               fProfV2Resolution[10];  //! resolution parameters for v2
305         TProfile*               fProfV3;                //! extracted v3
306         TProfile*               fProfV3Cumulant;        //! v3 cumulant
307         TProfile*               fProfV3Resolution[10];  //! resolution parameters for v3
308         // qa histograms for accepted pico tracks
309         TH1F*                   fHistPicoTrackPt[10];   //! pt of all charged tracks
310         TH1F*                   fHistPicoTrackMult[10]; //! multiplicity of accepted pico tracks
311         TH2F*                   fHistPicoCat1[10];      //! pico tracks spd hit and refit
312         TH2F*                   fHistPicoCat2[10];      //! pico tracks wo spd hit w refit, constrained
313         TH2F*                   fHistPicoCat3[10];      //! pico tracks wo spd hit wo refit, constrained
314         // qa histograms for accepted emcal clusters
315         TH1F*                   fHistClusterPt[10];     //! pt emcal clusters
316         TH2F*                   fHistClusterEtaPhi[10]; //! eta phi emcal clusters
317         TH2F*                   fHistClusterEtaPhiWeighted[10]; //! eta phi emcal clusters, pt weighted
318         // qa event planes
319         TProfile*               fHistPsiControl;        //! event plane control histogram
320         TProfile*               fHistPsiSpread;         //! event plane spread histogram
321         TH1F*                   fHistPsiVZEROA;         //! psi 2 from vzero a
322         TH1F*                   fHistPsiVZEROC;         //! psi 2 from vzero c
323         TH1F*                   fHistPsiVZERO;          //! psi 2 from combined vzero
324         TH1F*                   fHistPsiTPC;            //! psi 2 from tpc
325         TH2F*                   fHistPsiVZEROAV0M;      //! psi 2 from vzero a
326         TH2F*                   fHistPsiVZEROCV0M;      //! psi 2 from vzero c
327         TH2F*                   fHistPsiVZEROVV0M;      //! psi 2 from combined vzero
328         TH2F*                   fHistPsiTPCiV0M;        //! psi 2 from tpc
329         TH2F*                   fHistPsiVZEROATRK;      //! psi 2 from vzero a
330         TH2F*                   fHistPsiVZEROCTRK;      //! psi 2 from vzero c
331         TH2F*                   fHistPsiVZEROTRK;       //! psi 2 from combined vzero
332         TH2F*                   fHistPsiTPCTRK;         //! psi 2 from tpc
333         // background
334         TH1F*                   fHistRhoPackage[10];    //! rho as estimated by emcal jet package
335         TH1F*                   fHistRho[10];           //! background
336         TH2F*                   fHistRhoVsMult;         //! rho versus multiplicity
337         TH2F*                   fHistRhoVsCent;         //! rho veruss centrality
338         TH2F*                   fHistRhoAVsMult;        //! rho * A vs multiplicity for all jets
339         TH2F*                   fHistRhoAVsCent;        //! rho * A vs centrality for all jets
340         // delta pt distributions
341         TH2F*                   fHistRCPhiEta[10];              //! random cone eta and phi
342         TH2F*                   fHistRhoVsRCPt[10];             //! rho * A vs rcpt
343         TH1F*                   fHistRCPt[10];                  //! rcpt
344         TH2F*                   fHistDeltaPtDeltaPhi2[10];      //! dpt vs dphi (psi2 - phi)
345         TH2F*                   fHistDeltaPtDeltaPhi2Rho0[10];  //! dpt vs dphi, rho_0
346         TH2F*                   fHistRCPhiEtaExLJ[10];          //! random cone eta and phi, excl leading jet
347         TH2F*                   fHistRhoVsRCPtExLJ[10];         //! rho * A vs rcpt, excl leading jet
348         TH1F*                   fHistRCPtExLJ[10];              //! rcpt, excl leading jet
349         TH2F*                   fHistDeltaPtDeltaPhi2ExLJ[10];  //! dpt vs dphi, excl leading jet
350         TH2F*                   fHistDeltaPtDeltaPhi2ExLJRho0[10];      //! dpt vs dphi, excl leading jet, rho_0
351         // jet histograms (after kinematic cuts)
352         TH1F*                   fHistJetPtRaw[10];              //! jet pt - no background subtraction
353         TH1F*                   fHistJetPt[10];                 //! pt of found jets (background subtracted)
354         TH2F*                   fHistJetEtaPhi[10];             //! eta and phi correlation
355         TH2F*                   fHistJetPtArea[10];             //! jet pt versus area
356         TH2F*                   fHistJetPtEta[10];              //! jet pt versus eta (temp control)
357         TH2F*                   fHistJetPtConstituents[10];     //! jet pt versus number of constituents
358         TH2F*                   fHistJetEtaRho[10];             //! jet eta versus rho
359         // in plane, out of plane jet spectra
360         TH2F*                   fHistJetPsi2Pt[10];             //! event plane dependence of jet pt
361         TH2F*                   fHistJetPsi2PtRho0[10];         //! event plane dependence of jet pt vs rho_0
362
363         AliAnalysisTaskJetV2(const AliAnalysisTaskJetV2&);                  // not implemented
364         AliAnalysisTaskJetV2& operator=(const AliAnalysisTaskJetV2&);       // not implemented
365
366         ClassDef(AliAnalysisTaskJetV2, 1);
367 };
368
369 #endif