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