1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
19 * this task is part of the emcal jet framework and should be run in the emcaljet train
20 * the following extensions to an accepted AliVEvent are expected:
22 * - background estimate rho
24 * aod's and esd's are handled transparently
25 * the task will attempt to estimate a phi-dependent background density rho
26 * by fitting vn harmonics to the dpt/dphi distribution
28 * author: Redmer Alexander Bertens, Utrecht Univeristy, Utrecht, Netherlands
29 * rbertens@cern.ch, rbertens@nikhef.nl, r.a.bertens@uu.nl
45 #include <AliAnalysisTask.h>
46 #include <AliAnalysisManager.h>
47 #include <AliCentrality.h>
48 #include <AliVVertex.h>
49 #include <AliVTrack.h>
50 #include <AliVVZERO.h>
51 #include <AliESDEvent.h>
52 #include <AliAODEvent.h>
53 #include <AliAODTrack.h>
54 #include <AliOADBContainer.h>
55 // emcal jet framework includes
56 #include <AliPicoTrack.h>
57 #include <AliEmcalJet.h>
58 #include <AliRhoParameter.h>
59 #include <AliLocalRhoParameter.h>
60 #include <AliAnalysisTaskJetV2.h>
61 #include <AliClusterContainer.h>
63 class AliAnalysisTaskJetV2;
66 ClassImp(AliAnalysisTaskJetV2)
68 AliAnalysisTaskJetV2::AliAnalysisTaskJetV2() : AliAnalysisTaskEmcalJet("AliAnalysisTaskJetV2", kTRUE),
69 fDebug(0), fRunToyMC(kFALSE), fLocalInit(0), fAttachToEvent(kTRUE), fFillHistograms(kTRUE), fFillQAHistograms(kTRUE), fReduceBinsXByFactor(-1.), fReduceBinsYByFactor(-1.), fNoEventWeightsForQC(kTRUE), fCentralityClasses(0), fExpectedRuns(0), fExpectedSemiGoodRuns(0), fUserSuppliedV2(0), fUserSuppliedV3(0), fUserSuppliedR2(0), fUserSuppliedR3(0), fTracksCont(0), fClusterCont(0), fJetsCont(0), fLeadingJet(0), fNAcceptedTracks(0), fNAcceptedTracksQCn(0), fFitModulationType(kNoFit), fFitGoodnessTest(kChi2Poisson), fQCRecovery(kTryFit), fUsePtWeight(kTRUE), fUsePtWeightErrorPropagation(kTRUE), fDetectorType(kVZEROComb), fAnalysisType( kCharged), fFitModulationOptions("QWLI"), fRunModeType(kGrid), fDataType(kESD), fCollisionType(kPbPb), fRandom(0), fRunNumber(-1), fMappedRunNumber(0), fInCentralitySelection(-1), fFitModulation(0), fFitControl(0), fMinPvalue(0.01), fMaxPvalue(1), fNameSmallRho(""), fCachedRho(0), fSoftTrackMinPt(0.15), fSoftTrackMaxPt(5.), fSemiGoodJetMinPhi(0.), fSemiGoodJetMaxPhi(4.), fSemiGoodTrackMinPhi(0.), fSemiGoodTrackMaxPhi(4.), fAbsVertexZ(10), fHistCentrality(0), fHistVertexz(0), fHistRunnumbersPhi(0), fHistRunnumbersEta(0), fHistPvalueCDFROOT(0), fHistPvalueCDFROOTCent(0), fHistChi2ROOTCent(0), fHistPChi2Root(0), fHistPvalueCDF(0), fHistPvalueCDFCent(0), fHistChi2Cent(0), fHistPChi2(0), fHistKolmogorovTest(0), fHistKolmogorovTestCent(0), fHistPKolmogorov(0), fHistRhoStatusCent(0), fHistUndeterminedRunQA(0), fMinDisanceRCtoLJ(0), fMaxCones(-1), fExcludeLeadingJetsFromFit(1.), fRebinSwapHistoOnTheFly(kTRUE), fPercentageOfFits(10.), fOutputList(0), fOutputListGood(0), fOutputListBad(0), fHistAnalysisSummary(0), fHistSwap(0), fProfV2(0), fProfV2Cumulant(0), fProfV3(0), fProfV3Cumulant(0), fHistPsiControl(0), fHistPsiSpread(0), fHistPsiVZEROA(0), fHistPsiVZEROC(0), fHistPsiVZERO(0), fHistPsiTPC(0), fHistPsiVZEROAV0M(0), fHistPsiVZEROCV0M(0), fHistPsiVZEROVV0M(0), fHistPsiTPCV0M(0), fHistPsiVZEROATRK(0), fHistPsiVZEROCTRK(0), fHistPsiVZEROTRK(0), fHistPsiTPCTRK(0), fHistRhoVsMult(0), fHistRhoVsCent(0), fHistRhoAVsMult(0), fHistRhoAVsCent(0), fVZEROgainEqualization(0x0), fVZEROgainEqualizationPerRing(kFALSE), fChi2A(0x0), fChi2C(0x0), fChi3A(0x0), fChi3C(0x0), fOADB(0x0)
71 for(Int_t i(0); i < 10; i++) {
72 fProfV2Resolution[i] = 0;
73 fProfV3Resolution[i] = 0;
74 fHistPicoTrackPt[i] = 0;
75 fHistPicoTrackMult[i] = 0;
79 fHistClusterPt[i] = 0;
80 fHistClusterEtaPhi[i] = 0;
81 fHistClusterEtaPhiWeighted[i] = 0;
82 fHistPsiTPCLeadingJet[i] = 0;
83 fHistPsiVZEROALeadingJet[i] = 0;
84 fHistPsiVZEROCLeadingJet[i] = 0;
85 fHistPsiVZEROCombLeadingJet[i] = 0;
86 fHistPsi2Correlation[i] = 0;
87 fHistRhoPackage[i] = 0;
90 fHistRhoVsRCPt[i] = 0;
92 fHistDeltaPtDeltaPhi2[i] = 0;
93 fHistDeltaPtDeltaPhi2Rho0[i] = 0;
94 fHistRCPhiEtaExLJ[i] = 0;
95 fHistRhoVsRCPtExLJ[i] = 0;
97 fHistDeltaPtDeltaPhi2ExLJ[i] = 0;
98 fHistDeltaPtDeltaPhi2ExLJRho0[i] = 0;
101 fHistJetEtaPhi[i] = 0;
102 fHistJetPtArea[i] = 0;
103 fHistJetPtEta[i] = 0;
104 fHistJetPtConstituents[i] = 0;
105 fHistJetEtaRho[i] = 0;
106 fHistJetPsi2Pt[i] = 0;
107 fHistJetPsi2PtRho0[i] = 0;
109 for(Int_t i(0); i < 9; i++) {
110 for(Int_t j(0); j < 2; j++) {
111 for(Int_t k(0); k < 2; k++) {
112 fMeanQ[i][j][k] = 0.;
113 fWidthQ[i][j][k] = 0.;
114 fMeanQv3[i][j][k] = 0.;
115 fWidthQv3[i][j][k] = 0.;
119 for(Int_t i(0); i < 4; i++) {
123 for(Int_t i(0); i < 8; i++) fUseVZERORing[i] = kTRUE;
124 // default constructor
126 //_____________________________________________________________________________
127 AliAnalysisTaskJetV2::AliAnalysisTaskJetV2(const char* name, runModeType type) : AliAnalysisTaskEmcalJet(name, kTRUE),
128 fDebug(0), fRunToyMC(kFALSE), fLocalInit(0), fAttachToEvent(kTRUE), fFillHistograms(kTRUE), fFillQAHistograms(kTRUE), fReduceBinsXByFactor(-1.), fReduceBinsYByFactor(-1.), fNoEventWeightsForQC(kTRUE), fCentralityClasses(0), fExpectedRuns(0), fExpectedSemiGoodRuns(0), fUserSuppliedV2(0), fUserSuppliedV3(0), fUserSuppliedR2(0), fUserSuppliedR3(0), fTracksCont(0), fClusterCont(0), fJetsCont(0), fLeadingJet(0), fNAcceptedTracks(0), fNAcceptedTracksQCn(0), fFitModulationType(kNoFit), fFitGoodnessTest(kChi2Poisson), fQCRecovery(kTryFit), fUsePtWeight(kTRUE), fUsePtWeightErrorPropagation(kTRUE), fDetectorType(kVZEROComb), fAnalysisType(kCharged), fFitModulationOptions("QWLI"), fRunModeType(type), fDataType(kESD), fCollisionType(kPbPb), fRandom(0), fRunNumber(-1), fMappedRunNumber(0), fInCentralitySelection(-1), fFitModulation(0), fFitControl(0), fMinPvalue(0.01), fMaxPvalue(1), fNameSmallRho(""), fCachedRho(0), fSoftTrackMinPt(0.15), fSoftTrackMaxPt(5.), fSemiGoodJetMinPhi(0.), fSemiGoodJetMaxPhi(4.), fSemiGoodTrackMinPhi(0.), fSemiGoodTrackMaxPhi(4.), fAbsVertexZ(10), fHistCentrality(0), fHistVertexz(0), fHistRunnumbersPhi(0), fHistRunnumbersEta(0), fHistPvalueCDFROOT(0), fHistPvalueCDFROOTCent(0), fHistChi2ROOTCent(0), fHistPChi2Root(0), fHistPvalueCDF(0), fHistPvalueCDFCent(0), fHistChi2Cent(0), fHistPChi2(0), fHistKolmogorovTest(0), fHistKolmogorovTestCent(0), fHistPKolmogorov(0), fHistRhoStatusCent(0), fHistUndeterminedRunQA(0), fMinDisanceRCtoLJ(0), fMaxCones(-1), fExcludeLeadingJetsFromFit(1.), fRebinSwapHistoOnTheFly(kTRUE), fPercentageOfFits(10.), fOutputList(0), fOutputListGood(0), fOutputListBad(0), fHistAnalysisSummary(0), fHistSwap(0), fProfV2(0), fProfV2Cumulant(0), fProfV3(0), fProfV3Cumulant(0), fHistPsiControl(0), fHistPsiSpread(0), fHistPsiVZEROA(0), fHistPsiVZEROC(0), fHistPsiVZERO(0), fHistPsiTPC(0), fHistPsiVZEROAV0M(0), fHistPsiVZEROCV0M(0), fHistPsiVZEROVV0M(0), fHistPsiTPCV0M(0), fHistPsiVZEROATRK(0), fHistPsiVZEROCTRK(0), fHistPsiVZEROTRK(0), fHistPsiTPCTRK(0), fHistRhoVsMult(0), fHistRhoVsCent(0), fHistRhoAVsMult(0), fHistRhoAVsCent(0), fVZEROgainEqualization(0x0), fVZEROgainEqualizationPerRing(kFALSE), fChi2A(0x0), fChi2C(0x0), fChi3A(0x0), fChi3C(0x0), fOADB(0x0)
130 for(Int_t i(0); i < 10; i++) {
131 fProfV2Resolution[i] = 0;
132 fProfV3Resolution[i] = 0;
133 fHistPicoTrackPt[i] = 0;
134 fHistPicoTrackMult[i] = 0;
135 fHistPicoCat1[i] = 0;
136 fHistPicoCat2[i] = 0;
137 fHistPicoCat3[i] = 0;
138 fHistClusterPt[i] = 0;
139 fHistClusterEtaPhi[i] = 0;
140 fHistClusterEtaPhiWeighted[i] = 0;
141 fHistPsiTPCLeadingJet[i] = 0;
142 fHistPsiVZEROALeadingJet[i] = 0;
143 fHistPsiVZEROCLeadingJet[i] = 0;
144 fHistPsiVZEROCombLeadingJet[i] = 0;
145 fHistPsi2Correlation[i] = 0;
146 fHistRhoPackage[i] = 0;
148 fHistRCPhiEta[i] = 0;
149 fHistRhoVsRCPt[i] = 0;
151 fHistDeltaPtDeltaPhi2[i] = 0;
152 fHistDeltaPtDeltaPhi2Rho0[i] = 0;
153 fHistRCPhiEtaExLJ[i] = 0;
154 fHistRhoVsRCPtExLJ[i] = 0;
155 fHistRCPtExLJ[i] = 0;
156 fHistDeltaPtDeltaPhi2ExLJ[i] = 0;
157 fHistDeltaPtDeltaPhi2ExLJRho0[i] = 0;
158 fHistJetPtRaw[i] = 0;
160 fHistJetEtaPhi[i] = 0;
161 fHistJetPtArea[i] = 0;
162 fHistJetPtEta[i] = 0;
163 fHistJetPtConstituents[i] = 0;
164 fHistJetEtaRho[i] = 0;
165 fHistJetPsi2Pt[i] = 0;
166 fHistJetPsi2PtRho0[i] = 0;
168 for(Int_t i(0); i < 9; i++) {
169 for(Int_t j(0); j < 2; j++) {
170 for(Int_t k(0); k < 2; k++) {
171 fMeanQ[i][j][k] = 0.;
172 fWidthQ[i][j][k] = 0.;
173 fMeanQv3[i][j][k] = 0.;
174 fWidthQv3[i][j][k] = 0.;
178 for(Int_t i(0); i < 4; i++) {
182 for(Int_t i(0); i < 8; i++) fUseVZERORing[i] = kTRUE;
185 DefineInput(0, TChain::Class());
186 DefineOutput(1, TList::Class());
187 switch (fRunModeType) {
189 gStyle->SetOptFit(1);
190 DefineOutput(2, TList::Class());
191 DefineOutput(3, TList::Class());
193 default: fDebug = -1; // suppress debug info explicitely when not running locally
195 switch (fCollisionType) {
197 fFitModulationType = kNoFit;
201 if(fLocalRhoName=="") fLocalRhoName = Form("LocalRhoFrom_%s", GetName());
203 //_____________________________________________________________________________
204 AliAnalysisTaskJetV2::~AliAnalysisTaskJetV2()
207 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
208 if(fOutputList) {delete fOutputList; fOutputList = 0x0;}
209 if(fOutputListGood) {delete fOutputListGood; fOutputListGood = 0x0;}
210 if(fOutputListBad) {delete fOutputListBad; fOutputListBad = 0x0;}
211 if(fFitModulation) {delete fFitModulation; fFitModulation = 0x0;}
212 if(fHistSwap) {delete fHistSwap; fHistSwap = 0x0;}
213 if(fCentralityClasses) {delete fCentralityClasses; fCentralityClasses = 0x0;}
214 if(fExpectedRuns) {delete fExpectedRuns; fExpectedRuns = 0x0;}
215 if(fExpectedSemiGoodRuns) {delete fExpectedSemiGoodRuns; fExpectedSemiGoodRuns = 0x0;}
216 if(fFitControl) {delete fFitControl; fFitControl = 0x0;}
217 if(fVZEROgainEqualization) {delete fVZEROgainEqualization; fVZEROgainEqualization = 0x0;}
218 if(fChi2A) {delete fChi2A; fChi2A = 0x0;}
219 if(fChi2C) {delete fChi2C; fChi2C = 0x0;}
220 if(fChi3A) {delete fChi3A; fChi3A = 0x0;}
221 if(fChi3C) {delete fChi3C; fChi3C = 0x0;}
222 if(fOADB && !fOADB->IsZombie()) {
223 fOADB->Close(); fOADB = 0x0;
224 } else if (fOADB) fOADB = 0x0;
226 //_____________________________________________________________________________
227 void AliAnalysisTaskJetV2::ExecOnce()
230 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
231 fLocalRho = new AliLocalRhoParameter(fLocalRhoName.Data(), 0);
233 if(!(InputEvent()->FindListObject(fLocalRho->GetName()))) {
234 InputEvent()->AddObject(fLocalRho);
236 AliFatal(Form("%s: Container with name %s already present. Aborting", GetName(), fLocalRho->GetName()));
239 AliAnalysisTaskEmcalJet::ExecOnce(); // init the base class
240 AliAnalysisTaskEmcalJet::SetVzRange(-1.*fAbsVertexZ, fAbsVertexZ);
241 if(!GetJetContainer()) AliFatal(Form("%s: Couldn't find jet container. Aborting !", GetName()));
243 //_____________________________________________________________________________
244 Bool_t AliAnalysisTaskJetV2::Notify()
246 // determine the run number to see if the track and jet cuts should be refreshed for semi-good TPC runs
247 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
248 if(fRunNumber != InputEvent()->GetRunNumber()) {
249 fRunNumber = InputEvent()->GetRunNumber(); // set the current run number
250 if(fDebug > 0) printf("__FUNC__ %s > NEW RUNNUMBER DETECTED \n ", __func__);
251 // check if this is 10h or 11h data
252 switch (fCollisionType) {
254 if(fDebug > 0) printf(" LHC10h data, assuming full acceptance, reading VZERO calibration DB \n ");
255 // for 10h data the vzero event plane calibration needs to be cached
256 ReadVZEROCalibration2010h();
257 // no need to change rho or acceptance for 10h, so we're done
264 if(fDebug > 0) printf(" checking runnumber to adjust acceptance on the fly \n");
267 // reset the cuts. should be a pointless operation except for the case where the run number changes
268 // from semi-good back to good on one node, which is not a likely scenario (unless trains will
269 // run as one masterjob)
270 switch (fAnalysisType) {
272 AliAnalysisTaskEmcalJet::SetJetPhiLimits(-10., 10.);
275 AliAnalysisTaskEmcalJet::SetJetPhiLimits(1.405 + GetJetRadius(), 3.135 - GetJetRadius());
278 AliAnalysisTaskEmcal::SetTrackPhiLimits(-10., 10.);
281 if(fCachedRho) { // if there's a cached rho, it's the default, so switch back
282 if(fDebug > 0) printf("__FUNC__ %s > replacing rho with cached rho \n ", __func__);
283 fRho = fCachedRho; // reset rho back to cached value. again, should be pointless
285 Bool_t flaggedAsSemiGood(kFALSE); // not flagged as anything
286 for(Int_t i(0); i < fExpectedSemiGoodRuns->GetSize(); i++) {
287 if(fExpectedSemiGoodRuns->At(i) == fRunNumber) { // run is semi-good
288 if(fDebug > 0) printf("__FUNC__ %s > semi-good tpc run detected, adjusting acceptance \n ", __func__);
289 flaggedAsSemiGood = kTRUE;
290 switch (fAnalysisType) {
291 // for full jets the jet acceptance does not have to be changed as emcal does not
292 // cover the tpc low voltage readout strips
294 AliAnalysisTaskEmcalJet::SetJetPhiLimits(fSemiGoodJetMinPhi, fSemiGoodJetMaxPhi); // just an acceptance cut, jets are obtained from full azimuth, so no edge effects
298 AliAnalysisTaskEmcal::SetTrackPhiLimits(fSemiGoodTrackMinPhi, fSemiGoodTrackMaxPhi); // only affects vn extraction, NOT jet finding
299 // for semi-good runs, also try to get the 'small rho' estimate, if it is available
300 AliRhoParameter* tempRho(dynamic_cast<AliRhoParameter*>(InputEvent()->FindListObject(fNameSmallRho.Data())));
302 if(fDebug > 0) printf("__FUNC__ %s > switching to small rho, caching normal rho \n ", __func__);
303 fHistAnalysisSummary->SetBinContent(54, 1.); // bookkeep the fact that small rho is used
304 fCachedRho = fRho; // cache the original rho ...
305 fRho = tempRho; // ... and use the small rho
309 if(!flaggedAsSemiGood) {
310 // in case the run is not a semi-good run, check if it is recognized as another run
311 // only done to catch unexpected runs
312 for(Int_t i(0); i < fExpectedRuns->GetSize(); i++) {
313 if(fExpectedRuns->At(i) == fRunNumber) break; // run is known, break the loop else store the number in a random bin
314 fHistUndeterminedRunQA->SetBinContent(TMath::Nint(10.*gRandom->Uniform(0.,.9))+1, fRunNumber);
316 fHistAnalysisSummary->SetBinContent(53, 1.); // bookkeep which rho estimate is used
321 //_____________________________________________________________________________
322 Bool_t AliAnalysisTaskJetV2::InitializeAnalysis()
324 // initialize the anaysis
325 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
326 // if not set, estimate the number of cones that would fit into the selected acceptance
327 if(fMaxCones <= 0) fMaxCones = TMath::CeilNint((TMath::Abs(GetJetContainer()->GetJetEtaMax()-GetJetContainer()->GetJetEtaMin())*TMath::Abs(GetJetContainer()->GetJetPhiMax()-GetJetContainer()->GetJetPhiMin()))/(TMath::Pi()*GetJetRadius()*GetJetRadius()));
328 // manually 'override' the default acceptance cuts of the emcal framework (use with caution)
329 if(fMinDisanceRCtoLJ==0) fMinDisanceRCtoLJ = GetJetRadius();
330 if(dynamic_cast<AliAODEvent*>(InputEvent())) fDataType = kAOD; // determine the datatype
331 else if(dynamic_cast<AliESDEvent*>(InputEvent())) fDataType = kESD;
332 fHistAnalysisSummary->SetBinContent(36, (int)fDataType);
333 if(!fRandom) fRandom = new TRandom3(0); // set randomizer and random seed
334 switch (fFitModulationType) {
335 case kNoFit : { SetModulationFit(new TF1("fix_kNoFit", "[0]", 0, TMath::TwoPi())); } break;
337 SetModulationFit(new TF1("fit_kV2", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
338 fFitModulation->SetParameter(0, 0.); // normalization
339 fFitModulation->SetParameter(3, 0.2); // v2
340 fFitModulation->FixParameter(1, 1.); // constant
341 fFitModulation->FixParameter(2, 2.); // constant
344 SetModulationFit(new TF1("fit_kV3", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
345 fFitModulation->SetParameter(0, 0.); // normalization
346 fFitModulation->SetParameter(3, 0.2); // v3
347 fFitModulation->FixParameter(1, 1.); // constant
348 fFitModulation->FixParameter(2, 3.); // constant
350 default : { // for the combined fit, the 'direct fourier series' or the user supplied vn values we use v2 and v3
351 SetModulationFit(new TF1("fit_kCombined", "[0]*([1]+[2]*([3]*TMath::Cos([2]*(x-[4]))+[7]*TMath::Cos([5]*(x-[6]))))", 0, TMath::TwoPi()));
352 fFitModulation->SetParameter(0, 0.); // normalization
353 fFitModulation->SetParameter(3, 0.2); // v2
354 fFitModulation->FixParameter(1, 1.); // constant
355 fFitModulation->FixParameter(2, 2.); // constant
356 fFitModulation->FixParameter(5, 3.); // constant
357 fFitModulation->SetParameter(7, 0.2); // v3
360 switch (fRunModeType) {
361 case kGrid : { fFitModulationOptions += "N0"; } break;
364 FillAnalysisSummaryHistogram();
367 //_____________________________________________________________________________
368 TH1F* AliAnalysisTaskJetV2::BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Int_t c, Bool_t append)
370 // book a TH1F and connect it to the output container
371 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
372 if(fReduceBinsXByFactor > 0 ) bins = TMath::Nint(bins/fReduceBinsXByFactor);
373 if(!fOutputList) return 0x0;
375 if(c!=-1) { // format centrality dependent histograms accordingly
376 name = Form("%s_%i", name, c);
377 title += Form("_%i-%i", (int)(fCentralityClasses->At(c)), (int)(fCentralityClasses->At((1+c))));
379 title += Form(";%s;[counts]", x);
380 TH1F* histogram = new TH1F(name, title.Data(), bins, min, max);
382 if(append) fOutputList->Add(histogram);
385 //_____________________________________________________________________________
386 TH2F* AliAnalysisTaskJetV2::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, Bool_t append)
388 // book a TH2F and connect it to the output container
389 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
390 if(fReduceBinsXByFactor > 0 ) binsx = TMath::Nint(binsx/fReduceBinsXByFactor);
391 if(fReduceBinsYByFactor > 0 ) binsy = TMath::Nint(binsy/fReduceBinsYByFactor);
392 if(!fOutputList) return 0x0;
394 if(c!=-1) { // format centrality dependent histograms accordingly
395 name = Form("%s_%i", name, c);
396 title += Form("_%i-%i", (int)fCentralityClasses->At(c), (int)(fCentralityClasses->At((1+c))));
398 title += Form(";%s;%s", x, y);
399 TH2F* histogram = new TH2F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy);
401 if(append) fOutputList->Add(histogram);
404 //_____________________________________________________________________________
405 TH3F* AliAnalysisTaskJetV2::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, Bool_t append)
407 // book a TH2F and connect it to the output container
408 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
409 if(fReduceBinsXByFactor > 0 ) {
410 binsx = TMath::Nint(binsx/fReduceBinsXByFactor);
411 binsy = TMath::Nint(binsy/fReduceBinsXByFactor);
412 binsz = TMath::Nint(binsz/fReduceBinsXByFactor);
414 if(!fOutputList) return 0x0;
416 if(c!=-1) { // format centrality dependent histograms accordingly
417 name = Form("%s_%i", name, c);
418 title += Form("_%i-%i", (int)fCentralityClasses->At(c), (int)(fCentralityClasses->At((1+c))));
420 title += Form(";%s;%s;%s", x, y, z);
421 TH3F* histogram = new TH3F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy, binsz, minz, maxz);
423 if(append) fOutputList->Add(histogram);
426 //_____________________________________________________________________________
427 void AliAnalysisTaskJetV2::UserCreateOutputObjects()
429 // create output objects. also initializes some default values in case they aren't
430 // loaded via the AddTask macro
431 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
432 fOutputList = new TList();
433 fOutputList->SetOwner(kTRUE);
434 if(!fCentralityClasses) { // classes must be defined at this point
435 Double_t c[] = {0., 20., 40., 60., 80., 100.};
436 fCentralityClasses = new TArrayD(sizeof(c)/sizeof(c[0]), c);
438 if(!fExpectedRuns) { // expected runs must be defined at this point
439 Int_t r[] = {167813, 167988, 168066, 168068, 168069, 168076, 168104, 168212, 168311, 168322, 168325, 168341, 168361, 168362, 168458, 168460, 168461, 168992, 169091, 169094, 169138, 169143, 169167, 169417, 169835, 169837, 169838, 169846, 169855, 169858, 169859, 169923, 169956, 170027, 170036, 170081, /* up till here original good TPC list */169975, 169981, 170038, 170040, 170083, 170084, 170085, 170088, 170089, 170091, 170152, 170155, 170159, 170163, 170193, 170195, 170203, 170204, 170205, 170228, 170230, 170264, 170268, 170269, 170270, 170306, 170308, 170309, /* original semi-good tpc list */169415, 169411, 169035, 168988, 168984, 168826, 168777, 168512, 168511, 168467, 168464, 168342, 168310, 168115, 168108, 168107, 167987, 167915, 167903, /*new runs, good according to RCT */ 169238, 169160, 169156, 169148, 169145, 169144 /* run swith missing OROC 8 but seem ok in QA */};
440 fExpectedRuns = new TArrayI(sizeof(r)/sizeof(r[0]), r);
442 // set default semi-good runs only for 11h data
443 switch (fCollisionType) {
444 case kPbPb10h : break;
446 if(!fExpectedSemiGoodRuns) {
447 Int_t r[] = {169975, 169981, 170038, 170040, 170083, 170084, 170085, 170088, 170089, 170091, 170152, 170155, 170159, 170163, 170193, 170195, 170203, 170204, 170205, 170228, 170230, 170264, 170268, 170269, 170270, 170306, 170308, 170309};
448 fExpectedSemiGoodRuns = new TArrayI(sizeof(r)/sizeof(r[0]), r);
454 fHistCentrality = BookTH1F("fHistCentrality", "centrality", 102, -2, 100);
455 fHistVertexz = BookTH1F("fHistVertexz", "vertex z (cm)", 100, -12, 12);
457 // for some histograms adjust the bounds according to analysis acceptance
458 Double_t etaMin(-1.), etaMax(1.), phiMin(0.), phiMax(TMath::TwoPi());
459 switch (fAnalysisType) {
469 // pico track and emcal cluster kinematics
470 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
471 fHistPicoTrackPt[i] = BookTH1F("fHistPicoTrackPt", "p_{t} [GeV/c]", 100, 0, 100, i);
472 fHistPicoTrackMult[i] = BookTH1F("fHistPicoTrackMult", "multiplicity", 100, 0, 5000, i);
473 if(fFillQAHistograms) {
474 fHistPicoCat1[i] = BookTH2F("fHistPicoCat1", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
475 fHistPicoCat2[i] = BookTH2F("fHistPicoCat2", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
476 fHistPicoCat3[i] = BookTH2F("fHistPicoCat3", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
477 if(fAnalysisType == AliAnalysisTaskJetV2::kFull) {
478 fHistClusterPt[i] = BookTH1F("fHistClusterPt", "p_{t} [GeV/c]", 100, 0, 100, i);
479 fHistClusterEtaPhi[i] = BookTH2F("fHistClusterEtaPhi", "#eta", "#phi", 100, etaMax, etaMax, 100, phiMin, phiMax, i);
480 fHistClusterEtaPhiWeighted[i] = BookTH2F("fHistClusterEtaPhiWeighted", "#eta", "#phi", 100, etaMin, etaMax, 100, phiMin, phiMax, i);
482 fHistPsiTPCLeadingJet[i] = BookTH3F("fHistPsiTPCLeadingJet", "p_{t} [GeV/c]", "#Psi_{TPC}", "#varphi_{jet}", 70, -100, 250, 50, -1.*TMath::Pi()/2., TMath::Pi()/2., 50, phiMin, phiMax, i);
483 fHistPsiVZEROALeadingJet[i] = BookTH3F("fHistPsiVZEROALeadingJet", "p_{t} [GeV/c]", "#Psi_{VZEROA}", "#varphi_{jet}", 70, -100, 250, 50, -1.*TMath::Pi()/2., TMath::Pi()/2., 50, phiMin, phiMax, i);
484 fHistPsiVZEROCLeadingJet[i] = BookTH3F("fHistPsiVZEROCLeadingJet", "p_{t} [GeV/c]", "#Psi_{VZEROC}", "#varphi_{jet}", 70, -100, 250, 50, -1.*TMath::Pi()/2., TMath::Pi()/2., 50, phiMin, phiMax, i);
485 fHistPsiVZEROCombLeadingJet[i] = BookTH3F("fHistPsiVZEROCombLeadingJet", "p_{t} [GeV/c]", "#Psi_{VZEROComb}", "#varphi_{jet}", 70, -100, 250, 50, -1.*TMath::Pi()/2., TMath::Pi()/2., 50, phiMin, phiMax, i);
486 fHistPsi2Correlation[i] = BookTH3F("fHistPsi2Correlation", "#Psi_{TPC}", "#Psi_{VZEROA}", "#Psi_{VZEROC}", 20, -1.*TMath::Pi()/2., TMath::Pi()/2., 20, -1.*TMath::Pi()/2., TMath::Pi()/2., 20, -1.*TMath::Pi()/2., TMath::Pi()/2., i);
490 if(fFillQAHistograms) {
491 // event plane estimates and quality
492 fHistPsiControl = new TProfile("fHistPsiControl", "fHistPsiControl", 10, 0, 10);
493 fHistPsiControl->Sumw2();
494 fHistPsiSpread = new TProfile("fHistPsiSpread", "fHistPsiSpread", 4, 0, 4);
495 fHistPsiSpread->Sumw2();
496 fHistPsiControl->GetXaxis()->SetBinLabel(1, "<#Psi_{2, VZEROA}>");
497 fHistPsiControl->GetXaxis()->SetBinLabel(2, "<#Psi_{2, VZEROC}>");
498 fHistPsiControl->GetXaxis()->SetBinLabel(3, "<#Psi_{2, TPC}>");
499 fHistPsiControl->GetXaxis()->SetBinLabel(4, "<#Psi_{2, TPC, #eta < 0}>");
500 fHistPsiControl->GetXaxis()->SetBinLabel(5, "<#Psi_{2, TPC, #eta > 0}>");
501 fHistPsiControl->GetXaxis()->SetBinLabel(6, "<#Psi_{3, VZEROA}>");
502 fHistPsiControl->GetXaxis()->SetBinLabel(7, "<#Psi_{3, VZEROC}>");
503 fHistPsiControl->GetXaxis()->SetBinLabel(8, "<#Psi_{3, TPC}>");
504 fHistPsiControl->GetXaxis()->SetBinLabel(9, "<#Psi_{3, TPC, #eta < 0}>");
505 fHistPsiControl->GetXaxis()->SetBinLabel(10, "<#Psi_{3, TPC, #eta > 0}>");
506 fHistPsiSpread->GetXaxis()->SetBinLabel(1, "<#Psi_{2, VZEROA} - #Psi_{2, VZEROC}>");
507 fHistPsiSpread->GetXaxis()->SetBinLabel(2, "<#Psi_{2, VZEROC} - #Psi_{2, TPC}>");
508 fHistPsiSpread->GetXaxis()->SetBinLabel(3, "<#Psi_{2, VZEROC} - #Psi_{2, TPC}>");
509 fHistPsiSpread->GetXaxis()->SetBinLabel(4, "<#Psi_{2, TPC, #eta < 0} - #Psi_{2, TPC, #eta > 0}>");
510 fOutputList->Add(fHistPsiControl);
511 fOutputList->Add(fHistPsiSpread);
512 fHistPsiVZEROA = BookTH1F("fHistPsiVZEROA", "#Psi_{VZEROA}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
513 fHistPsiVZEROC = BookTH1F("fHistPsiVZEROC", "#Psi_{VZEROC}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
514 fHistPsiVZERO = BookTH1F("fHistPsiVZERO", "#Psi_{VZERO}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
515 fHistPsiTPC = BookTH1F("fHistPsiTPC", "#Psi_{TPC}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
516 fHistPsiVZEROAV0M = BookTH2F("fHistPsiVZEROAV0M", "V0M", "#Psi_{2, VZEROA}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
517 fHistPsiVZEROCV0M = BookTH2F("fHistPsiVZEROCV0M", "V0M", "#Psi_{2, VZEROC}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
518 fHistPsiVZEROVV0M = BookTH2F("fHistPsiVZEROV0M", "V0M", "#Psi_{2, VZERO}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
519 fHistPsiTPCV0M = BookTH2F("fHistPsiTPCV0M", "V0M", "#Psi_{2, TRK}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
520 fHistPsiVZEROATRK = BookTH2F("fHistPsiVZEROATRK", "TRK", "#Psi_{2, VZEROA}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
521 fHistPsiVZEROCTRK = BookTH2F("fHistPsiVZEROCTRK", "TRK", "#Psi_{2, VZEROC}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
522 fHistPsiVZEROTRK = BookTH2F("fHistPsiVZEROTRK", "TRK", "#Psi_{2, VZERO}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
523 fHistPsiTPCTRK = BookTH2F("fHistPsiTPCTRK", "TRK", "#Psi_{2, TRK}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
526 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i ++) {
527 fHistRhoPackage[i] = BookTH1F("fHistRhoPackage", "#rho [GeV/c]", 100, 0, 150, i);
528 fHistRho[i] = BookTH1F("fHistRho", "#rho [GeV/c]", 100, 0, 150, i);
530 fHistRhoVsMult = BookTH2F("fHistRhoVsMult", "multiplicity", "#rho [GeV/c]", 100, 0, 4000, 100, 0, 250);
531 fHistRhoVsCent = BookTH2F("fHistRhoVsCent", "centrality", "#rho [GeV/c]", 100, 0, 100, 100, 0, 250);
532 fHistRhoAVsMult = BookTH2F("fHistRhoAVsMult", "multiplicity", "#rho * A (jet) [GeV/c]", 100, 0, 4000, 100, 0, 50);
533 fHistRhoAVsCent = BookTH2F("fHistRhoAVsCent", "centrality", "#rho * A (jet) [GeV/c]", 100, 0, 100, 100, 0, 50);
535 TString detector("");
536 switch (fDetectorType) {
537 case kTPC : detector+="TPC";
539 case kVZEROA : detector+="VZEROA";
541 case kVZEROC : detector+="VZEROC";
543 case kVZEROComb : detector+="VZEROComb";
545 case kFixedEP : detector+="FixedEP";
549 // delta pt distributions
550 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i ++) {
551 if(fFillQAHistograms) fHistRCPhiEta[i] = BookTH2F("fHistRCPhiEta", "#phi (RC)", "#eta (RC)", 40, phiMin, phiMax, 40, etaMin, etaMax, i);
552 fHistRhoVsRCPt[i] = BookTH2F("fHistRhoVsRCPt", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i);
553 fHistRCPt[i] = BookTH1F("fHistRCPt", "p_{t} (RC) [GeV/c]", 130, -20, 150, i);
554 if(fFillQAHistograms) fHistRCPhiEtaExLJ[i] = BookTH2F("fHistRCPhiEtaExLJ", "#phi (RC)", "#eta (RC)", 40, phiMin, phiMax, 40, etaMin, etaMax, i);
555 fHistDeltaPtDeltaPhi2[i] = BookTH2F("fHistDeltaPtDeltaPhi2", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
556 fHistDeltaPtDeltaPhi2Rho0[i] = BookTH2F("fHistDeltaPtDeltaPhi2Rho0", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
557 fHistRhoVsRCPtExLJ[i] = BookTH2F("fHistRhoVsRCPtExLJ", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i);
558 fHistRCPtExLJ[i] = BookTH1F("fHistRCPtExLJ", "p_{t} (RC) [GeV/c]", 130, -20, 150, i);
559 fHistDeltaPtDeltaPhi2ExLJ[i] = BookTH2F("fHistDeltaPtDeltaPhi2ExLJ", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
560 fHistDeltaPtDeltaPhi2ExLJRho0[i] = BookTH2F("fHistDeltaPtDeltaPhi2ExLJRho0", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
561 // jet histograms (after kinematic cuts)
562 fHistJetPtRaw[i] = BookTH1F("fHistJetPtRaw", "p_{t, jet} RAW [GeV/c]", 200, -50, 150, i);
563 fHistJetPt[i] = BookTH1F("fHistJetPt", "p_{t, jet} [GeV/c]", 350, -100, 250, i);
564 if(fFillQAHistograms) fHistJetEtaPhi[i] = BookTH2F("fHistJetEtaPhi", "#eta", "#phi", 100, etaMin, etaMax, 100, phiMin, phiMax, i);
565 fHistJetPtArea[i] = BookTH2F("fHistJetPtArea", "p_{t, jet} [GeV/c]", "Area", 175, -100, 250, 30, 0, 0.9, i);
566 fHistJetPtEta[i] = BookTH2F("fHistJetPtEta", "p_{t, jet} [GeV/c]", "Eta", 175, -100, 250, 30, etaMin, etaMax, i);
567 fHistJetPtConstituents[i] = BookTH2F("fHistJetPtConstituents", "p_{t, jet} [GeV/c]", "no. of constituents", 350, -100, 250, 60, 0, 150, i);
568 fHistJetEtaRho[i] = BookTH2F("fHistJetEtaRho", "#eta", "#rho", 100, etaMin, etaMax, 100, 0, 300, i);
569 // in plane and out of plane spectra
570 fHistJetPsi2Pt[i] = BookTH2F("fHistJetPsi2Pt", Form("#phi_{jet} - #Psi_{2, %s}", detector.Data()), "p_{t, jet} [GeV/c]", 40, 0., TMath::Pi(), 350, -100, 250, i);
571 fHistJetPsi2PtRho0[i] = BookTH2F("fHistJetPsi2PtRho0", Form("#phi_{jet} - #Psi_{2, %s}", detector.Data()), "p_{t, jet} [GeV/c]", 40, 0., TMath::Pi(), 350, -100, 250, i);
572 // profiles for all correlator permutations which are necessary to calculate each second and third order event plane resolution
573 fProfV2Resolution[i] = new TProfile(Form("fProfV2Resolution_%i", i), Form("fProfV2Resolution_%i", i), 11, -0.5, 10.5);
574 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(3, "<cos(2(#Psi_{VZEROA} - #Psi_{VZEROC}))>");
575 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(4, "<cos(2(#Psi_{VZEROC} - #Psi_{VZEROA}))>");
576 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(5, "<cos(2(#Psi_{VZEROA} - #Psi_{TPC}))>");
577 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(6, "<cos(2(#Psi_{TPC} - #Psi_{VZEROA}))>");
578 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(7, "<cos(2(#Psi_{VZEROC} - #Psi_{TPC}))>");
579 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(8, "<cos(2(#Psi_{TPC} - #Psi_{VZEROC}))>");
580 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(9, "<cos(2(#Psi_{VZERO} - #Psi_{TPC_A}))>");
581 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(10, "<cos(2(#Psi_{VZERO} - #Psi_{TPC_B}))>");
582 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(11, "<cos(2(#Psi_{TPC_A} - #Psi_{TPC_B}))>");
583 fOutputList->Add(fProfV2Resolution[i]);
584 fProfV3Resolution[i] = new TProfile(Form("fProfV3Resolution_%i", i), Form("fProfV3Resolution_%i", i), 11, -0.5, 10.5);
585 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(3, "<cos(3(#Psi_{VZEROA} - #Psi_{VZEROC}))>");
586 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(4, "<cos(3(#Psi_{VZEROC} - #Psi_{VZEROA}))>");
587 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(5, "<cos(3(#Psi_{VZEROA} - #Psi_{TPC}))>");
588 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(6, "<cos(3(#Psi_{TPC} - #Psi_{VZEROA}))>");
589 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(7, "<cos(3(#Psi_{VZEROC} - #Psi_{TPC}))>");
590 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(8, "<cos(3(#Psi_{TPC} - #Psi_{VZEROC}))>");
591 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(9, "<cos(3(#Psi_{VZERO} - #Psi_{TPC_A}))>");
592 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(10, "<cos(3(#Psi_{VZERO} - #Psi_{TPC_B}))>");
593 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(11, "<cos(3(#Psi_{TPC_A} - #Psi_{TPC_B}))>");
594 fOutputList->Add(fProfV3Resolution[i]);
597 Float_t temp[fCentralityClasses->GetSize()];
598 for(Int_t i(0); i < fCentralityClasses->GetSize(); i++) temp[i] = fCentralityClasses->At(i);
599 fProfV2 = new TProfile("fProfV2", "fProfV2", fCentralityClasses->GetSize()-1, temp);
600 fProfV3 = new TProfile("fProfV3", "fProfV3", fCentralityClasses->GetSize()-1, temp);
601 fOutputList->Add(fProfV2);
602 fOutputList->Add(fProfV3);
603 switch (fFitModulationType) {
605 fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
606 fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
607 fOutputList->Add(fProfV2Cumulant);
608 fOutputList->Add(fProfV3Cumulant);
611 fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
612 fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
613 fOutputList->Add(fProfV2Cumulant);
614 fOutputList->Add(fProfV3Cumulant);
618 // for the histograms initialized below, binning is fixed to runnumbers or flags
619 fReduceBinsXByFactor = 1;
620 fReduceBinsYByFactor = 1;
621 if(fFillQAHistograms) {
622 fHistRunnumbersEta = new TH2F("fHistRunnumbersEta", "fHistRunnumbersEta", fExpectedRuns->GetSize()+1, -.5, fExpectedRuns->GetSize()+.5, 100, -1.1, 1.1);
623 fHistRunnumbersEta->Sumw2();
624 fOutputList->Add(fHistRunnumbersEta);
625 fHistRunnumbersPhi = new TH2F("fHistRunnumbersPhi", "fHistRunnumbersPhi", fExpectedRuns->GetSize()+1, -.5, fExpectedRuns->GetSize()+.5, 100, -0.2, TMath::TwoPi()+0.2);
626 fHistRunnumbersPhi->Sumw2();
627 fOutputList->Add(fHistRunnumbersPhi);
628 for(Int_t i(0); i < fExpectedRuns->GetSize(); i++) {
629 fHistRunnumbersPhi->GetXaxis()->SetBinLabel(i+1, Form("%i", fExpectedRuns->At(i)));
630 fHistRunnumbersEta->GetXaxis()->SetBinLabel(i+1, Form("%i", fExpectedRuns->At(i)));
632 fHistRunnumbersPhi->GetXaxis()->SetBinLabel(fExpectedRuns->GetSize()+1, "undetermined");
633 fHistRunnumbersEta->GetXaxis()->SetBinLabel(fExpectedRuns->GetSize()+1, "undetermined");
635 fHistAnalysisSummary = BookTH1F("fHistAnalysisSummary", "flag", 54, -0.5, 54.5);
636 fHistSwap = new TH1F("fHistSwap", "fHistSwap", 20, 0, TMath::TwoPi());
637 if(fUsePtWeight) fHistSwap->Sumw2();
639 if(fUserSuppliedV2) fOutputList->Add(fUserSuppliedV2);
640 if(fUserSuppliedV3) fOutputList->Add(fUserSuppliedV3);
641 if(fUserSuppliedR2) fOutputList->Add(fUserSuppliedR2);
642 if(fUserSuppliedR3) fOutputList->Add(fUserSuppliedR3);
643 // increase readability of output list
645 // cdf and pdf of chisquare distribution
646 fHistPvalueCDF = BookTH1F("fHistPvalueCDF", "CDF #chi^{2}", 50, 0, 1);
647 fHistPvalueCDFCent = BookTH2F("fHistPvalueCDFCent", "centrality", "p-value", 40, 0, 100, 40, 0, 1);
648 fHistChi2Cent = BookTH2F("fHistChi2Cent", "centrality", "#tilde{#chi^{2}}", 100, 0, 100, 100, 0, 5);
649 fHistPChi2 = BookTH2F("fHistPChi2", "p-value", "#tilde{#chi^{2}}", 1000, 0, 1, 100, 0, 5);
650 fHistKolmogorovTest = BookTH1F("fHistKolmogorovTest", "KolmogorovTest", 50, 0, 1);
651 fHistKolmogorovTestCent = BookTH2F("fHistKolmogorovTestCent", "centrality", "Kolmogorov p", 40, 0, 100, 45, 0, 1);
652 fHistPvalueCDFROOT = BookTH1F("fHistPvalueCDFROOT", "CDF #chi^{2} ROOT", 50, 0, 1);
653 fHistPvalueCDFROOTCent = BookTH2F("fHistPvalueCDFROOTCent", "centrality", "p-value ROOT", 40, 0, 100, 45, 0, 1);
654 fHistChi2ROOTCent = BookTH2F("fHistChi2ROOTCent", "centrality", "#tilde{#chi^{2}}", 40, 0, 100, 45, 0, 5);
655 fHistPChi2Root = BookTH2F("fHistPChi2Root", "p-value", "#tilde{#chi^{2}} ROOT", 1000, 0, 1, 100, 0, 5);
656 fHistPKolmogorov = BookTH2F("fHistPKolmogorov", "p-value", "kolmogorov p",40, 0, 1, 40, 0, 1);
657 fHistRhoStatusCent = BookTH2F("fHistRhoStatusCent", "centrality", "status [-1=lin was better, 0=ok, 1 = failed]", 101, -1, 100, 3, -1.5, 1.5);
658 fHistUndeterminedRunQA = BookTH1F("fHistUndeterminedRunQA", "runnumber", 10, 0, 10);
660 PostData(1, fOutputList);
662 switch (fRunModeType) {
664 fOutputListGood = new TList();
665 fOutputListGood->SetOwner(kTRUE);
666 fOutputListBad = new TList();
667 fOutputListBad->SetOwner(kTRUE);
668 PostData(2, fOutputListGood);
669 PostData(3, fOutputListBad);
674 // get the containers
675 fTracksCont = GetParticleContainer("Tracks");
676 fClusterCont = GetClusterContainer(0); // get the default cluster container
677 fJetsCont = GetJetContainer("Jets");
679 //_____________________________________________________________________________
680 Bool_t AliAnalysisTaskJetV2::Run()
682 // called for each accepted event (call made from user exec of parent class)
683 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
684 if(!fTracks||!fJets||!fRho) {
685 if(!fTracks) printf(" > Failed to retrieve fTracks ! < \n");
686 if(!fJets) printf(" > Failed to retrieve fJets ! < \n");
687 if(!fRho) printf(" > Failed to retrieve fRho ! < \n");
690 if(!fLocalInit) fLocalInit = InitializeAnalysis();
691 // reject the event if expected data is missing
692 if(!PassesCuts(InputEvent())) return kFALSE;
693 // cache the leading jet within acceptance
694 fLeadingJet = GetLeadingJet();
696 fLocalRho->SetVal(fRho->GetVal());
697 // place holder arrays for the event planes
699 // [0][0] psi2a [1,0] psi2c
700 // [0][1] psi3a [1,1] psi3c
701 Double_t vzero[2][2];
702 /* for the combined vzero event plane
704 * not fully implmemented yet, use with caution ! */
705 Double_t vzeroComb[2];
708 // evaluate the actual event planes
709 switch (fDetectorType) {
711 // for fixed, fix all ep's to default values
712 tpc[0] = 0.; tpc[1] = 1.;
713 vzero[0][0] = 0.; vzero[0][1] = 1.;
714 vzero[1][0] = 0.; vzero[1][1] = 1.;
715 vzeroComb[0] = 0.; vzeroComb[1] = 1.;
718 // else grab the actual data
719 CalculateEventPlaneVZERO(vzero);
720 CalculateEventPlaneCombinedVZERO(vzeroComb);
721 CalculateEventPlaneTPC(tpc);
724 Double_t psi2(-1), psi3(-1);
725 // arrays which will hold the fit parameters
726 switch (fDetectorType) { // determine the detector type for the rho fit
727 case kTPC : { psi2 = tpc[0]; psi3 = tpc[1]; } break;
728 case kVZEROA : { psi2 = vzero[0][0]; psi3 = vzero[0][1]; } break;
729 case kVZEROC : { psi2 = vzero[1][0]; psi3 = vzero[1][1]; } break;
730 case kVZEROComb : { psi2 = vzeroComb[0]; psi3 = vzeroComb[1];} break;
731 case kFixedEP : { psi2 = 0.; psi3 = 1.;} break;
734 switch (fFitModulationType) { // do the fits
736 switch (fCollisionType) {
737 case kPythia : { // background is zero for pp jets
738 fFitModulation->FixParameter(0, 0);
739 fLocalRho->SetVal(0);
742 fFitModulation->FixParameter(0, fLocalRho->GetVal());
746 case kV2 : { // only v2
747 if(CorrectRho(psi2, psi3)) {
748 fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
749 if(fUserSuppliedR2) {
750 Double_t r(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
751 if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
753 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
756 case kV3 : { // only v3
757 if(CorrectRho(psi2, psi3)) {
758 if(fUserSuppliedR3) {
759 Double_t r(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
760 if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
762 fProfV3->Fill(fCent, fFitModulation->GetParameter(3));
763 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
766 case kQC2 : { // qc2 analysis
767 if(CorrectRho(psi2, psi3)) {
768 if(fUserSuppliedR2 && fUserSuppliedR3) {
769 // note for the qc method, resolution is REVERSED to go back to v2obs
770 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
771 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
772 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
773 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)*r3);
775 if (fUsePtWeight) { // use weighted weights
776 Double_t dQCnM11 = (fNoEventWeightsForQC) ? 1. : QCnM11();
777 fProfV2->Fill(fCent, fFitModulation->GetParameter(3), dQCnM11);
778 fProfV3->Fill(fCent, fFitModulation->GetParameter(7), dQCnM11);
780 Double_t dQCnM = (fNoEventWeightsForQC) ? 2. : QCnM();
781 fProfV2->Fill(fCent, fFitModulation->GetParameter(3), dQCnM*(dQCnM-1));
782 fProfV3->Fill(fCent, fFitModulation->GetParameter(7), dQCnM*(dQCnM-1));
784 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
788 if(CorrectRho(psi2, psi3)) {
789 if(fUserSuppliedR2 && fUserSuppliedR3) {
790 // note for the qc method, resolution is REVERSED to go back to v2obs
791 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
792 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
793 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
794 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)*r3);
796 if (fUsePtWeight) { // use weighted weights
797 fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM1111()*/);
798 fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM1111()*/);
800 fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
801 fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
804 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
807 if(CorrectRho(psi2, psi3)) {
808 if(fUserSuppliedR2 && fUserSuppliedR3) {
809 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
810 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
811 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r2);
812 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)/r3);
814 fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
815 fProfV3->Fill(fCent, fFitModulation->GetParameter(7));
816 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
820 // if all went well, update the local rho parameter
821 fLocalRho->SetLocalRho(fFitModulation);
822 // fill a number of histograms. event qa needs to be filled first as it also determines the runnumber for the track qa
823 if(fFillQAHistograms) FillQAHistograms(InputEvent());
824 if(fFillHistograms) FillHistogramsAfterSubtraction(psi2, vzero, vzeroComb, tpc);
825 // send the output to the connected output container
826 PostData(1, fOutputList);
827 switch (fRunModeType) {
829 PostData(2, fOutputListGood);
830 PostData(3, fOutputListBad);
836 //_____________________________________________________________________________
837 void AliAnalysisTaskJetV2::Exec(Option_t* c)
839 // for stand alone, avoid framework event setup
840 switch (fCollisionType) {
842 // need to call ExecOnce as it is not loaded otherwise
843 if(!fLocalRho) AliAnalysisTaskJetV2::ExecOnce();
844 AliAnalysisTaskJetV2::Run();
847 AliAnalysisTaskSE::Exec(c);
851 //_____________________________________________________________________________
852 Double_t AliAnalysisTaskJetV2::CalculateEventPlaneChi(Double_t res)
854 // return chi for given resolution to combine event plane estimates from two subevents
855 // see Phys. Rev. C no. CS6346 (http://arxiv.org/abs/nucl-ex/9805001)
856 Double_t chi(2.), delta(1.), con((TMath::Sqrt(TMath::Pi()))/(2.*TMath::Sqrt(2)));
857 for (Int_t i(0); i < 15; i++) {
858 chi = ((con*chi*TMath::Exp(-chi*chi/4.)*(TMath::BesselI0(chi*chi/4.)+TMath::BesselI1(chi*chi/4.))) < res) ? chi + delta : chi - delta;
863 //_____________________________________________________________________________
864 void AliAnalysisTaskJetV2::CalculateEventPlaneVZERO(Double_t vzero[2][2]) const
866 // get the vzero event plane (a and c separately)
867 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
868 switch (fCollisionType) {
870 // for 10h data, get the calibrated q-vector from the database
871 Double_t QA2[] = {-999., -999.};
872 Double_t QA3[] = {-999., -999.};
873 Double_t QC2[] = {-999., -999.};
874 Double_t QC3[] = {-999., -999.};
875 CalculateQvectorVZERO(QA2, QA3, QC2, QC3);
876 vzero[0][0] = .5*TMath::ATan2(QA2[1], QA2[0]);
877 vzero[1][0] = .5*TMath::ATan2(QC2[1], QC2[0]);
878 vzero[0][1] = (1./3.)*TMath::ATan2(QA3[1], QA3[0]);
879 vzero[1][1] = (1./3.)*TMath::ATan2(QC3[1], QC3[0]);
882 // by default use the ep from the event header (make sure EP selection task is enabeled!)
883 Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0), h(0);
884 vzero[0][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 2, a, b);
885 vzero[1][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 2, c, d);
886 vzero[0][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 3, e, f);
887 vzero[1][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 3, g, h);
892 //_____________________________________________________________________________
893 void AliAnalysisTaskJetV2::CalculateEventPlaneCombinedVZERO(Double_t* comb) const
895 // return the combined vzero event plane
896 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
897 switch (fCollisionType) {
898 // for 10h data call calibration info
900 // get the calibrated q-vectors
901 Double_t Q2[] = {-999., -999.};
902 Double_t Q3[] = {-999., -999.};
903 // return if something isn't ok from the calibration side
904 CalculateQvectorCombinedVZERO(Q2, Q3);
905 comb[0] = .5*TMath::ATan2(Q2[1], Q2[0]);
906 comb[1] = (1./3.)*TMath::ATan2(Q3[1], Q3[0]);
909 // for all other types use calibrated event plane from the event header
910 Double_t a(0), b(0), c(0), d(0);
911 comb[0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 10, 2, a, b);
912 comb[1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 10, 3, c, d);
916 //_____________________________________________________________________________
917 void AliAnalysisTaskJetV2::CalculateEventPlaneTPC(Double_t* tpc)
919 // grab the TPC event plane
920 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
921 fNAcceptedTracks = 0; // reset the track counter
922 Double_t qx2(0), qy2(0); // for psi2
923 Double_t qx3(0), qy3(0); // for psi3
925 Float_t excludeInEta = -999;
926 if(fExcludeLeadingJetsFromFit > 0 ) { // remove the leading jet from ep estimate
927 if(fLeadingJet) excludeInEta = fLeadingJet->Eta();
929 for(Int_t iTPC(0); iTPC < fTracksCont->GetNEntries(); iTPC++) {
930 AliVParticle* track = fTracksCont->GetParticle(iTPC);
931 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
932 if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta) < GetJetContainer()->GetJetRadius()*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - GetJetContainer()->GetJetRadius() - GetJetContainer()->GetJetEtaMax() ) > 0 )) continue;
934 qx2+= TMath::Cos(2.*track->Phi());
935 qy2+= TMath::Sin(2.*track->Phi());
936 qx3+= TMath::Cos(3.*track->Phi());
937 qy3+= TMath::Sin(3.*track->Phi());
940 tpc[0] = .5*TMath::ATan2(qy2, qx2);
941 tpc[1] = (1./3.)*TMath::ATan2(qy3, qx3);
943 //_____________________________________________________________________________
944 void AliAnalysisTaskJetV2::CalculateEventPlaneResolution(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc)
946 // fill the profiles for the resolution parameters
947 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
948 fProfV2Resolution[fInCentralitySelection]->Fill(2., TMath::Cos(2.*(vzero[0][0] - vzero[1][0])));
949 fProfV2Resolution[fInCentralitySelection]->Fill(3., TMath::Cos(2.*(vzero[1][0] - vzero[0][0])));
950 fProfV2Resolution[fInCentralitySelection]->Fill(4., TMath::Cos(2.*(vzero[0][0] - tpc[0])));
951 fProfV2Resolution[fInCentralitySelection]->Fill(5., TMath::Cos(2.*(tpc[0] - vzero[0][0])));
952 fProfV2Resolution[fInCentralitySelection]->Fill(6., TMath::Cos(2.*(vzero[1][0] - tpc[0])));
953 fProfV2Resolution[fInCentralitySelection]->Fill(7., TMath::Cos(2.*(tpc[0] - vzero[1][0])));
954 fProfV3Resolution[fInCentralitySelection]->Fill(2., TMath::Cos(3.*(vzero[0][0] - vzero[1][0])));
955 fProfV3Resolution[fInCentralitySelection]->Fill(3., TMath::Cos(3.*(vzero[1][0] - vzero[0][0])));
956 fProfV3Resolution[fInCentralitySelection]->Fill(4., TMath::Cos(3.*(vzero[0][0] - tpc[0])));
957 fProfV3Resolution[fInCentralitySelection]->Fill(5., TMath::Cos(3.*(tpc[0] - vzero[0][0])));
958 fProfV3Resolution[fInCentralitySelection]->Fill(6., TMath::Cos(3.*(vzero[1][0] - tpc[0])));
959 fProfV3Resolution[fInCentralitySelection]->Fill(7., TMath::Cos(3.*(tpc[0] - vzero[1][0])));
960 // for the resolution of the combined vzero event plane, use two tpc halves as uncorrelated subdetectors
961 Double_t qx2a(0), qy2a(0); // for psi2a, negative eta
962 Double_t qx3a(0), qy3a(0); // for psi3a, negative eta
963 Double_t qx2b(0), qy2b(0); // for psi2a, positive eta
964 Double_t qx3b(0), qy3b(0); // for psi3a, positive eta
966 Int_t iTracks(fTracks->GetEntriesFast());
967 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
968 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
969 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
970 if(track->Eta() < 0 ) {
971 qx2a+= TMath::Cos(2.*track->Phi());
972 qy2a+= TMath::Sin(2.*track->Phi());
973 qx3a+= TMath::Cos(3.*track->Phi());
974 qy3a+= TMath::Sin(3.*track->Phi());
975 } else if (track->Eta() > 0) {
976 qx2b+= TMath::Cos(2.*track->Phi());
977 qy2b+= TMath::Sin(2.*track->Phi());
978 qx3b+= TMath::Cos(3.*track->Phi());
979 qy3b+= TMath::Sin(3.*track->Phi());
983 Double_t tpca2(.5*TMath::ATan2(qy2a, qx2a));
984 Double_t tpca3((1./3.)*TMath::ATan2(qy3a, qx3a));
985 Double_t tpcb2(.5*TMath::ATan2(qy2b, qx2b));
986 Double_t tpcb3((1./3.)*TMath::ATan2(qy3b, qx3b));
987 fProfV2Resolution[fInCentralitySelection]->Fill(8., TMath::Cos(2.*(vzeroComb[0] - tpca2)));
988 fProfV2Resolution[fInCentralitySelection]->Fill(9., TMath::Cos(2.*(vzeroComb[0] - tpcb2)));
989 fProfV2Resolution[fInCentralitySelection]->Fill(10., TMath::Cos(2.*(tpca2 - tpcb2)));
990 fProfV3Resolution[fInCentralitySelection]->Fill(8., TMath::Cos(3.*(vzeroComb[1] - tpca3)));
991 fProfV3Resolution[fInCentralitySelection]->Fill(9., TMath::Cos(3.*(vzeroComb[1] - tpcb3)));
992 fProfV3Resolution[fInCentralitySelection]->Fill(10., TMath::Cos(3.*(tpca3 - tpcb3)));
994 //_____________________________________________________________________________
995 void AliAnalysisTaskJetV2::CalculateQvectorVZERO(Double_t Qa2[2], Double_t Qc2[2], Double_t Qa3[2], Double_t Qc3[2]) const
997 // return the calibrated 2nd and 3rd order q-vectors for vzeroa and vzeroc
998 // function takes arrays as arguments, which correspond to vzero info in the following way
1000 // Qa2[0] = Qx2 for vzero A Qa2[1] = Qy2 for vzero A (etc)
1002 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1003 // placeholders for geometric information
1004 Double_t phi(-999.), weight(-999.);
1005 // reset placeholders for Q-vector components
1006 Qa2[0] = 0.; Qc2[0] = 0.; Qa3[0] = 0.; Qc3[0] = 0.;
1007 Qa2[1] = 0.; Qc2[1] = 0.; Qa3[1] = 0.; Qc3[1] = 0.;
1009 for(Int_t i(0); i < 64; i++) {
1010 // loop over all scintillators, construct Q-vectors in the same loop
1011 phi = TMath::PiOver4()*(0.5+i%8);
1013 // note that disabled rings have already been excluded in ReadVZEROCalibration2010h
1014 if(i<32) { // v0c side
1015 if(i < 8) weight = InputEvent()->GetVZEROData()->GetMultiplicity(i)*fVZEROCpol[0]/fVZEROgainEqualization->GetBinContent(1+i);
1016 else if (i < 16 ) weight = InputEvent()->GetVZEROData()->GetMultiplicity(i)*fVZEROCpol[1]/fVZEROgainEqualization->GetBinContent(1+i);
1017 else if (i < 24 ) weight = InputEvent()->GetVZEROData()->GetMultiplicity(i)*fVZEROCpol[2]/fVZEROgainEqualization->GetBinContent(1+i);
1018 else if (i < 32 ) weight = InputEvent()->GetVZEROData()->GetMultiplicity(i)*fVZEROCpol[3]/fVZEROgainEqualization->GetBinContent(1+i);
1019 // fill Q-vectors for v0c side
1020 Qc2[0]+=weight*TMath::Cos(2.*phi);
1021 Qc3[0]+=weight*TMath::Cos(3.*phi);
1022 Qc2[1]+=weight*TMath::Sin(2.*phi);
1023 Qc3[1]+=weight*TMath::Sin(3.*phi);
1024 } else { // v0a side
1025 if( i < 40) weight = InputEvent()->GetVZEROData()->GetMultiplicity(i)*fVZEROApol[0]/fVZEROgainEqualization->GetBinContent(1+i);
1026 else if ( i < 48 ) weight = InputEvent()->GetVZEROData()->GetMultiplicity(i)*fVZEROApol[1]/fVZEROgainEqualization->GetBinContent(1+i);
1027 else if ( i < 56 ) weight = InputEvent()->GetVZEROData()->GetMultiplicity(i)*fVZEROApol[2]/fVZEROgainEqualization->GetBinContent(1+i);
1028 else if ( i < 64 ) weight = InputEvent()->GetVZEROData()->GetMultiplicity(i)*fVZEROApol[3]/fVZEROgainEqualization->GetBinContent(1+i);
1029 // fill Q-vectors for v0a side
1030 Qa2[0]+=weight*TMath::Cos(2.*phi);
1031 Qa3[0]+=weight*TMath::Cos(3.*phi);
1032 Qa2[1]+=weight*TMath::Sin(2.*phi);
1033 Qa3[1]+=weight*TMath::Sin(3.*phi);
1036 // get the cache index and read the correction terms from the cache
1037 Int_t VZEROcentralityBin(GetVZEROCentralityBin());
1038 Double_t Qx2amean = fMeanQ[VZEROcentralityBin][1][0];
1039 Double_t Qx2arms = fWidthQ[VZEROcentralityBin][1][0];
1040 Double_t Qy2amean = fMeanQ[VZEROcentralityBin][1][1];
1041 Double_t Qy2arms = fWidthQ[VZEROcentralityBin][1][1];
1043 Double_t Qx2cmean = fMeanQ[VZEROcentralityBin][0][0];
1044 Double_t Qx2crms = fWidthQ[VZEROcentralityBin][0][0];
1045 Double_t Qy2cmean = fMeanQ[VZEROcentralityBin][0][1];
1046 Double_t Qy2crms = fWidthQ[VZEROcentralityBin][0][1];
1048 Double_t Qx3amean = fMeanQv3[VZEROcentralityBin][1][0];
1049 Double_t Qx3arms = fWidthQv3[VZEROcentralityBin][1][0];
1050 Double_t Qy3amean = fMeanQv3[VZEROcentralityBin][1][1];
1051 Double_t Qy3arms = fWidthQv3[VZEROcentralityBin][1][1];
1053 Double_t Qx3cmean = fMeanQv3[VZEROcentralityBin][0][0];
1054 Double_t Qx3crms = fWidthQv3[VZEROcentralityBin][0][0];
1055 Double_t Qy3cmean = fMeanQv3[VZEROcentralityBin][0][1];
1056 Double_t Qy3crms = fWidthQv3[VZEROcentralityBin][0][1];
1058 // update the weighted q-vectors with the re-centered values
1059 Qa2[0] = (Qa2[0] - Qx2amean)/Qx2arms;
1060 Qa2[1] = (Qa2[1] - Qy2amean)/Qy2arms;
1061 Qc2[0] = (Qc2[0] - Qx2cmean)/Qx2crms;
1062 Qc2[1] = (Qc2[1] - Qy2cmean)/Qy2crms;
1064 Qa3[0] = (Qa3[0] - Qx3amean)/Qx3arms;
1065 Qa3[1] = (Qa3[1] - Qy3amean)/Qy3arms;
1066 Qc3[0] = (Qc3[0] - Qx3cmean)/Qx3crms;
1067 Qc3[1] = (Qc3[0] - Qy3cmean)/Qy3crms;
1069 //_____________________________________________________________________________
1070 void AliAnalysisTaskJetV2::CalculateQvectorCombinedVZERO(Double_t Q2[2], Double_t Q3[2]) const
1072 // calculate calibrated q-vector of the combined vzeroa, vzeroc system
1073 // this is somewhat ugly as CalculateQvectorCombinedVZERO is called more than once per event
1074 // but for now it will have to do ...
1075 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1077 // first step: retrieve the q-vectors component-wise per vzero detector
1078 Double_t QA2[] = {-999., -999.};
1079 Double_t QA3[] = {-999., -999.};
1080 Double_t QC2[] = {-999., -999.};
1081 Double_t QC3[] = {-999., -999.};
1082 CalculateQvectorVZERO(QA2, QA3, QC2, QC3);
1084 // get cache index and retrieve the chi weights for this centrality
1085 Int_t VZEROcentralityBin(GetVZEROCentralityBin());
1086 Double_t chi2A(fChi2A->At(VZEROcentralityBin));
1087 Double_t chi2C(fChi2C->At(VZEROcentralityBin));
1088 Double_t chi3A(fChi3A->At(VZEROcentralityBin));
1089 Double_t chi3C(fChi3C->At(VZEROcentralityBin));
1091 // combine the vzera and vzeroc signal
1092 Q2[0] = chi2A*chi2A*QA2[0]+chi2C*chi2C*QC2[0];
1093 Q2[1] = chi2A*chi2A*QA2[1]+chi2C*chi2C*QC2[1];
1094 Q3[0] = chi3A*chi3A*QA3[0]+chi3C*chi3C*QC3[0];
1095 Q3[1] = chi3A*chi3A*QC3[1]+chi3C*chi3C*QC3[1];
1097 //_____________________________________________________________________________
1098 void AliAnalysisTaskJetV2::CalculateRandomCone(Float_t &pt, Float_t &eta, Float_t &phi,
1099 AliParticleContainer* tracksCont, AliClusterContainer* clusterCont, AliEmcalJet* jet) const
1101 // get a random cone
1102 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1103 pt = 0; eta = 0; phi = 0;
1104 Float_t etaJet(999), phiJet(999), dJet(999); // no jet: same as jet very far away
1105 if(jet) { // if a leading jet is given, use its kinematic properties to exclude it
1106 etaJet = jet->Eta();
1107 phiJet = jet->Phi();
1109 // the random cone acceptance has to equal the jet acceptance
1110 // this also insures safety when runnnig on the semi-good tpc runs for 11h data,
1111 // where jet acceptance is adjusted to reduced acceptance - hence random cone acceptance as well
1112 Float_t minPhi(GetJetContainer()->GetJetPhiMin()), maxPhi(GetJetContainer()->GetJetPhiMax());
1113 if(maxPhi > TMath::TwoPi()) maxPhi = TMath::TwoPi();
1114 if(minPhi < 0 ) minPhi = 0.;
1115 // construct a random cone and see if it's far away enough from the leading jet
1116 Int_t attempts(1000);
1119 eta = gRandom->Uniform(GetJetContainer()->GetJetEtaMin(), GetJetContainer()->GetJetEtaMax());
1120 phi = gRandom->Uniform(minPhi, maxPhi);
1122 dJet = TMath::Sqrt((etaJet-eta)*(etaJet-eta)+(phiJet-phi)*(phiJet-phi));
1123 if(dJet > fMinDisanceRCtoLJ) break;
1124 else if (attempts == 0) {
1125 printf(" > No random cone after 1000 tries, giving up ... !\n");
1129 // get the charged energy (if tracks are provided)
1131 AliVParticle* track = tracksCont->GetNextAcceptParticle(0);
1133 Float_t etaTrack(track->Eta()), phiTrack(track->Phi());
1134 // get distance from cone
1135 if(TMath::Abs(phiTrack-phi) > TMath::Abs(phiTrack - phi + TMath::TwoPi())) phiTrack+=TMath::TwoPi();
1136 if(TMath::Abs(phiTrack-phi) > TMath::Abs(phiTrack - phi - TMath::TwoPi())) phiTrack-=TMath::TwoPi();
1137 if(TMath::Sqrt(TMath::Abs((etaTrack-eta)*(etaTrack-eta)+(phiTrack-phi)*(phiTrack-phi))) <= GetJetRadius()) pt += track->Pt();
1138 track = tracksCont->GetNextAcceptParticle();
1141 // get the neutral energy (if clusters are provided)
1143 TLorentzVector momentum;
1144 AliVCluster* cluster = clusterCont->GetNextAcceptCluster(0);
1146 cluster->GetMomentum(momentum, const_cast<Double_t*>(fVertex));
1147 Float_t etaClus(momentum.Eta()), phiClus(momentum.Phi());
1148 // get distance from cone
1149 if(TMath::Abs(phiClus-phi) > TMath::Abs(phiClus - phi + TMath::TwoPi())) phiClus+=TMath::TwoPi();
1150 if(TMath::Abs(phiClus-phi) > TMath::Abs(phiClus - phi - TMath::TwoPi())) phiClus-=TMath::TwoPi();
1151 if(TMath::Sqrt(TMath::Abs((etaClus-eta)*(etaClus-eta)+(phiClus-phi)*(phiClus-phi))) <= GetJetRadius()) pt += momentum.Pt();
1152 cluster = clusterCont->GetNextAcceptCluster();
1156 //_____________________________________________________________________________
1157 Double_t AliAnalysisTaskJetV2::CalculateQC2(Int_t harm) {
1158 // get the second order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho
1159 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1160 Double_t reQ(0), imQ(0), modQ(0), M11(0), M(0);
1161 if(fUsePtWeight) { // for the weighted 2-nd order q-cumulant
1162 QCnQnk(harm, 1, reQ, imQ); // get the weighted 2-nd order q-vectors
1163 modQ = reQ*reQ+imQ*imQ; // get abs Q-squared
1164 M11 = QCnM11(); // equals S2,1 - S1,2
1165 return (M11 > 0) ? ((modQ - QCnS(1,2))/M11) : -999;
1166 } // else return the non-weighted 2-nd order q-cumulant
1167 QCnQnk(harm, 0, reQ, imQ); // get the non-weighted 2-nd order q-vectors
1168 modQ = reQ*reQ+imQ*imQ; // get abs Q-squared
1170 return (M > 1) ? (modQ - M)/(M*(M-1)) : -999;
1172 //_____________________________________________________________________________
1173 Double_t AliAnalysisTaskJetV2::CalculateQC4(Int_t harm) {
1174 // get the fourth order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho
1175 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1176 Double_t reQn1(0), imQn1(0), reQ2n2(0), imQ2n2(0), reQn3(0), imQn3(0), M1111(0), M(0);
1177 Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0); // terms of the calculation
1178 if(fUsePtWeight) { // for the weighted 4-th order q-cumulant
1179 QCnQnk(harm, 1, reQn1, imQn1);
1180 QCnQnk(harm*2, 2, reQ2n2, imQ2n2);
1181 QCnQnk(harm, 3, reQn3, imQn3);
1182 // fill in the terms ...
1183 a = (reQn1*reQn1+imQn1*imQn1)*(reQn1*reQn1+imQn1*imQn1);
1184 b = reQ2n2*reQ2n2 + imQ2n2*imQ2n2;
1185 c = -2.*(reQ2n2*reQn1*reQn1-reQ2n2*imQn1*imQn1+2.*imQ2n2*reQn1*imQn1);
1186 d = 8.*(reQn3*reQn1+imQn3*imQn1);
1187 e = -4.*QCnS(1,2)*(reQn1*reQn1+imQn1*imQn1);
1191 return (M1111 > 0) ? (a+b+c+d+e+f+g)/M1111 : -999;
1192 } // else return the unweighted case
1193 Double_t reQn(0), imQn(0), reQ2n(0), imQ2n(0);
1194 QCnQnk(harm, 0, reQn, imQn);
1195 QCnQnk(harm*2, 0, reQ2n, imQ2n);
1196 // fill in the terms ...
1198 if(M < 4) return -999;
1199 a = (reQn*reQn+imQn*imQn)*(reQn*reQn+imQn*imQn);
1200 b = reQ2n*reQ2n + imQ2n*imQ2n;
1201 c = -2.*(reQ2n*reQn*reQn-reQ2n*imQn*imQn+2.*imQ2n*reQn*imQn);
1202 e = -4.*(M-2)*(reQn*reQn+imQn*imQn);
1204 return (a+b+c+e+f)/(M*(M-1)*(M-2)*(M-3));
1206 //_____________________________________________________________________________
1207 void AliAnalysisTaskJetV2::QCnQnk(Int_t n, Int_t k, Double_t &reQ, Double_t &imQ) {
1208 // get the weighted n-th order q-vector, pass real and imaginary part as reference
1209 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1210 if(!fTracks) return;
1211 fNAcceptedTracksQCn = 0;
1212 Int_t iTracks(fTracks->GetEntriesFast());
1213 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
1214 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
1215 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
1216 fNAcceptedTracksQCn++;
1217 // for the unweighted case, k equals zero and the weight doesn't contribute to the equation below
1218 reQ += TMath::Power(track->Pt(), k) * TMath::Cos(((double)n)*track->Phi());
1219 imQ += TMath::Power(track->Pt(), k) * TMath::Sin(((double)n)*track->Phi());
1222 //_____________________________________________________________________________
1223 void AliAnalysisTaskJetV2::QCnDiffentialFlowVectors(
1224 TClonesArray* pois, TArrayD* ptBins, Bool_t vpart, Double_t* repn, Double_t* impn,
1225 Double_t *mp, Double_t *reqn, Double_t *imqn, Double_t* mq, Int_t n)
1227 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1228 // get unweighted differential flow vectors
1229 Int_t iPois(pois->GetEntriesFast());
1231 for(Int_t i(0); i < iPois; i++) {
1232 for(Int_t ptBin(0); ptBin < ptBins->GetSize()-1; ptBin++) {
1233 AliVTrack* poi = static_cast<AliVTrack*>(pois->At(i));
1234 if(PassesCuts(poi)) {
1235 if(poi->Pt() >= ptBins->At(ptBin) && poi->Pt() < ptBins->At(ptBin+1)) {
1236 // fill the flow vectors assuming that all poi's are in the rp selection (true by design)
1237 repn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
1238 impn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
1240 reqn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
1241 imqn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
1248 for(Int_t i(0); i < iPois; i++) {
1249 for(Int_t ptBin(0); ptBin < ptBins->GetSize()-1; ptBin++) {
1250 AliEmcalJet* poi = static_cast<AliEmcalJet*>(pois->At(i));
1251 if(PassesCuts(poi)) {
1252 Double_t pt(poi->Pt()-poi->Area()*fLocalRho->GetLocalVal(poi->Phi(), GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1253 if(pt >= ptBins->At(ptBin) && pt < ptBins->At(ptBin+1)) {
1254 repn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
1255 impn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
1256 mp[ptBin]++; // qn isn't filled, no overlap between poi's and rp's
1263 //_____________________________________________________________________________
1264 Double_t AliAnalysisTaskJetV2::QCnS(Int_t i, Int_t j) {
1265 // get the weighted ij-th order autocorrelation correction
1266 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1267 if(!fTracks || i <= 0 || j <= 0) return -999;
1268 Int_t iTracks(fTracks->GetEntriesFast());
1270 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
1271 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
1272 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
1273 Sij+=TMath::Power(track->Pt(), j);
1275 return TMath::Power(Sij, i);
1277 //_____________________________________________________________________________
1278 Double_t AliAnalysisTaskJetV2::QCnM() {
1279 // get multiplicity for unweighted q-cumulants. function QCnQnk should be called first
1280 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1281 return (Double_t) fNAcceptedTracksQCn;
1283 //_____________________________________________________________________________
1284 Double_t AliAnalysisTaskJetV2::QCnM11() {
1285 // get multiplicity weights for the weighted two particle cumulant
1286 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1287 return (QCnS(2,1) - QCnS(1,2));
1289 //_____________________________________________________________________________
1290 Double_t AliAnalysisTaskJetV2::QCnM1111() {
1291 // get multiplicity weights for the weighted four particle cumulant
1292 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1293 return (QCnS(4,1)-6*QCnS(1,2)*QCnS(2,1)+8*QCnS(1,3)*QCnS(1,1)+3*QCnS(2,2)-6*QCnS(1,4));
1295 //_____________________________________________________________________________
1296 Bool_t AliAnalysisTaskJetV2::QCnRecovery(Double_t psi2, Double_t psi3) {
1297 // decides how to deal with the situation where c2 or c3 is negative
1298 // returns kTRUE depending on whether or not a modulated rho is used for the jet background
1299 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1300 if(TMath::AreEqualAbs(fFitModulation->GetParameter(3), .0, 1e-10) && TMath::AreEqualAbs(fFitModulation->GetParameter(7), .0,1e-10)) {
1301 fFitModulation->SetParameter(7, 0);
1302 fFitModulation->SetParameter(3, 0);
1303 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1304 return kTRUE; // v2 and v3 have physical null values
1306 switch (fQCRecovery) {
1307 case kFixedRho : { // roll back to the original rho
1308 fFitModulation->SetParameter(7, 0);
1309 fFitModulation->SetParameter(3, 0);
1310 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1311 return kFALSE; // rho is forced to be fixed
1313 case kNegativeVn : {
1314 Double_t c2(fFitModulation->GetParameter(3));
1315 Double_t c3(fFitModulation->GetParameter(7));
1316 if( c2 < 0 ) c2 = -1.*TMath::Sqrt(-1.*c2);
1317 if( c3 < 0 ) c3 = -1.*TMath::Sqrt(-1.*c3);
1318 fFitModulation->SetParameter(3, c2);
1319 fFitModulation->SetParameter(7, c3);
1320 return kTRUE; // is this a physical quantity ?
1323 fitModulationType tempType(fFitModulationType); // store temporarily
1324 fFitModulationType = kCombined;
1325 fFitModulation->SetParameter(7, 0);
1326 fFitModulation->SetParameter(3, 0);
1327 Bool_t pass(CorrectRho(psi2, psi3)); // do the fit and all quality checks
1328 fFitModulationType = tempType; // roll back for next event
1331 default : return kFALSE;
1335 //_____________________________________________________________________________
1336 Bool_t AliAnalysisTaskJetV2::CorrectRho(Double_t psi2, Double_t psi3)
1338 // get rho' -> rho(phi)
1339 // two routines are available, both can be used with or without pt weights
1340 // [1] get vn from q-cumulants or as an integrated value from a user supplied histogram
1341 // in case of cumulants, both cumulants and vn values are stored. in both cases, v2 and v3
1342 // are expected. a check is performed to see if rho has no negative local minimum
1343 // for full description, see Phys. Rev. C 83, 044913
1344 // since the cn distribution has negative values, vn = sqrt(cn) can be imaginary sometimes
1345 // in this case one can either roll back to the 'original' rixed rho, do a fit for vn or take use
1346 // vn = - sqrt(|cn|)
1347 // [2] fitting a fourier expansion to the de/dphi distribution
1348 // the fit can be done with either v2, v3 or a combination.
1349 // in all cases, a cut can be made on the p-value of the chi-squared value of the fit
1350 // and a check can be performed to see if rho has no negative local minimum
1351 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1352 Int_t freeParams(2); // free parameters of the fit (for NDF)
1353 switch (fFitModulationType) { // for approaches where no fitting is required
1355 fFitModulation->FixParameter(4, psi2);
1356 fFitModulation->FixParameter(6, psi3);
1357 fFitModulation->FixParameter(3, CalculateQC2(2)); // set here with cn, vn = sqrt(cn)
1358 fFitModulation->FixParameter(7, CalculateQC2(3));
1359 // first fill the histos of the raw cumulant distribution
1360 if (fUsePtWeight) { // use weighted weights
1361 Double_t dQCnM11 = (fNoEventWeightsForQC) ? 1. : QCnM11();
1362 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), dQCnM11);
1363 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), dQCnM11);
1365 Double_t dQCnM = (fNoEventWeightsForQC) ? 2. : QCnM();
1366 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), dQCnM*(dQCnM-1));
1367 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), dQCnM*(dQCnM-1));
1369 // then see if one of the cn value is larger than zero and vn is readily available
1370 if(fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) {
1371 fFitModulation->FixParameter(3, TMath::Sqrt(fFitModulation->GetParameter(3)));
1372 fFitModulation->FixParameter(7, TMath::Sqrt(fFitModulation->GetParameter(7)));
1373 } else if (!QCnRecovery(psi2, psi3)) return kFALSE; // try to recover the cumulant, this will set v2 and v3
1374 if(fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { // general check
1375 fFitModulation->SetParameter(7, 0);
1376 fFitModulation->SetParameter(3, 0);
1377 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1383 fFitModulation->FixParameter(4, psi2);
1384 fFitModulation->FixParameter(6, psi3);
1385 fFitModulation->FixParameter(3, CalculateQC4(2)); // set here with cn, vn = sqrt(cn)
1386 fFitModulation->FixParameter(7, CalculateQC4(3));
1387 // first fill the histos of the raw cumulant distribution
1388 if (fUsePtWeight) { // use weighted weights
1389 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
1390 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
1392 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
1393 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
1395 // then see if one of the cn value is larger than zero and vn is readily available
1396 if(fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) {
1397 fFitModulation->FixParameter(3, TMath::Sqrt(fFitModulation->GetParameter(3)));
1398 fFitModulation->FixParameter(7, TMath::Sqrt(fFitModulation->GetParameter(7)));
1399 } else if (!QCnRecovery(psi2, psi3)) return kFALSE; // try to recover the cumulant, this will set v2 and v3
1400 if(fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { // general check
1401 fFitModulation->SetParameter(7, 0);
1402 fFitModulation->SetParameter(3, 0);
1403 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1407 case kIntegratedFlow : {
1408 // use v2 and v3 values from an earlier iteration over the data
1409 fFitModulation->FixParameter(3, fUserSuppliedV2->GetBinContent(fUserSuppliedV2->GetXaxis()->FindBin(fCent)));
1410 fFitModulation->FixParameter(4, psi2);
1411 fFitModulation->FixParameter(6, psi3);
1412 fFitModulation->FixParameter(7, fUserSuppliedV3->GetBinContent(fUserSuppliedV3->GetXaxis()->FindBin(fCent)));
1413 if(fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) {
1414 fFitModulation->SetParameter(7, 0);
1415 fFitModulation->SetParameter(3, 0);
1416 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1423 TString detector("");
1424 switch (fDetectorType) {
1425 case kTPC : detector+="TPC";
1427 case kVZEROA : detector+="VZEROA";
1429 case kVZEROC : detector+="VZEROC";
1431 case kVZEROComb : detector+="VZEROComb";
1433 case kFixedEP : detector+="FixedEP";
1437 Int_t iTracks(fTracks->GetEntriesFast());
1438 Double_t excludeInEta = -999;
1439 Double_t excludeInPhi = -999;
1440 Double_t excludeInPt = -999;
1441 if(iTracks <= 0 || fLocalRho->GetVal() <= 0 ) return kFALSE; // no use fitting an empty event ...
1442 if(fExcludeLeadingJetsFromFit > 0 ) {
1444 excludeInEta = fLeadingJet->Eta();
1445 excludeInPhi = fLeadingJet->Phi();
1446 excludeInPt = fLeadingJet->Pt();
1449 // check the acceptance of the track selection that will be used
1450 // if one uses e.g. semi-good tpc tracks, accepance in phi is reduced to 0 < phi < 4
1451 // the defaults (-10 < phi < 10) which accept all, are then overwritten
1452 Double_t lowBound(0.), upBound(TMath::TwoPi()); // bounds for fit
1453 if(GetParticleContainer()->GetParticlePhiMin() > lowBound) lowBound = GetParticleContainer()->GetParticlePhiMin();
1454 if(GetParticleContainer()->GetParticlePhiMax() < upBound) upBound = GetParticleContainer()->GetParticlePhiMax();
1455 fHistSwap->Reset(); // clear the histogram
1456 TH1F _tempSwap; // on stack for quick access
1457 TH1F _tempSwapN; // on stack for quick access, bookkeeping histogram
1458 if(fRebinSwapHistoOnTheFly) {
1459 if(fNAcceptedTracks < 49) fNAcceptedTracks = 49; // avoid aliasing effects
1460 _tempSwap = TH1F("_tempSwap", "_tempSwap", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), lowBound, upBound);
1461 if(fUsePtWeightErrorPropagation) _tempSwapN = TH1F("_tempSwapN", "_tempSwapN", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), lowBound, upBound);
1462 if(fUsePtWeight) _tempSwap.Sumw2();
1464 else _tempSwap = *fHistSwap; // now _tempSwap holds the desired histo
1465 // non poissonian error when using pt weights
1466 Double_t totalpts(0.), totalptsquares(0.), totalns(0.);
1467 for(Int_t i(0); i < iTracks; i++) {
1468 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1469 if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta) < GetJetContainer()->GetJetRadius()*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - GetJetContainer()->GetJetRadius() - GetJetContainer()->GetJetEtaMax() ) > 0 )) continue;
1470 if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
1472 _tempSwap.Fill(track->Phi(), track->Pt());
1473 if(fUsePtWeightErrorPropagation) {
1474 totalpts += track->Pt();
1475 totalptsquares += track->Pt()*track->Pt();
1477 _tempSwapN.Fill(track->Phi());
1480 else _tempSwap.Fill(track->Phi());
1482 if(fUsePtWeight && fUsePtWeightErrorPropagation) {
1483 // in the case of pt weights overwrite the poissonian error estimate which is assigned by root by a more sophisticated appraoch
1484 // the assumption here is that the bin error will be dominated by the uncertainty in the mean pt in a bin and in the uncertainty
1485 // of the number of tracks in a bin, the first of which will be estimated from the sample standard deviation of all tracks in the
1486 // event, for the latter use a poissonian estimate. the two contrubitions are assumed to be uncorrelated
1487 if(totalns < 2) return kFALSE; // not one track passes the cuts > 2 avoids possible division by 0 later on
1488 for(Int_t l = 0; l < _tempSwap.GetNbinsX(); l++) {
1489 if(_tempSwapN.GetBinContent(l+1) == 0) {
1490 _tempSwap.SetBinContent(l+1,0);
1491 _tempSwap.SetBinError(l+1,0);
1494 Double_t vartimesnsq = totalptsquares*totalns - totalpts*totalpts;
1495 Double_t variance = vartimesnsq/(totalns*(totalns-1.));
1496 Double_t SDOMSq = variance / _tempSwapN.GetBinContent(l+1);
1497 Double_t SDOMSqOverMeanSq = SDOMSq * _tempSwapN.GetBinContent(l+1) * _tempSwapN.GetBinContent(l+1) / (_tempSwapN.GetBinContent(l+1) * _tempSwapN.GetBinContent(l+1));
1498 Double_t poissonfrac = 1./_tempSwapN.GetBinContent(l+1);
1499 Double_t vartotalfrac = SDOMSqOverMeanSq + poissonfrac;
1500 Double_t vartotal = vartotalfrac * _tempSwap.GetBinContent(l+1) * _tempSwap.GetBinContent(l+1);
1501 if(vartotal > 0.0001) _tempSwap.SetBinError(l+1,TMath::Sqrt(vartotal));
1503 _tempSwap.SetBinContent(l+1,0);
1504 _tempSwap.SetBinError(l+1,0);
1509 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1510 switch (fFitModulationType) {
1512 fFitModulation->FixParameter(0, fLocalRho->GetVal() );
1516 fFitModulation->FixParameter(4, psi2);
1520 fFitModulation->FixParameter(4, psi3);
1524 fFitModulation->FixParameter(4, psi2);
1525 fFitModulation->FixParameter(6, psi3);
1528 case kFourierSeries : {
1529 // in this approach, an explicit calculation will be made of vn = sqrt(xn^2+yn^2)
1530 // where x[y] = Integrate[r(phi)cos[sin](n phi)dphi, 0, 2pi]
1531 Double_t cos2(0), sin2(0), cos3(0), sin3(0), sumPt(0);
1532 for(Int_t i(0); i < iTracks; i++) {
1533 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1534 if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
1535 sumPt += track->Pt();
1536 cos2 += track->Pt()*TMath::Cos(2*PhaseShift(track->Phi()-psi2));
1537 sin2 += track->Pt()*TMath::Sin(2*PhaseShift(track->Phi()-psi2));
1538 cos3 += track->Pt()*TMath::Cos(3*PhaseShift(track->Phi()-psi3));
1539 sin3 += track->Pt()*TMath::Sin(3*PhaseShift(track->Phi()-psi3));
1541 fFitModulation->SetParameter(3, TMath::Sqrt(cos2*cos2+sin2*sin2)/fLocalRho->GetVal());
1542 fFitModulation->SetParameter(4, psi2);
1543 fFitModulation->SetParameter(6, psi3);
1544 fFitModulation->SetParameter(7, TMath::Sqrt(cos3*cos3+sin3*sin3)/fLocalRho->GetVal());
1549 // toy mc, just here to check procedure, azimuthal profile is filled from hypothesis so p-value distribution should be flat
1550 Int_t _bins = _tempSwap.GetXaxis()->GetNbins();
1551 TF1* _tempFit = new TF1("temp_fit_kCombined", "[0]*([1]+[2]*([3]*TMath::Cos([2]*(x-[4]))+[7]*TMath::Cos([5]*(x-[6]))))", 0, TMath::TwoPi());
1552 _tempFit->SetParameter(0, fFitModulation->GetParameter(0)); // normalization
1553 _tempFit->SetParameter(3, 0.1); // v2
1554 _tempFit->FixParameter(1, 1.); // constant
1555 _tempFit->FixParameter(2, 2.); // constant
1556 _tempFit->FixParameter(5, 3.); // constant
1557 _tempFit->FixParameter(4, fFitModulation->GetParameter(4));
1558 _tempFit->FixParameter(6, fFitModulation->GetParameter(6));
1559 _tempFit->SetParameter(7, 0.1); // v3
1560 _tempSwap.Reset(); // rese bin content
1561 for(int _binsI = 0; _binsI < _bins*_bins; _binsI++) _tempSwap.Fill(_tempFit->GetRandom());
1563 _tempSwap.Fit(fFitModulation, fFitModulationOptions.Data(), "", lowBound, upBound);
1564 // the quality of the fit is evaluated from 1 - the cdf of the chi square distribution
1565 // three methods are available, all with their drawbacks. all are stored, one is selected to do the cut
1566 Int_t NDF(_tempSwap.GetXaxis()->GetNbins()-freeParams);
1567 if(NDF == 0 || (float)NDF <= 0.) return kFALSE;
1568 Double_t CDF(1.-ChiSquareCDF(NDF, ChiSquare(_tempSwap, fFitModulation)));
1569 Double_t CDFROOT(1.-ChiSquareCDF(NDF, fFitModulation->GetChisquare()));
1570 Double_t CDFKolmogorov(KolmogorovTest(_tempSwap, fFitModulation));
1571 // fill the values and centrality correlation (redundant but easy on the eyes)
1572 fHistPvalueCDF->Fill(CDF);
1573 fHistPvalueCDFCent->Fill(fCent, CDF);
1574 fHistPvalueCDFROOT->Fill(CDFROOT);
1575 fHistPvalueCDFROOTCent->Fill(fCent, CDFROOT);
1576 fHistKolmogorovTest->Fill(CDFKolmogorov);
1577 fHistChi2ROOTCent->Fill(fCent, fFitModulation->GetChisquare()/((float)NDF));
1578 fHistChi2Cent->Fill(fCent, ChiSquare(_tempSwap, fFitModulation)/((float)NDF));
1579 fHistKolmogorovTestCent->Fill(fCent, CDFKolmogorov);
1580 fHistPChi2Root->Fill(CDFROOT, fFitModulation->GetChisquare()/((float)NDF));
1581 fHistPChi2->Fill(CDF, ChiSquare(_tempSwap, fFitModulation)/((float)NDF));
1582 fHistPKolmogorov->Fill(CDF, CDFKolmogorov);
1584 // variable CDF is used for making cuts, so we fill it with the selected p-value
1585 switch (fFitGoodnessTest) {
1589 case kChi2Poisson : break; // CDF is already CDF
1590 case kKolmogorov : {
1591 CDF = CDFKolmogorov;
1597 // as an additional quality check, see if fitting a control fit has a higher significance
1598 _tempSwap.Fit(fFitControl, fFitModulationOptions.Data(), "", lowBound, upBound);
1599 Double_t CDFControl(-1.);
1600 switch (fFitGoodnessTest) {
1602 CDFControl = 1.-ChiSquareCDF(fFitControl->GetNDF(), fFitModulation->GetChisquare());
1604 case kChi2Poisson : {
1605 CDFControl = 1.-ChiSquareCDF(fFitControl->GetNDF(), ChiSquare(_tempSwap, fFitModulation));
1607 case kKolmogorov : {
1608 CDFControl = KolmogorovTest(_tempSwap, fFitControl);
1612 if(CDFControl > CDF) {
1613 CDF = -1.; // control fit is more significant, so throw out the 'old' fit
1614 fHistRhoStatusCent->Fill(fCent, -1);
1617 if(CDF >= fMinPvalue && CDF <= fMaxPvalue && ( fFitModulation->GetMinimum(0, TMath::TwoPi()) > 0)) {
1618 // fit quality. not that although with limited acceptance the fit is performed on just
1619 // part of phase space, the requirement that energy desntiy is larger than zero is applied
1620 // to the FULL spectrum
1621 fHistRhoStatusCent->Fill(fCent, 0.);
1622 // for LOCAL didactic purposes, save the best and the worst fits
1623 // this routine can produce a lot of output histograms (it's not memory 'safe') and will not work on GRID
1624 // since the output will become unmergeable (i.e. different nodes may produce conflicting output)
1625 switch (fRunModeType) {
1627 if(fRandom->Uniform(0, 100) > fPercentageOfFits) break;
1628 static Int_t didacticCounterBest(0);
1629 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
1630 TF1* didacticFit = (TF1*)fFitModulation->Clone(Form("fit_%i_CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
1631 switch(fFitModulationType) {
1633 // to make a nice picture also plot the separate components (v2 and v3) of the fit
1634 // only done for cobined fit where there are actually components to split ...
1635 TF1* v0(new TF1("dfit_kV2", "[0]", 0, TMath::TwoPi()));
1636 v0->SetParameter(0, didacticFit->GetParameter(0)); // normalization
1637 v0->SetLineColor(kMagenta);
1638 v0->SetLineStyle(7);
1639 didacticProfile->GetListOfFunctions()->Add(v0);
1640 TF1* v2(new TF1("dfit_kV2", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
1641 v2->SetParameter(0, didacticFit->GetParameter(0)); // normalization
1642 v2->SetParameter(3, didacticFit->GetParameter(3)); // v2
1643 v2->FixParameter(1, 1.); // constant
1644 v2->FixParameter(2, 2.); // constant
1645 v2->FixParameter(4, didacticFit->GetParameter(4)); // psi2
1646 v2->SetLineColor(kGreen);
1647 didacticProfile->GetListOfFunctions()->Add(v2);
1648 TF1* v3(new TF1("dfit_kV3", "[0]*([1]+[2]*[3]*TMath::Cos([5]*(x-[4])))", 0, TMath::TwoPi()));
1649 v3->SetParameter(0, didacticFit->GetParameter(0)); // normalization
1650 v3->SetParameter(3, didacticFit->GetParameter(7)); // v3
1651 v3->FixParameter(1, 1.); // constant
1652 v3->FixParameter(2, 2.); // constant
1653 v3->FixParameter(4, didacticFit->GetParameter(6)); // psi3
1654 v3->FixParameter(5, 3.); // constant
1655 v3->SetLineColor(kCyan);
1656 didacticProfile->GetListOfFunctions()->Add(v3);
1660 didacticProfile->GetListOfFunctions()->Add(didacticFit);
1661 didacticProfile->GetYaxis()->SetTitle("#frac{d #sum #it{p}_{T}}{d #varphi} [GeV/#it{c}]");
1662 didacticProfile->GetXaxis()->SetTitle("#varphi");
1663 fOutputListGood->Add(didacticProfile);
1664 didacticCounterBest++;
1665 TH2F* didacticSurface = BookTH2F(Form("surface_%s", didacticProfile->GetName()), "#phi", "#eta", 50, 0, TMath::TwoPi(), 50, -1, 1, -1, kFALSE);
1666 for(Int_t i(0); i < iTracks; i++) {
1667 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1668 if(PassesCuts(track)) {
1669 if(fUsePtWeight) didacticSurface->Fill(track->Phi(), track->Eta(), track->Pt());
1670 else didacticSurface->Fill(track->Phi(), track->Eta());
1673 if(fExcludeLeadingJetsFromFit) { // visualize the excluded region
1674 TF2 *f2 = new TF2(Form("%s_LJ", didacticSurface->GetName()),"[0]*TMath::Gaus(x,[1],[2])*TMath::Gaus(y,[3],[4])", 0, TMath::TwoPi(), -1, 1);
1675 f2->SetParameters(excludeInPt/3.,excludeInPhi,.1,excludeInEta,.1);
1676 didacticSurface->GetListOfFunctions()->Add(f2);
1678 fOutputListGood->Add(didacticSurface);
1682 } else { // if the fit is of poor quality revert to the original rho estimate
1683 switch (fRunModeType) { // again see if we want to save the fit
1685 static Int_t didacticCounterWorst(0);
1686 if(fRandom->Uniform(0, 100) > fPercentageOfFits) break;
1687 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data() ));
1688 TF1* didacticFit = (TF1*)fFitModulation->Clone(Form("fit_%i_p_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data()));
1689 didacticProfile->GetListOfFunctions()->Add(didacticFit);
1690 fOutputListBad->Add(didacticProfile);
1691 didacticCounterWorst++;
1695 switch (fFitModulationType) {
1696 case kNoFit : break; // nothing to do
1697 case kCombined : fFitModulation->SetParameter(7, 0); // no break
1698 case kFourierSeries : fFitModulation->SetParameter(7, 0); // no break
1699 default : { // needs to be done if there was a poor fit
1700 fFitModulation->SetParameter(3, 0);
1701 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1704 if(CDF > -.5) fHistRhoStatusCent->Fill(fCent, 1.);
1705 return kFALSE; // return false if the fit is rejected
1709 //_____________________________________________________________________________
1710 Bool_t AliAnalysisTaskJetV2::PassesCuts(AliVEvent* event)
1713 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1714 switch (fCollisionType) {
1716 fInCentralitySelection = 0;
1721 if(!event || !AliAnalysisTaskEmcal::IsEventSelected()) return kFALSE;
1722 if(TMath::Abs(InputEvent()->GetPrimaryVertex()->GetZ()) > 10.) return kFALSE;
1723 // aod and esd specific checks
1724 switch (fDataType) {
1726 AliESDEvent* esdEvent = static_cast<AliESDEvent*>(InputEvent());
1727 if( (!esdEvent) || (TMath::Abs(esdEvent->GetPrimaryVertexSPD()->GetZ() - esdEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE;
1730 AliAODEvent* aodEvent = static_cast<AliAODEvent*>(InputEvent());
1731 if( (!aodEvent) || (TMath::Abs(aodEvent->GetPrimaryVertexSPD()->GetZ() - aodEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE;
1735 fCent = InputEvent()->GetCentrality()->GetCentralityPercentile("V0M");
1736 if(fCent <= fCentralityClasses->At(0) || fCent >= fCentralityClasses->At(fCentralityClasses->GetSize()-1) || TMath::Abs(fCent-InputEvent()->GetCentrality()->GetCentralityPercentile("TRK")) > 5.) return kFALSE;
1737 // determine centrality class
1738 fInCentralitySelection = -1;
1739 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
1740 if(fCent >= fCentralityClasses->At(i) && fCent <= fCentralityClasses->At(1+i)) {
1741 fInCentralitySelection = i;
1744 if(fInCentralitySelection<0) return kFALSE; // should be null op
1745 // see if input containers are filled
1746 if(fTracks->GetEntries() < 1) return kFALSE;
1747 if(fRho->GetVal() <= 0 ) return kFALSE;
1748 if(fAnalysisType == AliAnalysisTaskJetV2::kFull && !fClusterCont) return kFALSE;
1751 //_____________________________________________________________________________
1752 void AliAnalysisTaskJetV2::FillHistogramsAfterSubtraction(Double_t psi2, Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc)
1755 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1756 FillTrackHistograms();
1757 if(fAnalysisType == AliAnalysisTaskJetV2::kFull) FillClusterHistograms();
1758 FillJetHistograms(psi2);
1759 if(fFillQAHistograms) FillEventPlaneHistograms(vzero, vzeroComb, tpc);
1760 FillRhoHistograms();
1761 FillDeltaPtHistograms(psi2);
1763 //_____________________________________________________________________________
1764 void AliAnalysisTaskJetV2::FillTrackHistograms() const
1766 // fill track histograms
1767 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1768 Int_t iTracks(fTracks->GetEntriesFast()), iAcceptedTracks(0);
1769 for(Int_t i(0); i < iTracks; i++) {
1770 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1771 if(!PassesCuts(track)) continue;
1773 fHistPicoTrackPt[fInCentralitySelection]->Fill(track->Pt());
1774 if(fFillQAHistograms) FillQAHistograms(track);
1776 fHistPicoTrackMult[fInCentralitySelection]->Fill(iAcceptedTracks);
1778 //_____________________________________________________________________________
1779 void AliAnalysisTaskJetV2::FillClusterHistograms() const
1781 // fill cluster histograms
1782 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1783 if(!fClusterCont) return;
1784 Int_t iClusters(fClusterCont->GetNClusters());
1785 TLorentzVector clusterLorentzVector;
1786 for(Int_t i(0); i < iClusters; i++) {
1787 AliVCluster* cluster = fClusterCont->GetCluster(i);
1788 if (!PassesCuts(cluster)) continue;
1789 cluster->GetMomentum(clusterLorentzVector, const_cast<Double_t*>(fVertex));
1790 fHistClusterPt[fInCentralitySelection]->Fill(clusterLorentzVector.Pt());
1791 fHistClusterEtaPhi[fInCentralitySelection]->Fill(clusterLorentzVector.Eta(), clusterLorentzVector.Phi());
1792 fHistClusterEtaPhiWeighted[fInCentralitySelection]->Fill(clusterLorentzVector.Eta(), clusterLorentzVector.Phi(), clusterLorentzVector.Pt());
1796 //_____________________________________________________________________________
1797 void AliAnalysisTaskJetV2::FillEventPlaneHistograms(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc) const
1799 // fill event plane histograms, only called in qa mode
1800 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1801 fHistPsiControl->Fill(0.5, vzero[0][0]); // vzero a psi2
1802 fHistPsiControl->Fill(1.5, vzero[1][0]); // vzero c psi2
1803 fHistPsiControl->Fill(2.5, tpc[0]); // tpc psi 2
1804 fHistPsiControl->Fill(5.5, vzero[0][1]); // vzero a psi3
1805 fHistPsiControl->Fill(6.5, vzero[1][1]); // vzero b psi3
1806 fHistPsiControl->Fill(7.5, tpc[1]); // tpc psi 3
1807 fHistPsiVZEROA->Fill(vzero[0][0]);
1808 fHistPsiVZEROC->Fill(vzero[1][0]);
1809 fHistPsiVZERO->Fill(vzeroComb[0]);
1810 fHistPsiTPC->Fill(tpc[0]);
1811 fHistPsiSpread->Fill(0.5, TMath::Abs(vzero[0][0]-vzero[1][0]));
1812 fHistPsiSpread->Fill(1.5, TMath::Abs(vzero[0][0]-tpc[0]));
1813 fHistPsiSpread->Fill(2.5, TMath::Abs(vzero[1][0]-tpc[0]));
1814 // event plane vs centrality QA histo's to check recentering
1815 Double_t TRK(InputEvent()->GetCentrality()->GetCentralityPercentile("TRK"));
1816 Double_t V0M(InputEvent()->GetCentrality()->GetCentralityPercentile("V0M"));
1817 fHistPsiVZEROAV0M->Fill(V0M, vzero[0][0]);
1818 fHistPsiVZEROCV0M->Fill(V0M, vzero[1][0]);
1819 fHistPsiVZEROVV0M->Fill(V0M, vzeroComb[0]);
1820 fHistPsiTPCV0M->Fill(V0M, tpc[0]);
1821 fHistPsiVZEROATRK->Fill(TRK, vzero[0][0]);
1822 fHistPsiVZEROCTRK->Fill(TRK, vzero[1][0]);
1823 fHistPsiVZEROTRK->Fill(TRK, vzeroComb[0]);
1824 fHistPsiTPCTRK->Fill(TRK, tpc[0]);
1825 // leading jet vs event plane bias
1827 Double_t rho(fLocalRho->GetLocalVal(fLeadingJet->Phi(), GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1828 Double_t pt(fLeadingJet->Pt() - fLeadingJet->Area()*rho);
1829 fHistPsiTPCLeadingJet[fInCentralitySelection]->Fill(pt, tpc[0], fLeadingJet->Phi());
1830 fHistPsiVZEROALeadingJet[fInCentralitySelection]->Fill(pt, vzero[0][0], fLeadingJet->Phi());
1831 fHistPsiVZEROCLeadingJet[fInCentralitySelection]->Fill(pt, vzero[1][0], fLeadingJet->Phi());
1832 fHistPsiVZEROCombLeadingJet[fInCentralitySelection]->Fill(pt, vzeroComb[0], fLeadingJet->Phi());
1834 // correlation of event planes
1835 fHistPsi2Correlation[fInCentralitySelection]->Fill(tpc[0], vzero[0][0], vzero[1][0]);
1837 //_____________________________________________________________________________
1838 void AliAnalysisTaskJetV2::FillRhoHistograms()
1840 // fill rho histograms
1841 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1842 fHistRhoPackage[fInCentralitySelection]->Fill(fLocalRho->GetVal()); // save the rho estimate from the emcal jet package
1843 // get multiplicity FIXME inefficient
1844 Int_t iJets(fJets->GetEntriesFast());
1845 Double_t rho(fLocalRho->GetLocalVal(TMath::Pi(), TMath::Pi(), fLocalRho->GetVal()));
1846 fHistRho[fInCentralitySelection]->Fill(rho);
1847 fHistRhoVsMult->Fill(fTracks->GetEntries(), rho);
1848 fHistRhoVsCent->Fill(fCent, rho);
1849 for(Int_t i(0); i < iJets; i++) {
1850 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
1851 if(!PassesCuts(jet)) continue;
1852 fHistRhoAVsMult->Fill(fTracks->GetEntries(), rho * jet->Area());
1853 fHistRhoAVsCent->Fill(fCent, rho * jet->Area());
1856 //_____________________________________________________________________________
1857 void AliAnalysisTaskJetV2::FillDeltaPtHistograms(Double_t psi2) const
1859 // fill delta pt histograms
1860 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1862 const Float_t areaRC = GetJetRadius()*GetJetRadius()*TMath::Pi();
1863 // we're retrieved the leading jet, now get a random cone
1864 for(i = 0; i < fMaxCones; i++) {
1865 Float_t pt(0), eta(0), phi(0);
1866 // get a random cone without constraints on leading jet position
1867 CalculateRandomCone(pt, eta, phi, fTracksCont, fClusterCont, 0x0);
1869 if(fFillQAHistograms) fHistRCPhiEta[fInCentralitySelection]->Fill(phi, eta);
1870 fHistRhoVsRCPt[fInCentralitySelection]->Fill(pt, fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal())*areaRC);
1871 fHistRCPt[fInCentralitySelection]->Fill(pt);
1872 fHistDeltaPtDeltaPhi2[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1873 fHistDeltaPtDeltaPhi2Rho0[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetVal());
1876 // get a random cone excluding leading jet area
1877 CalculateRandomCone(pt, eta, phi, fTracksCont, fClusterCont, fLeadingJet);
1879 if(fFillQAHistograms) fHistRCPhiEtaExLJ[fInCentralitySelection]->Fill(phi, eta);
1880 fHistRhoVsRCPtExLJ[fInCentralitySelection]->Fill(pt, fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal())*areaRC);
1881 fHistRCPtExLJ[fInCentralitySelection]->Fill(pt);
1882 fHistDeltaPtDeltaPhi2ExLJ[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1883 fHistDeltaPtDeltaPhi2ExLJRho0[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetVal());
1887 //_____________________________________________________________________________
1888 void AliAnalysisTaskJetV2::FillJetHistograms(Double_t psi2)
1890 // fill jet histograms
1891 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1892 Int_t iJets(fJets->GetEntriesFast());
1893 for(Int_t i(0); i < iJets; i++) {
1894 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
1895 if(PassesCuts(jet)) {
1896 Double_t pt(jet->Pt()), area(jet->Area()), eta(jet->Eta()), phi(jet->Phi());
1897 Double_t rho(fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1898 fHistJetPtRaw[fInCentralitySelection]->Fill(pt);
1899 fHistJetPt[fInCentralitySelection]->Fill(pt-area*rho);
1900 if(fFillQAHistograms) fHistJetEtaPhi[fInCentralitySelection]->Fill(eta, phi);
1901 fHistJetPtArea[fInCentralitySelection]->Fill(pt-area*rho, area);
1902 fHistJetPtEta[fInCentralitySelection]->Fill(pt-area*rho, eta);
1903 fHistJetPsi2Pt[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt-area*rho);
1904 fHistJetPsi2PtRho0[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt-area*fLocalRho->GetVal());
1905 fHistJetPtConstituents[fInCentralitySelection]->Fill(pt-area*rho, jet->GetNumberOfConstituents());
1906 fHistJetEtaRho[fInCentralitySelection]->Fill(eta, pt/area);
1910 //_____________________________________________________________________________
1911 void AliAnalysisTaskJetV2::FillQAHistograms(AliVTrack* vtrack) const
1913 // fill qa histograms for pico tracks
1914 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1916 AliPicoTrack* track = static_cast<AliPicoTrack*>(vtrack);
1917 fHistRunnumbersPhi->Fill(fMappedRunNumber, track->Phi());
1918 fHistRunnumbersEta->Fill(fMappedRunNumber, track->Eta());
1919 Int_t type((int)(track->GetTrackType()));
1922 fHistPicoCat1[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1925 fHistPicoCat2[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1928 fHistPicoCat3[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1933 //_____________________________________________________________________________
1934 void AliAnalysisTaskJetV2::FillQAHistograms(AliVEvent* vevent)
1936 // fill qa histograms for events
1937 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1939 fHistVertexz->Fill(vevent->GetPrimaryVertex()->GetZ());
1940 fHistCentrality->Fill(fCent);
1941 Int_t runNumber(InputEvent()->GetRunNumber());
1942 for(fMappedRunNumber = 0; fMappedRunNumber < fExpectedRuns->GetSize(); fMappedRunNumber++) {
1943 if(fExpectedRuns->At(fMappedRunNumber) == runNumber) return;
1945 if(fDebug > 0) printf("\n > TASK %s CANNOT IDENTIFY RUN - CONFIGURATION COULD BE INCORRECT < \n", GetName());
1947 //_____________________________________________________________________________
1948 void AliAnalysisTaskJetV2::FillAnalysisSummaryHistogram() const
1950 // fill the analysis summary histrogram, saves all relevant analysis settigns
1951 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1952 fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fJetRadius");
1953 fHistAnalysisSummary->SetBinContent(2, GetJetContainer()->GetJetRadius());
1954 fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fJetEtaMin");
1955 fHistAnalysisSummary->SetBinContent(3, GetJetContainer()->GetJetEtaMin());
1956 fHistAnalysisSummary->GetXaxis()->SetBinLabel(4, "fJetEtaMax");
1957 fHistAnalysisSummary->SetBinContent(4, GetJetContainer()->GetJetEtaMax());
1958 fHistAnalysisSummary->GetXaxis()->SetBinLabel(5, "fJetPhiMin");
1959 fHistAnalysisSummary->SetBinContent(5, GetJetContainer()->GetJetPhiMin());
1960 fHistAnalysisSummary->GetXaxis()->SetBinLabel(6, "fJetPhiMax");
1961 fHistAnalysisSummary->SetBinContent(6, GetJetContainer()->GetJetPhiMin());
1962 fHistAnalysisSummary->GetXaxis()->SetBinLabel(16, "fForceBeamType");
1963 fHistAnalysisSummary->SetBinContent(16, fForceBeamType);
1964 fHistAnalysisSummary->GetXaxis()->SetBinLabel(17, "fMinCent");
1965 fHistAnalysisSummary->SetBinContent(17, fMinCent);
1966 fHistAnalysisSummary->GetXaxis()->SetBinLabel(18, "fMaxCent");
1967 fHistAnalysisSummary->SetBinContent(18, fMaxCent);
1968 fHistAnalysisSummary->GetXaxis()->SetBinLabel(19, "fMinVz");
1969 fHistAnalysisSummary->SetBinContent(19, fMinVz);
1970 fHistAnalysisSummary->GetXaxis()->SetBinLabel(20, "fMaxVz");
1971 fHistAnalysisSummary->SetBinContent(20, fMaxVz);
1972 fHistAnalysisSummary->GetXaxis()->SetBinLabel(21, "fOffTrigger");
1973 fHistAnalysisSummary->SetBinContent(21, fOffTrigger);
1974 fHistAnalysisSummary->GetXaxis()->SetBinLabel(34, "fitModulationType");
1975 fHistAnalysisSummary->SetBinContent(34, (int)fFitModulationType);
1976 fHistAnalysisSummary->GetXaxis()->SetBinLabel(35, "runModeType");
1977 fHistAnalysisSummary->SetBinContent(35, (int)fRunModeType);
1978 fHistAnalysisSummary->GetXaxis()->SetBinLabel(36, "data type");
1979 fHistAnalysisSummary->SetBinContent(36, (int)fDataType);
1980 fHistAnalysisSummary->GetXaxis()->SetBinLabel(37, "iterator");
1981 fHistAnalysisSummary->SetBinContent(37, 1.);
1982 fHistAnalysisSummary->GetXaxis()->SetBinLabel(38, "fMinPvalue");
1983 fHistAnalysisSummary->SetBinContent(38, fMinPvalue);
1984 fHistAnalysisSummary->GetXaxis()->SetBinLabel(39, "fMaxPvalue");
1985 fHistAnalysisSummary->SetBinContent(39, fMaxPvalue);
1986 fHistAnalysisSummary->GetXaxis()->SetBinLabel(40, "fExcludeLeadingJetsFromFit");
1987 fHistAnalysisSummary->SetBinContent(40, fExcludeLeadingJetsFromFit);
1988 fHistAnalysisSummary->GetXaxis()->SetBinLabel(41, "fRebinSwapHistoOnTheFly");
1989 fHistAnalysisSummary->SetBinContent(41, (int)fRebinSwapHistoOnTheFly);
1990 fHistAnalysisSummary->GetXaxis()->SetBinLabel(42, "fUsePtWeight");
1991 fHistAnalysisSummary->SetBinContent(42, (int)fUsePtWeight);
1992 fHistAnalysisSummary->GetXaxis()->SetBinLabel(44, "fSoftTrackMinPt");
1993 fHistAnalysisSummary->SetBinContent(44, fSoftTrackMinPt);
1994 fHistAnalysisSummary->GetXaxis()->SetBinLabel(45, "fSoftTrackMaxPt");
1995 fHistAnalysisSummary->SetBinContent(45, fSoftTrackMaxPt);
1996 fHistAnalysisSummary->GetXaxis()->SetBinLabel(46, "fMaxCones");
1997 fHistAnalysisSummary->SetBinContent(46, fMaxCones);
1998 fHistAnalysisSummary->GetXaxis()->SetBinLabel(47, "used rho");
1999 fHistAnalysisSummary->GetXaxis()->SetBinLabel(48, "used small rho");
2001 //_____________________________________________________________________________
2002 void AliAnalysisTaskJetV2::Terminate(Option_t *)
2005 switch (fRunModeType) {
2007 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
2008 AliAnalysisTaskJetV2::Dump();
2009 for(Int_t i(0); i < fHistAnalysisSummary->GetXaxis()->GetNbins(); i++) printf( " > flag: %s \t content %.2f \n", fHistAnalysisSummary->GetXaxis()->GetBinLabel(1+i), fHistAnalysisSummary->GetBinContent(1+i));
2014 //_____________________________________________________________________________
2015 void AliAnalysisTaskJetV2::SetModulationFit(TF1* fit)
2017 // set modulation fit
2018 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
2019 if (fFitModulation) delete fFitModulation;
2020 fFitModulation = fit;
2022 //_____________________________________________________________________________
2023 void AliAnalysisTaskJetV2::SetUseControlFit(Bool_t c)
2026 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
2027 if (fFitControl) delete fFitControl;
2029 fFitControl = new TF1("controlFit", "pol0", 0, TMath::TwoPi());
2030 } else fFitControl = 0x0;
2032 //_____________________________________________________________________________
2033 TH1F* AliAnalysisTaskJetV2::GetResolutionFromOuptutFile(detectorType det, Int_t h, TArrayD* cen)
2035 // INTERFACE METHOD FOR OUTPUTFILE
2036 // get the detector resolution, user has ownership of the returned histogram
2037 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
2039 printf(" > Please add fOutputList first < \n");
2043 (cen) ? r = new TH1F("R", "R", cen->GetSize()-1, cen->GetArray()) : r = new TH1F("R", "R", 10, 0, 10);
2044 if(!cen) r->GetXaxis()->SetTitle("number of centrality bin");
2045 r->GetYaxis()->SetTitle(Form("Resolution #Psi_{%i}", h));
2046 for(Int_t i(0); i < 10; i++) {
2047 TProfile* temp((TProfile*)fOutputList->FindObject(Form("fProfV%iResolution_%i", h, i)));
2049 Double_t a(temp->GetBinContent(3)), b(temp->GetBinContent(5)), c(temp->GetBinContent(7));
2050 Double_t d(temp->GetBinContent(9)), e(temp->GetBinContent(10)), f(temp->GetBinContent(11));
2051 Double_t _a(temp->GetBinError(3)), _b(temp->GetBinError(5)), _c(temp->GetBinError(7));
2052 Double_t _d(temp->GetBinError(9)), _e(temp->GetBinError(10)), _f(temp->GetBinError(11));
2053 if(a <= 0 || b <= 0 || c <= 0 || d <= 0 || e <= 0 || f <= 0) continue;
2056 r->SetBinContent(1+i, TMath::Sqrt((a*b)/c));
2057 if(i==0) r->SetNameTitle("VZEROA resolution", "VZEROA resolution");
2058 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
2061 r->SetBinContent(1+i, TMath::Sqrt((a*c)/b));
2062 if(i==0) r->SetNameTitle("VZEROC resolution", "VZEROC resolution");
2063 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
2066 r->SetBinContent(1+i, TMath::Sqrt((b*c)/a));
2067 if(i==0) r->SetNameTitle("TPC resolution", "TPC resolution");
2068 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
2071 r->SetBinContent(1+i, TMath::Sqrt((d*e)/f));
2072 if(i==0) r->SetNameTitle("VZEROComb resolution", "VZEROComb resolution");
2073 r->SetBinError(1+i, TMath::Sqrt(_d*_d+_e*_e+_f*_f));
2080 //_____________________________________________________________________________
2081 TH1F* AliAnalysisTaskJetV2::CorrectForResolutionDiff(TH1F* v, detectorType det, TArrayD* cen, Int_t c, Int_t h)
2083 // INTERFACE METHOD FOR OUTPUT FILE
2084 // correct the supplied differential vn histogram v for detector resolution
2085 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
2086 TH1F* r(GetResolutionFromOuptutFile(det, h, cen));
2088 printf(" > Couldn't find resolution < \n");
2091 Double_t res(1./r->GetBinContent(1+r->FindBin(c)));
2092 TF1* line = new TF1("line", "pol0", 0, 200);
2093 line->SetParameter(0, res);
2097 //_____________________________________________________________________________
2098 TH1F* AliAnalysisTaskJetV2::CorrectForResolutionInt(TH1F* v, detectorType det, TArrayD* cen, Int_t h)
2100 // INTERFACE METHOD FOR OUTPUT FILE
2101 // correct the supplied intetrated vn histogram v for detector resolution
2102 // integrated vn must have the same centrality binning as the resolotion correction
2103 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
2104 TH1F* r(GetResolutionFromOuptutFile(det, h, cen));
2108 //_____________________________________________________________________________
2109 TH1F* AliAnalysisTaskJetV2::GetDifferentialQC(TProfile* refCumulants, TProfile* diffCumlants, TArrayD* ptBins, Int_t h)
2111 // get differential QC
2112 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
2113 Double_t r(refCumulants->GetBinContent(h-1)); // v2 reference flow
2114 if(r > 0) r = TMath::Sqrt(r);
2115 TH1F* qc = new TH1F(Form("QC2v%i", h), Form("QC2v%i", h), ptBins->GetSize()-1, ptBins->GetArray());
2116 Double_t a(0), b(0), c(0); // dummy variables
2117 for(Int_t i(0); i < ptBins->GetSize(); i++) {
2119 a = diffCumlants->GetBinContent(1+i);
2120 b = diffCumlants->GetBinError(1+i);
2122 qc->SetBinContent(1+i, c);
2123 (a <= 0 || b <= 0) ? qc->SetBinError(1+i, b) : qc->SetBinError(1+i, TMath::Sqrt(c*c*b*b/(a*a)));
2129 //_____________________________________________________________________________
2130 void AliAnalysisTaskJetV2::ReadVZEROCalibration2010h()
2132 // necessary for calibration of 10h vzero event plane. code copied from flow package
2133 // (duplicate, but i didn't want to introduce an ulgy dependency )
2134 // this function is only called when the runnumber changes
2135 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
2137 // 1) check if the proper chi weights for merging vzero a and vzero c ep are present
2138 // if not, use sane defaults. centrality binning is equal to that given in the fVZEROcentralityBin snippet
2140 // chi values can be calculated using the static helper function
2141 // AliAnalysisTaskJetV2::CalculateEventPlaneChi(Double_t res) where res is the event plane
2142 // resolution in a given centrality bin
2144 // the resolutions that were used for these defaults are
2145 // this might need a bit of updating as they were read 'by-eye' from a performance plot ..
2146 // Double_t R2VZEROA[] = {.35, .40, .48, .50, .48, .45, .38, .26, .16};
2147 // Double_t R2VZEROC[] = {.45, .60, .70, .73, .68, .60, .40, .36, .17};
2148 // Double_t R3VZEROA[] = {.22, .23, .22, .19, .15, .12, .08, .00, .00};
2149 // Double_t R3VZEROC[] = {.30, .30, .28, .25, .22, .17, .11, .00, .00};
2151 Double_t chiC2[] = {0.771423, 1.10236, 1.38116, 1.48077, 1.31964, 1.10236, 0.674622, 0.600403, 0.273865};
2152 Double_t chiA2[] = {0.582214, 0.674622, 0.832214, 0.873962, 0.832214, 0.771423, 0.637146, 0.424255, 0.257385};
2153 Double_t chiC3[] = {0.493347, 0.493347, 0.458557, 0.407166, 0.356628, 0.273865, 0.176208, 6.10352e-05, 6.10352e-05};
2154 Double_t chiA3[] = {0.356628, 0.373474, 0.356628, 0.306702, 0.24115, 0.192322, 0.127869, 6.10352e-05, 6.10352e-05};
2156 if(!fChi2A) fChi2A = new TArrayD(9, chiA2);
2157 if(!fChi2C) fChi2C = new TArrayD(9, chiC2);
2158 if(!fChi3A) fChi3A = new TArrayD(9, chiA3);
2159 if(!fChi3C) fChi3C = new TArrayD(9, chiC3);
2161 // 2) open database file
2162 fOADB = TFile::Open("$ALICE_ROOT/OADB/PWGCF/VZERO/VZEROcalibEP.root");
2163 if(fOADB->IsZombie()){
2164 printf("OADB file $ALICE_ROOT/OADB/PWGCF/VZERO/VZEROcalibEP.root cannot be opened, CALIBRATION FAILED !");
2168 AliOADBContainer *cont = (AliOADBContainer*) fOADB->Get("hMultV0BefCorr");
2170 // see if database is readable
2171 printf("OADB object hMultV0BefCorr is not available in the file\n");
2174 Int_t run(fRunNumber);
2175 if(!(cont->GetObject(run))){
2176 // if the run isn't recognized fall back to a default run
2177 printf("OADB object hMultVZEROBefCorr is not available for run %i (used default run 137366)\n",run);
2180 // step 3) get the proper multiplicity weights from the vzero signal
2181 fVZEROgainEqualization = ((TH2F*)cont->GetObject(run))->ProfileX();
2182 if(!fVZEROgainEqualization) {
2183 AliFatal(Form("%s: Fatal error, couldn't read fVZEROgainEqualization from OADB object < \n", GetName()));
2187 TF1* fpol0 = new TF1("fpol0","pol0");
2188 if(fVZEROgainEqualizationPerRing) {
2189 // do the calibration per ring
2190 // start with the vzero c rings (segments 0 through 31)
2191 fVZEROgainEqualization->Fit(fpol0, "", "", 0, 8);
2192 (fUseVZERORing[0]) ? SetVZEROCpol(0, fpol0->GetParameter(0)) : SetVZEROCpol(0, 0.);
2193 fVZEROgainEqualization->Fit(fpol0, "", "", 8, 16);
2194 (fUseVZERORing[1]) ? SetVZEROCpol(1, fpol0->GetParameter(0)) : SetVZEROCpol(1, 0.);
2195 fVZEROgainEqualization->Fit(fpol0, "", "", 16, 24);
2196 (fUseVZERORing[2]) ? SetVZEROCpol(2, fpol0->GetParameter(0)) : SetVZEROCpol(2, 0.);
2197 fVZEROgainEqualization->Fit(fpol0, "", "", 24, 32);
2198 (fUseVZERORing[3]) ? SetVZEROCpol(3, fpol0->GetParameter(0)) : SetVZEROCpol(3, 0.);
2199 // same thing for vero A
2200 fVZEROgainEqualization->Fit(fpol0, "", "", 32, 40);
2201 (fUseVZERORing[4]) ? SetVZEROApol(0, fpol0->GetParameter(0)) : SetVZEROApol(0, 0.);
2202 fVZEROgainEqualization->Fit(fpol0, "", "", 40, 48);
2203 (fUseVZERORing[5]) ? SetVZEROApol(1, fpol0->GetParameter(0)) : SetVZEROApol(1, 0.);
2204 fVZEROgainEqualization->Fit(fpol0, "", "", 48, 56);
2205 (fUseVZERORing[6]) ? SetVZEROApol(2, fpol0->GetParameter(0)) : SetVZEROApol(2, 0.);
2206 fVZEROgainEqualization->Fit(fpol0, "", "", 56, 64);
2207 (fUseVZERORing[7]) ? SetVZEROApol(3, fpol0->GetParameter(0)) : SetVZEROApol(3, 0.);
2209 // do the calibration in one go. the calibration will still be
2210 // stored per ring, but each ring has the same weight now
2211 // this should be the default for the analysis as the database is tuned to this configuration
2212 fVZEROgainEqualization->Fit(fpol0,"","",0,31);
2213 for(Int_t i(0); i < 4; i++) SetVZEROCpol(i, fpol0->GetParameter(0));
2214 fVZEROgainEqualization->Fit(fpol0,"","",32,64);
2215 for(Int_t i(0); i < 4; i++) SetVZEROApol(i, fpol0->GetParameter(0));
2218 // step 4) extract the information to re-weight the q-vectors
2219 for(Int_t iside=0;iside<2;iside++){
2220 for(Int_t icoord=0;icoord<2;icoord++){
2221 for(Int_t i=0;i < 9;i++){
2223 if(iside==0 && icoord==0)
2224 snprintf(namecont,100,"hQxc2_%i",i);
2225 else if(iside==1 && icoord==0)
2226 snprintf(namecont,100,"hQxa2_%i",i);
2227 else if(iside==0 && icoord==1)
2228 snprintf(namecont,100,"hQyc2_%i",i);
2229 else if(iside==1 && icoord==1)
2230 snprintf(namecont,100,"hQya2_%i",i);
2232 cont = (AliOADBContainer*) fOADB->Get(namecont);
2234 printf("OADB object %s is not available in the file\n",namecont);
2238 if(!(cont->GetObject(run))){
2239 printf("OADB object %s is not available for run %i (used run 137366)\n",namecont,run);
2243 // store info for all centralities to cache
2244 fMeanQ[i][iside][icoord] = ((TH1F *) cont->GetObject(run))->GetMean();
2245 fWidthQ[i][iside][icoord] = ((TH1F *) cont->GetObject(run))->GetRMS();
2248 if(iside==0 && icoord==0)
2249 snprintf(namecont,100,"hQxc3_%i",i);
2250 else if(iside==1 && icoord==0)
2251 snprintf(namecont,100,"hQxa3_%i",i);
2252 else if(iside==0 && icoord==1)
2253 snprintf(namecont,100,"hQyc3_%i",i);
2254 else if(iside==1 && icoord==1)
2255 snprintf(namecont,100,"hQya3_%i",i);
2257 cont = (AliOADBContainer*) fOADB->Get(namecont);
2259 printf("OADB object %s is not available in the file\n",namecont);
2263 if(!(cont->GetObject(run))){
2264 printf("OADB object %s is not available for run %i (used run 137366)\n",namecont,run);
2267 // store info for all centralities to cache
2268 fMeanQv3[i][iside][icoord] = ((TH1F *) cont->GetObject(run))->GetMean();
2269 fWidthQv3[i][iside][icoord] = ((TH1F *) cont->GetObject(run))->GetRMS();
2273 // cleanup. the opened file is closed in the destructor, otherwise fVZEROgainEqualization is no longer available
2276 //_____________________________________________________________________________
2277 Int_t AliAnalysisTaskJetV2::GetVZEROCentralityBin() const
2279 // return cache index number corresponding to the event centrality
2280 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
2281 Float_t v0Centr(InputEvent()->GetCentrality()->GetCentralityPercentile("V0M"));
2282 if(v0Centr < 5) return 0;
2283 else if(v0Centr < 10) return 1;
2284 else if(v0Centr < 20) return 2;
2285 else if(v0Centr < 30) return 3;
2286 else if(v0Centr < 40) return 4;
2287 else if(v0Centr < 50) return 5;
2288 else if(v0Centr < 60) return 6;
2289 else if(v0Centr < 70) return 7;
2292 //_____________________________________________________________________________