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