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