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