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 SetEventPlaneWeights(TH1F* ep) {fEventPlaneWeights = ep; }
108 void SetNameRhoSmall(TString name) {fNameSmallRho = name; }
109 void SetRandomSeed(TRandom3* r) {if (fRandom) delete fRandom; fRandom = r; }
110 void SetModulationFit(TF1* fit);
111 void SetUseControlFit(Bool_t c);
112 void SetModulationFitMinMaxP(Float_t m, Float_t n) {fMinPvalue = m; fMaxPvalue = n; }
113 void SetModulationFitType(fitModulationType type) {fFitModulationType = type; }
114 void SetGoodnessTest(fitGoodnessTest test) {fFitGoodnessTest = test; }
115 void SetQCnRecoveryType(qcRecovery type) {fQCRecovery = type; }
116 void SetModulationFitOptions(TString opt) {fFitModulationOptions = opt; }
117 void SetReferenceDetector(detectorType type) {fDetectorType = type; }
118 void SetAnalysisType(analysisType type) {fAnalysisType = type; }
119 void SetCollisionType(collisionType type) {fCollisionType = type; }
120 void SetUsePtWeight(Bool_t w) {
122 if(!fUsePtWeight) fUsePtWeightErrorPropagation = kFALSE; }
123 void SetUsePtWeightErrorPropagation(Bool_t w) {fUsePtWeightErrorPropagation = w; }
124 void SetRunModeType(runModeType type) {fRunModeType = type; }
125 void SetAbsVertexZ(Float_t v) {fAbsVertexZ = v; }
126 void SetMinDistanceRctoLJ(Float_t m) {fMinDisanceRCtoLJ = m; }
127 void SetMaxNoRandomCones(Int_t m) {fMaxCones = m; }
128 void SetExcludeLeadingJetsFromFit(Float_t n) {fExcludeLeadingJetsFromFit = n; }
129 void SetRebinSwapHistoOnTheFly(Bool_t r) {fRebinSwapHistoOnTheFly = r; }
130 void SetSaveThisPercentageOfFits(Float_t p) {fPercentageOfFits = p; }
131 // setters specific to the vzero calibration for 10h data
132 void SetVZEROApol(Int_t ring, Float_t f) {fVZEROApol[ring]=f;}
133 void SetVZEROCpol(Int_t ring, Float_t f) {fVZEROCpol[ring]=f;}
134 void SetVZEROgainEqualizationPerRing(Bool_t s) {fVZEROgainEqualizationPerRing = s;}
135 void SetUseVZERORing(Int_t i, Bool_t u) {
136 // exclude vzero rings: 0 through 7 can be excluded by calling this setter multiple times
137 // 0 corresponds to segment ID 0 through 7, etc
138 fUseVZERORing[i] = u;
139 fVZEROgainEqualizationPerRing = kTRUE; // must be true for this option
142 void SetChi2VZEROA(TArrayD* a) { fChi2A = a;}
143 void SetChi2VZEROC(TArrayD* a) { fChi2C = a;}
144 void SetChi3VZEROA(TArrayD* a) { fChi3A = a;}
145 void SetChi3VZEROC(TArrayD* a) { fChi3C = a;}
148 TString GetJetsName() const {return GetJetContainer()->GetArrayName(); }
149 TString GetTracksName() const {return GetParticleContainer()->GetArrayName(); }
150 TString GetLocalRhoName() const {return fLocalRhoName; }
151 TArrayD* GetCentralityClasses() const {return fCentralityClasses;}
152 TProfile* GetResolutionParameters(Int_t h, Int_t c) const {return (h==2) ? fProfV2Resolution[c] : fProfV3Resolution[c];}
153 TList* GetOutputList() const {return fOutputList;}
154 AliLocalRhoParameter* GetLocalRhoParameter() const {return fLocalRho;}
155 Double_t GetJetRadius() const {return GetJetContainer()->GetJetRadius();}
156 AliEmcalJet* GetLeadingJet(AliLocalRhoParameter* localRho = 0x0);
157 void ExecMe() {ExecOnce();}
159 AliAnalysisTaskJetV2* ReturnMe() {return this;}
161 void SetSoftTrackMinMaxPt(Float_t min, Float_t max) {fSoftTrackMinPt = min; fSoftTrackMaxPt = max;}
162 void SetSemiGoodJetMinMaxPhi(Double_t a, Double_t b) {fSemiGoodJetMinPhi = a; fSemiGoodJetMaxPhi = b;}
163 void SetSemiGoodTrackMinMaxPhi(Double_t a, Double_t b) {fSemiGoodTrackMinPhi = a; fSemiGoodTrackMaxPhi = b;}
164 // numerical evaluations
165 static Double_t CalculateEventPlaneChi(Double_t res);
166 void CalculateEventPlaneVZERO(Double_t vzero[2][2]) const;
167 void CalculateEventPlaneCombinedVZERO(Double_t* comb) const;
168 void CalculateEventPlaneTPC(Double_t* tpc);
169 void CalculateEventPlaneResolution(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc);
170 void CalculateQvectorVZERO(Double_t Qa2[2], Double_t Qc2[2], Double_t Qa3[2], Double_t Qc3[2]) const;
171 void CalculateQvectorCombinedVZERO(Double_t Q2[2], Double_t Q3[2]) const;
172 void CalculateRandomCone(
176 AliParticleContainer* tracksCont,
177 AliClusterContainer* clusterCont = 0x0,
178 AliEmcalJet* jet = 0x0
180 Double_t CalculateQC2(Int_t harm);
181 Double_t CalculateQC4(Int_t harm);
182 // helper calculations for the q-cumulant analysis
183 void QCnQnk(Int_t n, Int_t k, Double_t &reQ, Double_t &imQ);
184 void QCnDiffentialFlowVectors(
185 TClonesArray* pois, TArrayD* ptBins, Bool_t vpart, Double_t* repn, Double_t* impn,
186 Double_t *mp, Double_t *reqn, Double_t *imqn, Double_t* mq, Int_t n);
187 Double_t QCnS(Int_t i, Int_t j);
191 Bool_t QCnRecovery(Double_t psi2, Double_t psi3);
193 Bool_t CorrectRho(Double_t psi2, Double_t psi3);
194 // event and track selection
195 /* inline */ Bool_t PassesCuts(AliVParticle* track) const { return AcceptTrack(track, 0); }
196 /* inline */ Bool_t PassesCuts(AliEmcalJet* jet) { return AcceptJet(jet, 0); }
197 /* inline */ Bool_t PassesCuts(AliVCluster* clus) const { return AcceptCluster(clus, 0); }
198 /* inline */ Bool_t PassesSimpleCuts(AliEmcalJet* jet) {
199 Float_t minPhi(GetJetContainer()->GetJetPhiMin()), maxPhi(GetJetContainer()->GetJetPhiMax());
200 Float_t minEta(GetJetContainer()->GetJetEtaMin()), maxEta(GetJetContainer()->GetJetEtaMax());
201 return (jet/* && jet->Pt() > 1.*/ && jet->Eta() > minEta && jet->Eta() < maxEta && jet->Phi() > minPhi && jet->Phi() < maxPhi && jet->Area() > .557*GetJetRadius()*GetJetRadius()*TMath::Pi());
203 Bool_t PassesCuts(AliVEvent* event);
204 Bool_t PassesCuts(const AliVCluster* track) const;
205 // filling histograms
206 void FillHistogramsAfterSubtraction(Double_t psi2, Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc);
207 void FillTrackHistograms() const;
208 void FillClusterHistograms() const;
209 void FillEventPlaneHistograms(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc) const;
210 void FillRhoHistograms();
211 void FillDeltaPtHistograms(Double_t psi2) const;
212 void FillJetHistograms(Double_t psi2);
213 void FillQAHistograms(AliVTrack* vtrack) const;
214 void FillQAHistograms(AliVEvent* vevent);
215 void FillWeightedTrackHistograms() const;
216 void FillWeightedClusterHistograms() const;
217 void FillWeightedEventPlaneHistograms(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc) const;
218 void FillWeightedRhoHistograms();
219 void FillWeightedDeltaPtHistograms(Double_t psi2) const;
220 void FillWeightedJetHistograms(Double_t psi2);
221 void FillWeightedQAHistograms(AliVTrack* vtrack) const;
222 void FillWeightedQAHistograms(AliVEvent* vevent);
223 void FillAnalysisSummaryHistogram() const;
224 virtual void Terminate(Option_t* option);
225 // interface methods for the output file
226 void SetOutputList(TList* l) {fOutputList = l;}
227 TH1F* GetResolutionFromOuptutFile(detectorType detector, Int_t h = 2, TArrayD* c = 0x0);
228 TH1F* CorrectForResolutionDiff(TH1F* v, detectorType detector, TArrayD* cen, Int_t c, Int_t h = 2);
229 TH1F* CorrectForResolutionInt(TH1F* v, detectorType detector, TArrayD* cen, Int_t h = 2);
230 TH1F* GetDifferentialQC(TProfile* refCumulants, TProfile* diffCumlants, TArrayD* ptBins, Int_t h);
231 void ReadVZEROCalibration2010h();
232 Int_t GetVZEROCentralityBin() const;
234 // analysis flags and settings
235 Int_t fDebug; // debug level (0 none, 1 fcn calls, 2 verbose)
236 Bool_t fRunToyMC; // run toy mc for fit routine
237 Bool_t fLocalInit; //! is the analysis initialized?
238 Bool_t fAttachToEvent; // attach local rho to the event
239 Bool_t fFillHistograms; // fill histograms
240 Bool_t fFillQAHistograms; // fill qa histograms
241 Float_t fReduceBinsXByFactor; // reduce the bins on x-axis of histo's by this much
242 Float_t fReduceBinsYByFactor; // reduce the bins on y-axis of histo's by this much
243 Bool_t fNoEventWeightsForQC; // don't store event weights for qc analysis
244 TArrayD* fCentralityClasses; //-> centrality classes (maximum 10)
245 TArrayI* fExpectedRuns; //-> array of expected run numbers, used for QA
246 TArrayI* fExpectedSemiGoodRuns; //-> array of expected semi-good runs, used for cuts and QA
247 TH1F* fUserSuppliedV2; // histo with integrated v2
248 TH1F* fUserSuppliedV3; // histo with integrated v3
249 TH1F* fUserSuppliedR2; // correct the extracted v2 with this r
250 TH1F* fUserSuppliedR3; // correct the extracted v3 with this r
251 TH1F* fEventPlaneWeights; // weight histo for the event plane
252 Float_t fEventPlaneWeight; //! the actual weight of an event
253 AliParticleContainer* fTracksCont; //! tracks
254 AliClusterContainer* fClusterCont; //! cluster container
255 AliJetContainer* fJetsCont; //! jets
256 AliEmcalJet* fLeadingJet; //! leading jet
257 AliEmcalJet* fLeadingJetAfterSub; //! leading jet after background subtraction
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
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
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
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
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* fHistPsiTPCV0M; //! 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 TH3F* fHistPsiTPCLeadingJet[10]; //! correlation tpc EP, LJ pt
353 TH3F* fHistPsiVZEROALeadingJet[10]; //! correlation vzeroa EP, LJ pt
354 TH3F* fHistPsiVZEROCLeadingJet[10]; //! correlation vzeroc EP, LJ pt
355 TH3F* fHistPsiVZEROCombLeadingJet[10];//! correlation vzerocomb EP, LJ pt
356 TH3F* fHistPsi2Correlation[10]; //! correlation of event planes
357 TH2F* fHistLeadingJetBackground[10]; //! geometric correlation of leading jet w/wo bkg subtraction
359 TH1F* fHistRhoPackage[10]; //! rho as estimated by emcal jet package
360 TH1F* fHistRho[10]; //! background
361 TH2F* fHistRhoVsMult; //! rho versus multiplicity
362 TH2F* fHistRhoVsCent; //! rho veruss centrality
363 TH2F* fHistRhoAVsMult; //! rho * A vs multiplicity for all jets
364 TH2F* fHistRhoAVsCent; //! rho * A vs centrality for all jets
365 // delta pt distributions
366 TH2F* fHistRCPhiEta[10]; //! random cone eta and phi
367 TH2F* fHistRhoVsRCPt[10]; //! rho * A vs rcpt
368 TH1F* fHistRCPt[10]; //! rcpt
369 TH2F* fHistDeltaPtDeltaPhi2[10]; //! dpt vs dphi (psi2 - phi)
370 TH2F* fHistDeltaPtDeltaPhi2Rho0[10]; //! dpt vs dphi, rho_0
371 TH2F* fHistRCPhiEtaExLJ[10]; //! random cone eta and phi, excl leading jet
372 TH2F* fHistRhoVsRCPtExLJ[10]; //! rho * A vs rcpt, excl leading jet
373 TH1F* fHistRCPtExLJ[10]; //! rcpt, excl leading jet
374 TH2F* fHistDeltaPtDeltaPhi2ExLJ[10]; //! dpt vs dphi, excl leading jet
375 TH2F* fHistDeltaPtDeltaPhi2ExLJRho0[10]; //! dpt vs dphi, excl leading jet, rho_0
376 // jet histograms (after kinematic cuts)
377 TH1F* fHistJetPtRaw[10]; //! jet pt - no background subtraction
378 TH1F* fHistJetPt[10]; //! pt of found jets (background subtracted)
379 TH2F* fHistJetEtaPhi[10]; //! eta and phi correlation
380 TH2F* fHistJetPtArea[10]; //! jet pt versus area
381 TH2F* fHistJetPtEta[10]; //! jet pt versus eta (temp control)
382 TH2F* fHistJetPtConstituents[10]; //! jet pt versus number of constituents
383 TH2F* fHistJetEtaRho[10]; //! jet eta versus rho
384 // in plane, out of plane jet spectra
385 TH2F* fHistJetPsi2Pt[10]; //! event plane dependence of jet pt
386 TH2F* fHistJetPsi2PtRho0[10]; //! event plane dependence of jet pt vs rho_0
387 // vzero event plane calibration cache for 10h data
388 Float_t fMeanQ[9][2][2]; //! recentering
389 Float_t fWidthQ[9][2][2]; //! recentering
390 Float_t fMeanQv3[9][2][2]; //! recentering
391 Float_t fWidthQv3[9][2][2]; //! recentering
392 TH1* fVZEROgainEqualization; //! equalization histo
393 Bool_t fVZEROgainEqualizationPerRing; // per ring vzero gain calibration
394 Float_t fVZEROApol[4]; //! calibration info per ring
395 Float_t fVZEROCpol[4]; //! calibration info per ring
396 Bool_t fUseVZERORing[8]; // kTRUE means the ring is included
397 TArrayD* fChi2A; // chi vs cent for vzero A ep_2
398 TArrayD* fChi2C; // chi vs cent for vzero C ep_2
399 TArrayD* fChi3A; // chi vs cent for vzero A ep_3
400 TArrayD* fChi3C; // chi vs cent for vzero C ep_3
401 TFile* fOADB; //! fOADB
404 AliAnalysisTaskJetV2(const AliAnalysisTaskJetV2&); // not implemented
405 AliAnalysisTaskJetV2& operator=(const AliAnalysisTaskJetV2&); // not implemented
407 ClassDef(AliAnalysisTaskJetV2, 3);