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
43 #include <AliAnalysisTask.h>
44 #include <AliAnalysisManager.h>
45 #include <AliCentrality.h>
46 #include <AliVVertex.h>
47 #include <AliVTrack.h>
48 #include <AliESDEvent.h>
49 #include <AliAODEvent.h>
50 #include <AliAODTrack.h>
51 // emcal jet framework includes
52 #include <AliPicoTrack.h>
53 #include <AliEmcalJet.h>
54 #include <AliRhoParameter.h>
55 #include <AliLocalRhoParameter.h>
56 #include <AliAnalysisTaskJetV2.h>
57 #include <AliClusterContainer.h>
59 class AliAnalysisTaskJetV2;
62 ClassImp(AliAnalysisTaskJetV2)
64 AliAnalysisTaskJetV2::AliAnalysisTaskJetV2() : AliAnalysisTaskEmcalJet("AliAnalysisTaskJetV2", kTRUE),
65 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(kTPC), 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.), fUseV0EventPlaneFromHeader(kTRUE), /*fExplicitOutlierCut(-1),*/ 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), fHistPsiTPCiV0M(0), fHistPsiVZEROATRK(0), fHistPsiVZEROCTRK(0), fHistPsiVZEROTRK(0), fHistPsiTPCTRK(0), fHistRhoVsMult(0), fHistRhoVsCent(0), fHistRhoAVsMult(0), fHistRhoAVsCent(0) {
66 for(Int_t i(0); i < 10; i++) {
67 fProfV2Resolution[i] = 0;
68 fProfV3Resolution[i] = 0;
69 fHistPicoTrackPt[i] = 0;
70 fHistPicoTrackMult[i] = 0;
74 fHistClusterPt[i] = 0;
75 fHistClusterEtaPhi[i] = 0;
76 fHistClusterEtaPhiWeighted[i] = 0;
77 fHistRhoPackage[i] = 0;
80 fHistRhoVsRCPt[i] = 0;
82 fHistDeltaPtDeltaPhi2[i] = 0;
83 fHistDeltaPtDeltaPhi2Rho0[i] = 0;
84 fHistRCPhiEtaExLJ[i] = 0;
85 fHistRhoVsRCPtExLJ[i] = 0;
87 fHistDeltaPtDeltaPhi2ExLJ[i] = 0;
88 fHistDeltaPtDeltaPhi2ExLJRho0[i] = 0;
91 fHistJetEtaPhi[i] = 0;
92 fHistJetPtArea[i] = 0;
94 fHistJetPtConstituents[i] = 0;
95 fHistJetEtaRho[i] = 0;
96 fHistJetPsi2Pt[i] = 0;
97 fHistJetPsi2PtRho0[i] = 0;
99 // default constructor
101 //_____________________________________________________________________________
102 AliAnalysisTaskJetV2::AliAnalysisTaskJetV2(const char* name, runModeType type) : AliAnalysisTaskEmcalJet(name, kTRUE),
103 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(kTPC), 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.), fUseV0EventPlaneFromHeader(kTRUE), /*fExplicitOutlierCut(-1),*/ 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), fHistPsiTPCiV0M(0), fHistPsiVZEROATRK(0), fHistPsiVZEROCTRK(0), fHistPsiVZEROTRK(0), fHistPsiTPCTRK(0), fHistRhoVsMult(0), fHistRhoVsCent(0), fHistRhoAVsMult(0), fHistRhoAVsCent(0) {
104 for(Int_t i(0); i < 10; i++) {
105 fProfV2Resolution[i] = 0;
106 fProfV3Resolution[i] = 0;
107 fHistPicoTrackPt[i] = 0;
108 fHistPicoTrackMult[i] = 0;
109 fHistPicoCat1[i] = 0;
110 fHistPicoCat2[i] = 0;
111 fHistPicoCat3[i] = 0;
112 fHistClusterPt[i] = 0;
113 fHistClusterEtaPhi[i] = 0;
114 fHistClusterEtaPhiWeighted[i] = 0;
115 fHistRhoPackage[i] = 0;
117 fHistRCPhiEta[i] = 0;
118 fHistRhoVsRCPt[i] = 0;
120 fHistDeltaPtDeltaPhi2[i] = 0;
121 fHistDeltaPtDeltaPhi2Rho0[i] = 0;
122 fHistRCPhiEtaExLJ[i] = 0;
123 fHistRhoVsRCPtExLJ[i] = 0;
124 fHistRCPtExLJ[i] = 0;
125 fHistDeltaPtDeltaPhi2ExLJ[i] = 0;
126 fHistDeltaPtDeltaPhi2ExLJRho0[i] = 0;
127 fHistJetPtRaw[i] = 0;
129 fHistJetEtaPhi[i] = 0;
130 fHistJetPtArea[i] = 0;
131 fHistJetPtEta[i] = 0;
132 fHistJetPtConstituents[i] = 0;
133 fHistJetEtaRho[i] = 0;
134 fHistJetPsi2Pt[i] = 0;
135 fHistJetPsi2PtRho0[i] = 0;
138 DefineInput(0, TChain::Class());
139 DefineOutput(1, TList::Class());
140 switch (fRunModeType) {
142 gStyle->SetOptFit(1);
143 DefineOutput(2, TList::Class());
144 DefineOutput(3, TList::Class());
146 default: fDebug = -1; // suppress debug info explicitely when not running locally
148 switch (fCollisionType) {
150 fFitModulationType = kNoFit;
154 if(fLocalRhoName=="") fLocalRhoName = Form("LocalRhoFrom_%s", GetName());
156 //_____________________________________________________________________________
157 AliAnalysisTaskJetV2::~AliAnalysisTaskJetV2()
160 if(fOutputList) {delete fOutputList; fOutputList = 0x0;}
161 if(fOutputListGood) {delete fOutputListGood; fOutputListGood = 0x0;}
162 if(fOutputListBad) {delete fOutputListBad; fOutputListBad = 0x0;}
163 if(fFitModulation) {delete fFitModulation; fFitModulation = 0x0;}
164 if(fHistSwap) {delete fHistSwap; fHistSwap = 0x0;}
165 if(fCentralityClasses) {delete fCentralityClasses; fCentralityClasses = 0x0;}
166 if(fExpectedRuns) {delete fExpectedRuns; fExpectedRuns = 0x0;}
167 if(fExpectedSemiGoodRuns) {delete fExpectedSemiGoodRuns; fExpectedSemiGoodRuns = 0x0;}
168 if(fFitControl) {delete fFitControl; fFitControl = 0x0;}
170 //_____________________________________________________________________________
171 void AliAnalysisTaskJetV2::ExecOnce()
174 fLocalRho = new AliLocalRhoParameter(fLocalRhoName.Data(), 0);
176 if(!(InputEvent()->FindListObject(fLocalRho->GetName()))) {
177 InputEvent()->AddObject(fLocalRho);
179 AliFatal(Form("%s: Container with name %s already present. Aborting", GetName(), fLocalRho->GetName()));
182 AliAnalysisTaskEmcalJet::ExecOnce(); // init the base class
183 AliAnalysisTaskEmcalJet::SetVzRange(-1.*fAbsVertexZ, fAbsVertexZ);
184 if(!GetJetContainer()) AliFatal(Form("%s: Couldn't find jet container. Aborting !", GetName()));
186 //_____________________________________________________________________________
187 Bool_t AliAnalysisTaskJetV2::Notify()
189 // determine the run number to see if the track and jet cuts should be refreshed for semi-good TPC runs
190 if(fRunNumber != InputEvent()->GetRunNumber()) {
191 fRunNumber = InputEvent()->GetRunNumber(); // set the current run number
192 if(fDebug > 0) printf("__FUNC__ %s > NEW RUNNUMBER DETECTED \n ", __func__);
193 // check if this is 10h or 11h data. for 10h we don't want to change the acceptance
194 switch (fCollisionType) {
200 // reset the cuts. should be a pointless operation except for the case where the run number changes
201 // from semi-good back to good on one node, which is not a likely scenario (unless trains will
202 // run as one masterjob)
203 AliAnalysisTaskEmcal::SetTrackPhiLimits(-10., 10.);
204 switch (fAnalysisType) {
206 AliAnalysisTaskEmcalJet::SetJetPhiLimits(-10., 10.);
209 AliAnalysisTaskEmcalJet::SetJetPhiLimits(1.405 + GetJetRadius(), 3.135 - GetJetRadius());
213 if(fCachedRho) { // if there's a cached rho, it's the default, so switch back
214 if(fDebug > 0) printf("__FUNC__ %s > replacing rho with cached rho \n ", __func__);
215 fRho = fCachedRho; // reset rho back to cached value. again, should be pointless
217 Bool_t flaggedAsSemiGood(kFALSE); // not flagged as anything
218 for(Int_t i(0); i < fExpectedSemiGoodRuns->GetSize(); i++) {
219 if(fExpectedSemiGoodRuns->At(i) == fRunNumber) { // run is semi-good
220 if(fDebug > 0) printf("__FUNC__ %s > semi-good tpc run detected, adjusting acceptance \n ", __func__);
221 flaggedAsSemiGood = kTRUE;
222 switch (fAnalysisType) {
223 // for full jets the jet acceptance does not have to be changed as emcal does not
224 // cover the tpc low voltage readout strips
226 AliAnalysisTaskEmcalJet::SetJetPhiLimits(fSemiGoodJetMinPhi, fSemiGoodJetMaxPhi); // just an acceptance cut, jets are obtained from full azimuth, so no edge effects
230 AliAnalysisTaskEmcal::SetTrackPhiLimits(fSemiGoodTrackMinPhi, fSemiGoodTrackMaxPhi); // only affects vn extraction, NOT jet finding
231 // for semi-good runs, also try to get the 'small rho' estimate, if it is available
232 AliRhoParameter* tempRho(dynamic_cast<AliRhoParameter*>(InputEvent()->FindListObject(fNameSmallRho.Data())));
234 if(fDebug > 0) printf("__FUNC__ %s > switching to small rho, caching normal rho \n ", __func__);
235 fHistAnalysisSummary->SetBinContent(54, 1.); // bookkeep the fact that small rho is used
236 fCachedRho = fRho; // cache the original rho ...
237 fRho = tempRho; // ... and use the small rho
241 if(!flaggedAsSemiGood) {
242 // in case the run is not a semi-good run, check if it is recognized as another run
243 // only done to catch unexpected runs
244 for(Int_t i(0); i < fExpectedRuns->GetSize(); i++) {
245 if(fExpectedRuns->At(i) == fRunNumber) break; // run is known, break the loop else store the number in a random bin
246 fHistUndeterminedRunQA->SetBinContent(TMath::Nint(10.*gRandom->Uniform(0.,.9))+1, fRunNumber);
248 fHistAnalysisSummary->SetBinContent(53, 1.); // bookkeep which rho estimate is used
253 //_____________________________________________________________________________
254 Bool_t AliAnalysisTaskJetV2::InitializeAnalysis()
256 // initialize the anaysis
257 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
258 // if not set, estimate the number of cones that would fit into the selected acceptance
259 if(fMaxCones <= 0) fMaxCones = TMath::CeilNint((TMath::Abs(GetJetContainer()->GetJetEtaMax()-GetJetContainer()->GetJetEtaMin())*TMath::Abs(GetJetContainer()->GetJetPhiMax()-GetJetContainer()->GetJetPhiMin()))/(TMath::Pi()*GetJetRadius()*GetJetRadius()));
260 // manually 'override' the default acceptance cuts of the emcal framework (use with caution)
261 if(fMinDisanceRCtoLJ==0) fMinDisanceRCtoLJ = GetJetRadius();
262 if(dynamic_cast<AliAODEvent*>(InputEvent())) fDataType = kAOD; // determine the datatype
263 else if(dynamic_cast<AliESDEvent*>(InputEvent())) fDataType = kESD;
264 fHistAnalysisSummary->SetBinContent(36, (int)fDataType);
265 if(!fRandom) fRandom = new TRandom3(0); // set randomizer and random seed
266 switch (fFitModulationType) {
267 case kNoFit : { SetModulationFit(new TF1("fix_kNoFit", "[0]", 0, TMath::TwoPi())); } break;
269 SetModulationFit(new TF1("fit_kV2", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
270 fFitModulation->SetParameter(0, 0.); // normalization
271 fFitModulation->SetParameter(3, 0.2); // v2
272 fFitModulation->FixParameter(1, 1.); // constant
273 fFitModulation->FixParameter(2, 2.); // constant
276 SetModulationFit(new TF1("fit_kV3", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
277 fFitModulation->SetParameter(0, 0.); // normalization
278 fFitModulation->SetParameter(3, 0.2); // v3
279 fFitModulation->FixParameter(1, 1.); // constant
280 fFitModulation->FixParameter(2, 3.); // constant
282 default : { // for the combined fit, the 'direct fourier series' or the user supplied vn values we use v2 and v3
283 SetModulationFit(new TF1("fit_kCombined", "[0]*([1]+[2]*([3]*TMath::Cos([2]*(x-[4]))+[7]*TMath::Cos([5]*(x-[6]))))", 0, TMath::TwoPi()));
284 fFitModulation->SetParameter(0, 0.); // normalization
285 fFitModulation->SetParameter(3, 0.2); // v2
286 fFitModulation->FixParameter(1, 1.); // constant
287 fFitModulation->FixParameter(2, 2.); // constant
288 fFitModulation->FixParameter(5, 3.); // constant
289 fFitModulation->SetParameter(7, 0.2); // v3
292 switch (fRunModeType) {
293 case kGrid : { fFitModulationOptions += "N0"; } break;
296 FillAnalysisSummaryHistogram();
299 //_____________________________________________________________________________
300 TH1F* AliAnalysisTaskJetV2::BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Int_t c, Bool_t append)
302 // book a TH1F and connect it to the output container
303 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
304 if(fReduceBinsXByFactor > 0 ) bins = TMath::Nint(bins/fReduceBinsXByFactor);
305 if(!fOutputList) return 0x0;
307 if(c!=-1) { // format centrality dependent histograms accordingly
308 name = Form("%s_%i", name, c);
309 title += Form("_%i-%i", (int)(fCentralityClasses->At(c)), (int)(fCentralityClasses->At((1+c))));
311 title += Form(";%s;[counts]", x);
312 TH1F* histogram = new TH1F(name, title.Data(), bins, min, max);
314 if(append) fOutputList->Add(histogram);
317 //_____________________________________________________________________________
318 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)
320 // book a TH2F and connect it to the output container
321 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
322 if(fReduceBinsXByFactor > 0 ) binsx = TMath::Nint(binsx/fReduceBinsXByFactor);
323 if(fReduceBinsYByFactor > 0 ) binsy = TMath::Nint(binsy/fReduceBinsYByFactor);
324 if(!fOutputList) return 0x0;
326 if(c!=-1) { // format centrality dependent histograms accordingly
327 name = Form("%s_%i", name, c);
328 title += Form("_%i-%i", (int)fCentralityClasses->At(c), (int)(fCentralityClasses->At((1+c))));
330 title += Form(";%s;%s", x, y);
331 TH2F* histogram = new TH2F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy);
333 if(append) fOutputList->Add(histogram);
336 //_____________________________________________________________________________
337 void AliAnalysisTaskJetV2::UserCreateOutputObjects()
339 // create output objects. also initializes some default values in case they aren't
340 // loaded via the AddTask macro
341 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
342 fOutputList = new TList();
343 fOutputList->SetOwner(kTRUE);
344 if(!fCentralityClasses) { // classes must be defined at this point
345 Double_t c[] = {0., 20., 40., 60., 80., 100.};
346 fCentralityClasses = new TArrayD(sizeof(c)/sizeof(c[0]), c);
348 if(!fExpectedRuns) { // expected runs must be defined at this point
349 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 */};
350 fExpectedRuns = new TArrayI(sizeof(r)/sizeof(r[0]), r);
352 // set default semi-good runs only for 11h data
353 switch (fCollisionType) {
354 case kPbPb10h : break;
356 if(!fExpectedSemiGoodRuns) {
357 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};
358 fExpectedSemiGoodRuns = new TArrayI(sizeof(r)/sizeof(r[0]), r);
364 fHistCentrality = BookTH1F("fHistCentrality", "centrality", 102, -2, 100);
365 fHistVertexz = BookTH1F("fHistVertexz", "vertex z (cm)", 100, -12, 12);
367 // pico track and emcal cluster kinematics
368 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
369 fHistPicoTrackPt[i] = BookTH1F("fHistPicoTrackPt", "p_{t} [GeV/c]", 100, 0, 100, i);
370 fHistPicoTrackMult[i] = BookTH1F("fHistPicoTrackMult", "multiplicity", 100, 0, 5000, i);
371 if(fFillQAHistograms) {
372 fHistPicoCat1[i] = BookTH2F("fHistPicoCat1", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
373 fHistPicoCat2[i] = BookTH2F("fHistPicoCat2", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
374 fHistPicoCat3[i] = BookTH2F("fHistPicoCat3", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
375 if(fAnalysisType == AliAnalysisTaskJetV2::kFull) {
376 fHistClusterPt[i] = BookTH1F("fHistClusterPt", "p_{t} [GeV/c]", 100, 0, 100, i);
377 fHistClusterEtaPhi[i] = BookTH2F("fHistClusterEtaPhi", "#eta", "#phi", 100, -1., 1., 100, 0, TMath::TwoPi(), i);
378 fHistClusterEtaPhiWeighted[i] = BookTH2F("fHistClusterEtaPhiWeighted", "#eta", "#phi", 100, -1., 1., 100, 0, TMath::TwoPi(), i);
383 if(fFillQAHistograms) {
384 // event plane estimates and quality
385 fHistPsiControl = new TProfile("fHistPsiControl", "fHistPsiControl", 10, 0, 10);
386 fHistPsiControl->Sumw2();
387 fHistPsiSpread = new TProfile("fHistPsiSpread", "fHistPsiSpread", 4, 0, 4);
388 fHistPsiSpread->Sumw2();
389 fHistPsiControl->GetXaxis()->SetBinLabel(1, "<#Psi_{2, VZEROA}>");
390 fHistPsiControl->GetXaxis()->SetBinLabel(2, "<#Psi_{2, VZEROC}>");
391 fHistPsiControl->GetXaxis()->SetBinLabel(3, "<#Psi_{2, TPC}>");
392 fHistPsiControl->GetXaxis()->SetBinLabel(4, "<#Psi_{2, TPC, #eta < 0}>");
393 fHistPsiControl->GetXaxis()->SetBinLabel(5, "<#Psi_{2, TPC, #eta > 0}>");
394 fHistPsiControl->GetXaxis()->SetBinLabel(6, "<#Psi_{3, VZEROA}>");
395 fHistPsiControl->GetXaxis()->SetBinLabel(7, "<#Psi_{3, VZEROC}>");
396 fHistPsiControl->GetXaxis()->SetBinLabel(8, "<#Psi_{3, TPC}>");
397 fHistPsiControl->GetXaxis()->SetBinLabel(9, "<#Psi_{3, TPC, #eta < 0}>");
398 fHistPsiControl->GetXaxis()->SetBinLabel(10, "<#Psi_{3, TPC, #eta > 0}>");
399 fHistPsiSpread->GetXaxis()->SetBinLabel(1, "<#Psi_{2, VZEROA} - #Psi_{2, VZEROC}>");
400 fHistPsiSpread->GetXaxis()->SetBinLabel(2, "<#Psi_{2, VZEROC} - #Psi_{2, TPC}>");
401 fHistPsiSpread->GetXaxis()->SetBinLabel(3, "<#Psi_{2, VZEROC} - #Psi_{2, TPC}>");
402 fHistPsiSpread->GetXaxis()->SetBinLabel(4, "<#Psi_{2, TPC, #eta < 0} - #Psi_{2, TPC, #eta > 0}>");
403 fOutputList->Add(fHistPsiControl);
404 fOutputList->Add(fHistPsiSpread);
405 fHistPsiVZEROA = BookTH1F("fHistPsiVZEROA", "#Psi_{VZEROA}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
406 fHistPsiVZEROC = BookTH1F("fHistPsiVZEROC", "#Psi_{VZEROC}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
407 fHistPsiVZERO = BookTH1F("fHistPsiVZERO", "#Psi_{VZERO}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
408 fHistPsiTPC = BookTH1F("fHistPsiTPC", "#Psi_{TPC}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
409 fHistPsiVZEROAV0M = BookTH2F("fHistPsiVZEROAV0M", "V0M", "#Psi_{2, VZEROA}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
410 fHistPsiVZEROCV0M = BookTH2F("fHistPsiVZEROCV0M", "V0M", "#Psi_{2, VZEROC}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
411 fHistPsiVZEROVV0M = BookTH2F("fHistPsiVZEROV0M", "V0M", "#Psi_{2, VZERO}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
412 fHistPsiTPCiV0M = BookTH2F("fHistPsiTPCV0M", "V0M", "#Psi_{2, TRK}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
413 fHistPsiVZEROATRK = BookTH2F("fHistPsiVZEROATRK", "TRK", "#Psi_{2, VZEROA}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
414 fHistPsiVZEROCTRK = BookTH2F("fHistPsiVZEROCTRK", "TRK", "#Psi_{2, VZEROC}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
415 fHistPsiVZEROTRK = BookTH2F("fHistPsiVZEROTRK", "TRK", "#Psi_{2, VZERO}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
416 fHistPsiTPCTRK = BookTH2F("fHistPsiTPCTRK", "TRK", "#Psi_{2, TRK}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
419 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i ++) {
420 fHistRhoPackage[i] = BookTH1F("fHistRhoPackage", "#rho [GeV/c]", 100, 0, 150, i);
421 fHistRho[i] = BookTH1F("fHistRho", "#rho [GeV/c]", 100, 0, 150, i);
423 fHistRhoVsMult = BookTH2F("fHistRhoVsMult", "multiplicity", "#rho [GeV/c]", 100, 0, 4000, 100, 0, 250);
424 fHistRhoVsCent = BookTH2F("fHistRhoVsCent", "centrality", "#rho [GeV/c]", 100, 0, 100, 100, 0, 250);
425 fHistRhoAVsMult = BookTH2F("fHistRhoAVsMult", "multiplicity", "#rho * A (jet) [GeV/c]", 100, 0, 4000, 100, 0, 50);
426 fHistRhoAVsCent = BookTH2F("fHistRhoAVsCent", "centrality", "#rho * A (jet) [GeV/c]", 100, 0, 100, 100, 0, 50);
428 TString detector("");
429 switch (fDetectorType) {
430 case kTPC : detector+="TPC";
432 case kVZEROA : detector+="VZEROA";
434 case kVZEROC : detector+="VZEROC";
436 case kVZEROComb : detector+="VZEROComb";
440 // delta pt distributions
441 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i ++) {
442 if(fFillQAHistograms) fHistRCPhiEta[i] = BookTH2F("fHistRCPhiEta", "#phi (RC)", "#eta (RC)", 40, 0, TMath::TwoPi(), 40, -1, 1, i);
443 fHistRhoVsRCPt[i] = BookTH2F("fHistRhoVsRCPt", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i);
444 fHistRCPt[i] = BookTH1F("fHistRCPt", "p_{t} (RC) [GeV/c]", 130, -20, 150, i);
445 if(fFillQAHistograms) fHistRCPhiEtaExLJ[i] = BookTH2F("fHistRCPhiEtaExLJ", "#phi (RC)", "#eta (RC)", 40, 0, TMath::TwoPi(), 40, -1, 1, i);
446 fHistDeltaPtDeltaPhi2[i] = BookTH2F("fHistDeltaPtDeltaPhi2", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
447 fHistDeltaPtDeltaPhi2Rho0[i] = BookTH2F("fHistDeltaPtDeltaPhi2Rho0", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
448 fHistRhoVsRCPtExLJ[i] = BookTH2F("fHistRhoVsRCPtExLJ", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i);
449 fHistRCPtExLJ[i] = BookTH1F("fHistRCPtExLJ", "p_{t} (RC) [GeV/c]", 130, -20, 150, i);
450 fHistDeltaPtDeltaPhi2ExLJ[i] = BookTH2F("fHistDeltaPtDeltaPhi2ExLJ", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
451 fHistDeltaPtDeltaPhi2ExLJRho0[i] = BookTH2F("fHistDeltaPtDeltaPhi2ExLJRho0", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
452 // jet histograms (after kinematic cuts)
453 fHistJetPtRaw[i] = BookTH1F("fHistJetPtRaw", "p_{t, jet} RAW [GeV/c]", 200, -50, 150, i);
454 fHistJetPt[i] = BookTH1F("fHistJetPt", "p_{t, jet} [GeV/c]", 350, -100, 250, i);
455 if(fFillQAHistograms) fHistJetEtaPhi[i] = BookTH2F("fHistJetEtaPhi", "#eta", "#phi", 100, -1, 1, 100, 0, TMath::TwoPi(), i);
456 fHistJetPtArea[i] = BookTH2F("fHistJetPtArea", "p_{t, jet} [GeV/c]", "Area", 175, -100, 250, 30, 0, 0.9, i);
457 fHistJetPtEta[i] = BookTH2F("fHistJetPtEta", "p_{t, jet} [GeV/c]", "Eta", 175, -100, 250, 30, -0.9, 0.9, i);
458 fHistJetPtConstituents[i] = BookTH2F("fHistJetPtConstituents", "p_{t, jet} [GeV/c]", "Area", 350, -100, 250, 60, 0, 150, i);
459 fHistJetEtaRho[i] = BookTH2F("fHistJetEtaRho", "#eta", "#rho", 100, -1, 1, 100, 0, 300, i);
460 // in plane and out of plane spectra
461 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);
462 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);
463 // profiles for all correlator permutations which are necessary to calculate each second and third order event plane resolution
464 fProfV2Resolution[i] = new TProfile(Form("fProfV2Resolution_%i", i), Form("fProfV2Resolution_%i", i), 11, -0.5, 10.5);
465 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(3, "<cos(2(#Psi_{VZEROA} - #Psi_{VZEROC}))>");
466 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(4, "<cos(2(#Psi_{VZEROC} - #Psi_{VZEROA}))>");
467 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(5, "<cos(2(#Psi_{VZEROA} - #Psi_{TPC}))>");
468 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(6, "<cos(2(#Psi_{TPC} - #Psi_{VZEROA}))>");
469 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(7, "<cos(2(#Psi_{VZEROC} - #Psi_{TPC}))>");
470 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(8, "<cos(2(#Psi_{TPC} - #Psi_{VZEROC}))>");
471 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(9, "<cos(2(#Psi_{VZERO} - #Psi_{TPC_A}))>");
472 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(10, "<cos(2(#Psi_{VZERO} - #Psi_{TPC_B}))>");
473 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(11, "<cos(2(#Psi_{TPC_A} - #Psi_{TPC_B}))>");
474 fOutputList->Add(fProfV2Resolution[i]);
475 fProfV3Resolution[i] = new TProfile(Form("fProfV3Resolution_%i", i), Form("fProfV3Resolution_%i", i), 11, -0.5, 10.5);
476 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(3, "<cos(3(#Psi_{VZEROA} - #Psi_{VZEROC}))>");
477 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(4, "<cos(3(#Psi_{VZEROC} - #Psi_{VZEROA}))>");
478 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(5, "<cos(3(#Psi_{VZEROA} - #Psi_{TPC}))>");
479 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(6, "<cos(3(#Psi_{TPC} - #Psi_{VZEROA}))>");
480 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(7, "<cos(3(#Psi_{VZEROC} - #Psi_{TPC}))>");
481 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(8, "<cos(3(#Psi_{TPC} - #Psi_{VZEROC}))>");
482 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(9, "<cos(3(#Psi_{VZERO} - #Psi_{TPC_A}))>");
483 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(10, "<cos(3(#Psi_{VZERO} - #Psi_{TPC_B}))>");
484 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(11, "<cos(3(#Psi_{TPC_A} - #Psi_{TPC_B}))>");
485 fOutputList->Add(fProfV3Resolution[i]);
488 Float_t temp[fCentralityClasses->GetSize()];
489 for(Int_t i(0); i < fCentralityClasses->GetSize(); i++) temp[i] = fCentralityClasses->At(i);
490 fProfV2 = new TProfile("fProfV2", "fProfV2", fCentralityClasses->GetSize()-1, temp);
491 fProfV3 = new TProfile("fProfV3", "fProfV3", fCentralityClasses->GetSize()-1, temp);
492 fOutputList->Add(fProfV2);
493 fOutputList->Add(fProfV3);
494 switch (fFitModulationType) {
496 fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
497 fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
498 fOutputList->Add(fProfV2Cumulant);
499 fOutputList->Add(fProfV3Cumulant);
502 fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
503 fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
504 fOutputList->Add(fProfV2Cumulant);
505 fOutputList->Add(fProfV3Cumulant);
509 // for the histograms initialized below, binning is fixed to runnumbers or flags
510 fReduceBinsXByFactor = 1;
511 fReduceBinsYByFactor = 1;
512 if(fFillQAHistograms) {
513 fHistRunnumbersEta = new TH2F("fHistRunnumbersEta", "fHistRunnumbersEta", fExpectedRuns->GetSize()+1, -.5, fExpectedRuns->GetSize()+.5, 100, -1.1, 1.1);
514 fHistRunnumbersEta->Sumw2();
515 fOutputList->Add(fHistRunnumbersEta);
516 fHistRunnumbersPhi = new TH2F("fHistRunnumbersPhi", "fHistRunnumbersPhi", fExpectedRuns->GetSize()+1, -.5, fExpectedRuns->GetSize()+.5, 100, -0.2, TMath::TwoPi()+0.2);
517 fHistRunnumbersPhi->Sumw2();
518 fOutputList->Add(fHistRunnumbersPhi);
519 for(Int_t i(0); i < fExpectedRuns->GetSize(); i++) {
520 fHistRunnumbersPhi->GetXaxis()->SetBinLabel(i+1, Form("%i", fExpectedRuns->At(i)));
521 fHistRunnumbersEta->GetXaxis()->SetBinLabel(i+1, Form("%i", fExpectedRuns->At(i)));
523 fHistRunnumbersPhi->GetXaxis()->SetBinLabel(fExpectedRuns->GetSize()+1, "undetermined");
524 fHistRunnumbersEta->GetXaxis()->SetBinLabel(fExpectedRuns->GetSize()+1, "undetermined");
526 fHistAnalysisSummary = BookTH1F("fHistAnalysisSummary", "flag", 54, -0.5, 54.5);
527 fHistSwap = new TH1F("fHistSwap", "fHistSwap", 20, 0, TMath::TwoPi());
528 if(fUsePtWeight) fHistSwap->Sumw2();
530 if(fUserSuppliedV2) fOutputList->Add(fUserSuppliedV2);
531 if(fUserSuppliedV3) fOutputList->Add(fUserSuppliedV3);
532 if(fUserSuppliedR2) fOutputList->Add(fUserSuppliedR2);
533 if(fUserSuppliedR3) fOutputList->Add(fUserSuppliedR3);
534 // increase readability of output list
536 // cdf and pdf of chisquare distribution
537 fHistPvalueCDF = BookTH1F("fHistPvalueCDF", "CDF #chi^{2}", 50, 0, 1);
538 fHistPvalueCDFCent = BookTH2F("fHistPvalueCDFCent", "centrality", "p-value", 40, 0, 100, 40, 0, 1);
539 fHistChi2Cent = BookTH2F("fHistChi2Cent", "centrality", "#tilde{#chi^{2}}", 100, 0, 100, 100, 0, 5);
540 fHistPChi2 = BookTH2F("fHistPChi2", "p-value", "#tilde{#chi^{2}}", 1000, 0, 1, 100, 0, 5);
541 fHistKolmogorovTest = BookTH1F("fHistKolmogorovTest", "KolmogorovTest", 50, 0, 1);
542 fHistKolmogorovTestCent = BookTH2F("fHistKolmogorovTestCent", "centrality", "Kolmogorov p", 40, 0, 100, 45, 0, 1);
543 fHistPvalueCDFROOT = BookTH1F("fHistPvalueCDFROOT", "CDF #chi^{2} ROOT", 50, 0, 1);
544 fHistPvalueCDFROOTCent = BookTH2F("fHistPvalueCDFROOTCent", "centrality", "p-value ROOT", 40, 0, 100, 45, 0, 1);
545 fHistChi2ROOTCent = BookTH2F("fHistChi2ROOTCent", "centrality", "#tilde{#chi^{2}}", 40, 0, 100, 45, 0, 5);
546 fHistPChi2Root = BookTH2F("fHistPChi2Root", "p-value", "#tilde{#chi^{2}} ROOT", 1000, 0, 1, 100, 0, 5);
547 fHistPKolmogorov = BookTH2F("fHistPKolmogorov", "p-value", "kolmogorov p",40, 0, 1, 40, 0, 1);
548 fHistRhoStatusCent = BookTH2F("fHistRhoStatusCent", "centrality", "status [-1=lin was better, 0=ok, 1 = failed]", 101, -1, 100, 3, -1.5, 1.5);
549 fHistUndeterminedRunQA = BookTH1F("fHistUndeterminedRunQA", "runnumber", 10, 0, 10);
551 PostData(1, fOutputList);
553 switch (fRunModeType) {
555 fOutputListGood = new TList();
556 fOutputListGood->SetOwner(kTRUE);
557 fOutputListBad = new TList();
558 fOutputListBad->SetOwner(kTRUE);
559 PostData(2, fOutputListGood);
560 PostData(3, fOutputListBad);
565 // get the containers
566 fTracksCont = GetParticleContainer("Tracks");
567 fClusterCont = GetClusterContainer(0); // get the default cluster container
568 fJetsCont = GetJetContainer("Jets");
570 //_____________________________________________________________________________
571 Bool_t AliAnalysisTaskJetV2::Run()
573 // user exec: execute once for each event
574 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
575 if(!fTracks||!fJets||!fRho) return kFALSE;
576 if(!fLocalInit) fLocalInit = InitializeAnalysis();
577 // reject the event if expected data is missing
578 if(!PassesCuts(InputEvent())) return kFALSE;
579 fLeadingJet = GetLeadingJet(); // store the leading jet
581 fLocalRho->SetVal(fRho->GetVal());
582 // [0][0] psi2a [1,0] psi2c
583 // [0][1] psi3a [1,1] psi3c
584 Double_t vzero[2][2];
585 CalculateEventPlaneVZERO(vzero);
586 /* for the combined vzero event plane
588 * not fully implmemented yet, use with caution ! */
589 Double_t vzeroComb[2];
590 CalculateEventPlaneCombinedVZERO(vzeroComb);
593 CalculateEventPlaneTPC(tpc);
594 Double_t psi2(-1), psi3(-1);
595 // arrays which will hold the fit parameters
596 switch (fDetectorType) { // determine the detector type for the rho fit
597 case kTPC : { psi2 = tpc[0]; psi3 = tpc[1]; } break;
598 case kVZEROA : { psi2 = vzero[0][0]; psi3 = vzero[0][1]; } break;
599 case kVZEROC : { psi2 = vzero[1][0]; psi3 = vzero[1][1]; } break;
600 case kVZEROComb : { psi2 = vzeroComb[0]; psi3 = vzeroComb[1];} break;
603 switch (fFitModulationType) { // do the fits
605 switch (fCollisionType) {
606 case kPythia : { // background is zero for pp jets
607 fFitModulation->FixParameter(0, 0);
608 fLocalRho->SetVal(0);
611 fFitModulation->FixParameter(0, fLocalRho->GetVal());
615 case kV2 : { // only v2
616 if(CorrectRho(psi2, psi3)) {
617 fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
618 if(fUserSuppliedR2) {
619 Double_t r(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
620 if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
622 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
625 case kV3 : { // only v3
626 if(CorrectRho(psi2, psi3)) {
627 if(fUserSuppliedR3) {
628 Double_t r(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
629 if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
631 fProfV3->Fill(fCent, fFitModulation->GetParameter(3));
632 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
635 case kQC2 : { // qc2 analysis
636 if(CorrectRho(psi2, psi3)) {
637 if(fUserSuppliedR2 && fUserSuppliedR3) {
638 // note for the qc method, resolution is REVERSED to go back to v2obs
639 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
640 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
641 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
642 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)*r3);
644 if (fUsePtWeight) { // use weighted weights
645 Double_t dQCnM11 = (fNoEventWeightsForQC) ? 1. : QCnM11();
646 fProfV2->Fill(fCent, fFitModulation->GetParameter(3), dQCnM11);
647 fProfV3->Fill(fCent, fFitModulation->GetParameter(7), dQCnM11);
649 Double_t dQCnM = (fNoEventWeightsForQC) ? 2. : QCnM();
650 fProfV2->Fill(fCent, fFitModulation->GetParameter(3), dQCnM*(dQCnM-1));
651 fProfV3->Fill(fCent, fFitModulation->GetParameter(7), dQCnM*(dQCnM-1));
653 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
657 if(CorrectRho(psi2, psi3)) {
658 if(fUserSuppliedR2 && fUserSuppliedR3) {
659 // note for the qc method, resolution is REVERSED to go back to v2obs
660 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
661 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
662 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
663 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)*r3);
665 if (fUsePtWeight) { // use weighted weights
666 fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM1111()*/);
667 fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM1111()*/);
669 fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
670 fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
673 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
676 if(CorrectRho(psi2, psi3)) {
677 if(fUserSuppliedR2 && fUserSuppliedR3) {
678 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
679 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
680 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r2);
681 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)/r3);
683 fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
684 fProfV3->Fill(fCent, fFitModulation->GetParameter(7));
685 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
689 // if all went well, update the local rho parameter
690 fLocalRho->SetLocalRho(fFitModulation);
691 // fill a number of histograms. event qa needs to be filled first as it also determines the runnumber for the track qa
692 if(fFillQAHistograms) FillQAHistograms(InputEvent());
693 if(fFillHistograms) FillHistogramsAfterSubtraction(psi2, vzero, vzeroComb, tpc);
694 // send the output to the connected output container
695 PostData(1, fOutputList);
696 switch (fRunModeType) {
698 PostData(2, fOutputListGood);
699 PostData(3, fOutputListBad);
706 //_____________________________________________________________________________
707 void AliAnalysisTaskJetV2::CalculateEventPlaneVZERO(Double_t vzero[2][2]) const
709 // get the vzero event plane
710 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
711 if(fUseV0EventPlaneFromHeader) { // use the vzero from the header
712 Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0), h(0);
713 vzero[0][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 2, a, b);
714 vzero[1][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 2, c, d);
715 vzero[0][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 3, e, f);
716 vzero[1][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 3, g, h);
719 // grab the vzero event plane without recentering
720 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
721 Double_t qxa2(0), qya2(0), qxc2(0), qyc2(0); // for psi2
722 Double_t qxa3(0), qya3(0), qxc3(0), qyc3(0); // for psi3
723 for(Int_t iVZERO(0); iVZERO < 64; iVZERO++) {
724 Double_t phi(TMath::PiOver4()*(.5+iVZERO%8)), /* eta(0), */ weight(InputEvent()->GetVZEROEqMultiplicity(iVZERO));
726 qxa2 += weight*TMath::Cos(2.*phi);
727 qya2 += weight*TMath::Sin(2.*phi);
728 qxa3 += weight*TMath::Cos(3.*phi);
729 qya3 += weight*TMath::Sin(3.*phi);
732 qxc2 += weight*TMath::Cos(2.*phi);
733 qyc2 += weight*TMath::Sin(2.*phi);
734 qxc3 += weight*TMath::Cos(3.*phi);
735 qyc3 += weight*TMath::Sin(3.*phi);
738 vzero[0][0] = .5*TMath::ATan2(qya2, qxa2);
739 vzero[1][0] = .5*TMath::ATan2(qyc2, qxc2);
740 vzero[0][1] = (1./3.)*TMath::ATan2(qya3, qxa3);
741 vzero[1][1] = (1./3.)*TMath::ATan2(qyc3, qxc3);
743 //_____________________________________________________________________________
744 void AliAnalysisTaskJetV2::CalculateEventPlaneTPC(Double_t* tpc)
746 // grab the TPC event plane
747 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
748 fNAcceptedTracks = 0; // reset the track counter
749 Double_t qx2(0), qy2(0); // for psi2
750 Double_t qx3(0), qy3(0); // for psi3
752 Float_t excludeInEta = -999;
753 if(fExcludeLeadingJetsFromFit > 0 ) { // remove the leading jet from ep estimate
754 if(fLeadingJet) excludeInEta = fLeadingJet->Eta();
756 for(Int_t iTPC(0); iTPC < fTracksCont->GetNEntries(); iTPC++) {
757 AliVParticle* track = fTracksCont->GetParticle(iTPC);
758 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
759 if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta) < GetJetContainer()->GetJetRadius()*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - GetJetContainer()->GetJetRadius() - GetJetContainer()->GetJetEtaMax() ) > 0 )) continue;
761 qx2+= TMath::Cos(2.*track->Phi());
762 qy2+= TMath::Sin(2.*track->Phi());
763 qx3+= TMath::Cos(3.*track->Phi());
764 qy3+= TMath::Sin(3.*track->Phi());
767 tpc[0] = .5*TMath::ATan2(qy2, qx2);
768 tpc[1] = (1./3.)*TMath::ATan2(qy3, qx3);
770 //_____________________________________________________________________________
771 void AliAnalysisTaskJetV2::CalculateEventPlaneCombinedVZERO(Double_t* comb) const
773 // grab the combined vzero event plane
774 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
775 Double_t a(0), b(0), c(0), d(0);
776 comb[0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 10, 2, a, b);
777 comb[1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 10, 3, c, d);
779 //_____________________________________________________________________________
780 void AliAnalysisTaskJetV2::CalculateEventPlaneResolution(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc)
782 // fill the profiles for the resolution parameters
783 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
784 fProfV2Resolution[fInCentralitySelection]->Fill(2., TMath::Cos(2.*(vzero[0][0] - vzero[1][0])));
785 fProfV2Resolution[fInCentralitySelection]->Fill(3., TMath::Cos(2.*(vzero[1][0] - vzero[0][0])));
786 fProfV2Resolution[fInCentralitySelection]->Fill(4., TMath::Cos(2.*(vzero[0][0] - tpc[0])));
787 fProfV2Resolution[fInCentralitySelection]->Fill(5., TMath::Cos(2.*(tpc[0] - vzero[0][0])));
788 fProfV2Resolution[fInCentralitySelection]->Fill(6., TMath::Cos(2.*(vzero[1][0] - tpc[0])));
789 fProfV2Resolution[fInCentralitySelection]->Fill(7., TMath::Cos(2.*(tpc[0] - vzero[1][0])));
790 fProfV3Resolution[fInCentralitySelection]->Fill(2., TMath::Cos(3.*(vzero[0][0] - vzero[1][0])));
791 fProfV3Resolution[fInCentralitySelection]->Fill(3., TMath::Cos(3.*(vzero[1][0] - vzero[0][0])));
792 fProfV3Resolution[fInCentralitySelection]->Fill(4., TMath::Cos(3.*(vzero[0][0] - tpc[0])));
793 fProfV3Resolution[fInCentralitySelection]->Fill(5., TMath::Cos(3.*(tpc[0] - vzero[0][0])));
794 fProfV3Resolution[fInCentralitySelection]->Fill(6., TMath::Cos(3.*(vzero[1][0] - tpc[0])));
795 fProfV3Resolution[fInCentralitySelection]->Fill(7., TMath::Cos(3.*(tpc[0] - vzero[1][0])));
796 // for the resolution of the combined vzero event plane, use two tpc halves as uncorrelated subdetectors
797 Double_t qx2a(0), qy2a(0); // for psi2a, negative eta
798 Double_t qx3a(0), qy3a(0); // for psi3a, negative eta
799 Double_t qx2b(0), qy2b(0); // for psi2a, positive eta
800 Double_t qx3b(0), qy3b(0); // for psi3a, positive eta
802 Int_t iTracks(fTracks->GetEntriesFast());
803 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
804 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
805 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
806 if(track->Eta() < 0 ) {
807 qx2a+= TMath::Cos(2.*track->Phi());
808 qy2a+= TMath::Sin(2.*track->Phi());
809 qx3a+= TMath::Cos(3.*track->Phi());
810 qy3a+= TMath::Sin(3.*track->Phi());
811 } else if (track->Eta() > 0) {
812 qx2b+= TMath::Cos(2.*track->Phi());
813 qy2b+= TMath::Sin(2.*track->Phi());
814 qx3b+= TMath::Cos(3.*track->Phi());
815 qy3b+= TMath::Sin(3.*track->Phi());
819 Double_t tpca2(.5*TMath::ATan2(qy2a, qx2a));
820 Double_t tpca3((1./3.)*TMath::ATan2(qy3a, qx3a));
821 Double_t tpcb2(.5*TMath::ATan2(qy2b, qx2b));
822 Double_t tpcb3((1./3.)*TMath::ATan2(qy3b, qx3b));
823 fProfV2Resolution[fInCentralitySelection]->Fill(8., TMath::Cos(2.*(vzeroComb[0] - tpca2)));
824 fProfV2Resolution[fInCentralitySelection]->Fill(9., TMath::Cos(2.*(vzeroComb[0] - tpcb2)));
825 fProfV2Resolution[fInCentralitySelection]->Fill(10., TMath::Cos(2.*(tpca2 - tpcb2)));
826 fProfV3Resolution[fInCentralitySelection]->Fill(8., TMath::Cos(3.*(vzeroComb[1] - tpca3)));
827 fProfV3Resolution[fInCentralitySelection]->Fill(9., TMath::Cos(3.*(vzeroComb[1] - tpcb3)));
828 fProfV3Resolution[fInCentralitySelection]->Fill(10., TMath::Cos(3.*(tpca3 - tpcb3)));
830 //_____________________________________________________________________________
831 void AliAnalysisTaskJetV2::CalculateRandomCone(Float_t &pt, Float_t &eta, Float_t &phi,
832 AliParticleContainer* tracksCont, AliClusterContainer* clusterCont, AliEmcalJet* jet) const
835 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
836 pt = 0; eta = 0; phi = 0;
837 Float_t etaJet(999), phiJet(999), dJet(999); // no jet: same as jet very far away
838 if(jet) { // if a leading jet is given, use its kinematic properties to exclude it
842 // the random cone acceptance has to equal the jet acceptance
843 // this also insures safety when runnnig on the semi-good tpc runs for 11h data,
844 // where jet acceptance is adjusted to reduced acceptance - hence random cone acceptance as well
845 Float_t minPhi(GetJetContainer()->GetJetPhiMin()), maxPhi(GetJetContainer()->GetJetPhiMax());
846 if(maxPhi > TMath::TwoPi()) maxPhi = TMath::TwoPi();
847 if(minPhi < 0 ) minPhi = 0.;
848 // construct a random cone and see if it's far away enough from the leading jet
849 Int_t attempts(1000);
852 eta = gRandom->Uniform(GetJetContainer()->GetJetEtaMin(), GetJetContainer()->GetJetEtaMax());
853 phi = gRandom->Uniform(minPhi, maxPhi);
855 dJet = TMath::Sqrt((etaJet-eta)*(etaJet-eta)+(phiJet-phi)*(phiJet-phi));
856 if(dJet > fMinDisanceRCtoLJ) break;
857 else if (attempts == 0) {
858 printf(" > No random cone after 1000 tries, giving up ... !\n");
862 // get the charged energy (if tracks are provided)
864 AliVParticle* track = tracksCont->GetNextAcceptParticle(0);
866 Float_t etaTrack(track->Eta()), phiTrack(track->Phi());
867 // get distance from cone
868 if(TMath::Abs(phiTrack-phi) > TMath::Abs(phiTrack - phi + TMath::TwoPi())) phiTrack+=TMath::TwoPi();
869 if(TMath::Abs(phiTrack-phi) > TMath::Abs(phiTrack - phi - TMath::TwoPi())) phiTrack-=TMath::TwoPi();
870 if(TMath::Sqrt(TMath::Abs((etaTrack-eta)*(etaTrack-eta)+(phiTrack-phi)*(phiTrack-phi))) <= GetJetRadius()) pt += track->Pt();
871 track = tracksCont->GetNextAcceptParticle();
874 // get the neutral energy (if clusters are provided)
876 AliVCluster* cluster = clusterCont->GetNextAcceptCluster(0);
878 TLorentzVector momentum;
879 cluster->GetMomentum(momentum, const_cast<Double_t*>(fVertex));
880 Float_t etaClus(momentum.Eta()), phiClus(momentum.Phi());
881 // get distance from cone
882 if(TMath::Abs(phiClus-phi) > TMath::Abs(phiClus - phi + TMath::TwoPi())) phiClus+=TMath::TwoPi();
883 if(TMath::Abs(phiClus-phi) > TMath::Abs(phiClus - phi - TMath::TwoPi())) phiClus-=TMath::TwoPi();
884 if(TMath::Sqrt(TMath::Abs((etaClus-eta)*(etaClus-eta)+(phiClus-phi)*(phiClus-phi))) <= GetJetRadius()) pt += momentum.Pt();
885 cluster = clusterCont->GetNextAcceptCluster();
889 //_____________________________________________________________________________
890 Double_t AliAnalysisTaskJetV2::CalculateQC2(Int_t harm) {
891 // get the second order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho
892 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
893 Double_t reQ(0), imQ(0), modQ(0), M11(0), M(0);
894 if(fUsePtWeight) { // for the weighted 2-nd order q-cumulant
895 QCnQnk(harm, 1, reQ, imQ); // get the weighted 2-nd order q-vectors
896 modQ = reQ*reQ+imQ*imQ; // get abs Q-squared
897 M11 = QCnM11(); // equals S2,1 - S1,2
898 return (M11 > 0) ? ((modQ - QCnS(1,2))/M11) : -999;
899 } // else return the non-weighted 2-nd order q-cumulant
900 QCnQnk(harm, 0, reQ, imQ); // get the non-weighted 2-nd order q-vectors
901 modQ = reQ*reQ+imQ*imQ; // get abs Q-squared
903 return (M > 1) ? (modQ - M)/(M*(M-1)) : -999;
905 //_____________________________________________________________________________
906 Double_t AliAnalysisTaskJetV2::CalculateQC4(Int_t harm) {
907 // get the fourth order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho
908 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
909 Double_t reQn1(0), imQn1(0), reQ2n2(0), imQ2n2(0), reQn3(0), imQn3(0), M1111(0), M(0);
910 Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0); // terms of the calculation
911 if(fUsePtWeight) { // for the weighted 4-th order q-cumulant
912 QCnQnk(harm, 1, reQn1, imQn1);
913 QCnQnk(harm*2, 2, reQ2n2, imQ2n2);
914 QCnQnk(harm, 3, reQn3, imQn3);
915 // fill in the terms ...
916 a = (reQn1*reQn1+imQn1*imQn1)*(reQn1*reQn1+imQn1*imQn1);
917 b = reQ2n2*reQ2n2 + imQ2n2*imQ2n2;
918 c = -2.*(reQ2n2*reQn1*reQn1-reQ2n2*imQn1*imQn1+2.*imQ2n2*reQn1*imQn1);
919 d = 8.*(reQn3*reQn1+imQn3*imQn1);
920 e = -4.*QCnS(1,2)*(reQn1*reQn1+imQn1*imQn1);
924 return (M1111 > 0) ? (a+b+c+d+e+f+g)/M1111 : -999;
925 } // else return the unweighted case
926 Double_t reQn(0), imQn(0), reQ2n(0), imQ2n(0);
927 QCnQnk(harm, 0, reQn, imQn);
928 QCnQnk(harm*2, 0, reQ2n, imQ2n);
929 // fill in the terms ...
931 if(M < 4) return -999;
932 a = (reQn*reQn+imQn*imQn)*(reQn*reQn+imQn*imQn);
933 b = reQ2n*reQ2n + imQ2n*imQ2n;
934 c = -2.*(reQ2n*reQn*reQn-reQ2n*imQn*imQn+2.*imQ2n*reQn*imQn);
935 e = -4.*(M-2)*(reQn*reQn+imQn*imQn);
937 return (a+b+c+e+f)/(M*(M-1)*(M-2)*(M-3));
939 //_____________________________________________________________________________
940 void AliAnalysisTaskJetV2::QCnQnk(Int_t n, Int_t k, Double_t &reQ, Double_t &imQ) {
941 // get the weighted n-th order q-vector, pass real and imaginary part as reference
942 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
944 fNAcceptedTracksQCn = 0;
945 Int_t iTracks(fTracks->GetEntriesFast());
946 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
947 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
948 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
949 fNAcceptedTracksQCn++;
950 // for the unweighted case, k equals zero and the weight doesn't contribute to the equation below
951 reQ += TMath::Power(track->Pt(), k) * TMath::Cos(((double)n)*track->Phi());
952 imQ += TMath::Power(track->Pt(), k) * TMath::Sin(((double)n)*track->Phi());
955 //_____________________________________________________________________________
956 void AliAnalysisTaskJetV2::QCnDiffentialFlowVectors(
957 TClonesArray* pois, TArrayD* ptBins, Bool_t vpart, Double_t* repn, Double_t* impn,
958 Double_t *mp, Double_t *reqn, Double_t *imqn, Double_t* mq, Int_t n)
960 // get unweighted differential flow vectors
961 Int_t iPois(pois->GetEntriesFast());
963 for(Int_t i(0); i < iPois; i++) {
964 for(Int_t ptBin(0); ptBin < ptBins->GetSize()-1; ptBin++) {
965 AliVTrack* poi = static_cast<AliVTrack*>(pois->At(i));
966 if(PassesCuts(poi)) {
967 if(poi->Pt() >= ptBins->At(ptBin) && poi->Pt() < ptBins->At(ptBin+1)) {
968 // fill the flow vectors assuming that all poi's are in the rp selection (true by design)
969 repn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
970 impn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
972 reqn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
973 imqn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
980 for(Int_t i(0); i < iPois; i++) {
981 for(Int_t ptBin(0); ptBin < ptBins->GetSize()-1; ptBin++) {
982 AliEmcalJet* poi = static_cast<AliEmcalJet*>(pois->At(i));
983 if(PassesCuts(poi)) {
984 Double_t pt(poi->Pt()-poi->Area()*fLocalRho->GetLocalVal(poi->Phi(), GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
985 if(pt >= ptBins->At(ptBin) && pt < ptBins->At(ptBin+1)) {
986 repn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
987 impn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
988 mp[ptBin]++; // qn isn't filled, no overlap between poi's and rp's
995 //_____________________________________________________________________________
996 Double_t AliAnalysisTaskJetV2::QCnS(Int_t i, Int_t j) {
997 // get the weighted ij-th order autocorrelation correction
998 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
999 if(!fTracks || i <= 0 || j <= 0) return -999;
1000 Int_t iTracks(fTracks->GetEntriesFast());
1002 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
1003 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
1004 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
1005 Sij+=TMath::Power(track->Pt(), j);
1007 return TMath::Power(Sij, i);
1009 //_____________________________________________________________________________
1010 Double_t AliAnalysisTaskJetV2::QCnM() {
1011 // get multiplicity for unweighted q-cumulants. function QCnQnk should be called first
1012 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1013 return (Double_t) fNAcceptedTracksQCn;
1015 //_____________________________________________________________________________
1016 Double_t AliAnalysisTaskJetV2::QCnM11() {
1017 // get multiplicity weights for the weighted two particle cumulant
1018 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1019 return (QCnS(2,1) - QCnS(1,2));
1021 //_____________________________________________________________________________
1022 Double_t AliAnalysisTaskJetV2::QCnM1111() {
1023 // get multiplicity weights for the weighted four particle cumulant
1024 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1025 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));
1027 //_____________________________________________________________________________
1028 Bool_t AliAnalysisTaskJetV2::QCnRecovery(Double_t psi2, Double_t psi3) {
1029 // decides how to deal with the situation where c2 or c3 is negative
1030 // returns kTRUE depending on whether or not a modulated rho is used for the jet background
1031 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1032 if(TMath::AreEqualAbs(fFitModulation->GetParameter(3), .0, 1e-10) && TMath::AreEqualAbs(fFitModulation->GetParameter(7), .0,1e-10)) {
1033 fFitModulation->SetParameter(7, 0);
1034 fFitModulation->SetParameter(3, 0);
1035 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1036 return kTRUE; // v2 and v3 have physical null values
1038 switch (fQCRecovery) {
1039 case kFixedRho : { // roll back to the original rho
1040 fFitModulation->SetParameter(7, 0);
1041 fFitModulation->SetParameter(3, 0);
1042 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1043 return kFALSE; // rho is forced to be fixed
1045 case kNegativeVn : {
1046 Double_t c2(fFitModulation->GetParameter(3));
1047 Double_t c3(fFitModulation->GetParameter(7));
1048 if( c2 < 0 ) c2 = -1.*TMath::Sqrt(-1.*c2);
1049 if( c3 < 0 ) c3 = -1.*TMath::Sqrt(-1.*c3);
1050 fFitModulation->SetParameter(3, c2);
1051 fFitModulation->SetParameter(7, c3);
1052 return kTRUE; // is this a physical quantity ?
1055 fitModulationType tempType(fFitModulationType); // store temporarily
1056 fFitModulationType = kCombined;
1057 fFitModulation->SetParameter(7, 0);
1058 fFitModulation->SetParameter(3, 0);
1059 Bool_t pass(CorrectRho(psi2, psi3)); // do the fit and all quality checks
1060 fFitModulationType = tempType; // roll back for next event
1063 default : return kFALSE;
1067 //_____________________________________________________________________________
1068 Bool_t AliAnalysisTaskJetV2::CorrectRho(Double_t psi2, Double_t psi3)
1070 // get rho' -> rho(phi)
1071 // two routines are available, both can be used with or without pt weights
1072 // [1] get vn from q-cumulants or as an integrated value from a user supplied histogram
1073 // in case of cumulants, both cumulants and vn values are stored. in both cases, v2 and v3
1074 // are expected. a check is performed to see if rho has no negative local minimum
1075 // for full description, see Phys. Rev. C 83, 044913
1076 // since the cn distribution has negative values, vn = sqrt(cn) can be imaginary sometimes
1077 // in this case one can either roll back to the 'original' rixed rho, do a fit for vn or take use
1078 // vn = - sqrt(|cn|)
1079 // [2] fitting a fourier expansion to the de/dphi distribution
1080 // the fit can be done with either v2, v3 or a combination.
1081 // in all cases, a cut can be made on the p-value of the chi-squared value of the fit
1082 // and a check can be performed to see if rho has no negative local minimum
1083 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1084 Int_t freeParams(2); // free parameters of the fit (for NDF)
1085 switch (fFitModulationType) { // for approaches where no fitting is required
1087 fFitModulation->FixParameter(4, psi2);
1088 fFitModulation->FixParameter(6, psi3);
1089 fFitModulation->FixParameter(3, CalculateQC2(2)); // set here with cn, vn = sqrt(cn)
1090 fFitModulation->FixParameter(7, CalculateQC2(3));
1091 // first fill the histos of the raw cumulant distribution
1092 if (fUsePtWeight) { // use weighted weights
1093 Double_t dQCnM11 = (fNoEventWeightsForQC) ? 1. : QCnM11();
1094 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), dQCnM11);
1095 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), dQCnM11);
1097 Double_t dQCnM = (fNoEventWeightsForQC) ? 2. : QCnM();
1098 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), dQCnM*(dQCnM-1));
1099 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), dQCnM*(dQCnM-1));
1101 // then see if one of the cn value is larger than zero and vn is readily available
1102 if(fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) {
1103 fFitModulation->FixParameter(3, TMath::Sqrt(fFitModulation->GetParameter(3)));
1104 fFitModulation->FixParameter(7, TMath::Sqrt(fFitModulation->GetParameter(7)));
1105 } else if (!QCnRecovery(psi2, psi3)) return kFALSE; // try to recover the cumulant, this will set v2 and v3
1106 if(fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { // general check
1107 fFitModulation->SetParameter(7, 0);
1108 fFitModulation->SetParameter(3, 0);
1109 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1115 fFitModulation->FixParameter(4, psi2);
1116 fFitModulation->FixParameter(6, psi3);
1117 fFitModulation->FixParameter(3, CalculateQC4(2)); // set here with cn, vn = sqrt(cn)
1118 fFitModulation->FixParameter(7, CalculateQC4(3));
1119 // first fill the histos of the raw cumulant distribution
1120 if (fUsePtWeight) { // use weighted weights
1121 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
1122 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
1124 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
1125 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
1127 // then see if one of the cn value is larger than zero and vn is readily available
1128 if(fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) {
1129 fFitModulation->FixParameter(3, TMath::Sqrt(fFitModulation->GetParameter(3)));
1130 fFitModulation->FixParameter(7, TMath::Sqrt(fFitModulation->GetParameter(7)));
1131 } else if (!QCnRecovery(psi2, psi3)) return kFALSE; // try to recover the cumulant, this will set v2 and v3
1132 if(fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { // general check
1133 fFitModulation->SetParameter(7, 0);
1134 fFitModulation->SetParameter(3, 0);
1135 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1139 case kIntegratedFlow : {
1140 // use v2 and v3 values from an earlier iteration over the data
1141 fFitModulation->FixParameter(3, fUserSuppliedV2->GetBinContent(fUserSuppliedV2->GetXaxis()->FindBin(fCent)));
1142 fFitModulation->FixParameter(4, psi2);
1143 fFitModulation->FixParameter(6, psi3);
1144 fFitModulation->FixParameter(7, fUserSuppliedV3->GetBinContent(fUserSuppliedV3->GetXaxis()->FindBin(fCent)));
1145 if(fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) {
1146 fFitModulation->SetParameter(7, 0);
1147 fFitModulation->SetParameter(3, 0);
1148 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1155 TString detector("");
1156 switch (fDetectorType) {
1157 case kTPC : detector+="TPC";
1159 case kVZEROA : detector+="VZEROA";
1161 case kVZEROC : detector+="VZEROC";
1163 case kVZEROComb : detector+="VZEROComb";
1167 Int_t iTracks(fTracks->GetEntriesFast());
1168 Double_t excludeInEta = -999;
1169 Double_t excludeInPhi = -999;
1170 Double_t excludeInPt = -999;
1171 if(iTracks <= 0 || fLocalRho->GetVal() <= 0 ) return kFALSE; // no use fitting an empty event ...
1172 if(fExcludeLeadingJetsFromFit > 0 ) {
1174 excludeInEta = fLeadingJet->Eta();
1175 excludeInPhi = fLeadingJet->Phi();
1176 excludeInPt = fLeadingJet->Pt();
1179 // check the acceptance of the track selection that will be used
1180 // if one uses e.g. semi-good tpc tracks, accepance in phi is reduced to 0 < phi < 4
1181 // the defaults (-10 < phi < 10) which accept all, are then overwritten
1182 Double_t lowBound(0.), upBound(TMath::TwoPi()); // bounds for fit
1183 if(GetParticleContainer()->GetParticlePhiMin() > lowBound) lowBound = GetParticleContainer()->GetParticlePhiMin();
1184 if(GetParticleContainer()->GetParticlePhiMax() < upBound) upBound = GetParticleContainer()->GetParticlePhiMax();
1186 fHistSwap->Reset(); // clear the histogram
1187 TH1F _tempSwap; // on stack for quick access
1188 TH1F _tempSwapN; // on stack for quick access, bookkeeping histogram
1189 if(fRebinSwapHistoOnTheFly) {
1190 if(fNAcceptedTracks < 49) fNAcceptedTracks = 49; // avoid aliasing effects
1191 _tempSwap = TH1F("_tempSwap", "_tempSwap", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), lowBound, upBound);
1192 if(fUsePtWeightErrorPropagation) _tempSwapN = TH1F("_tempSwapN", "_tempSwapN", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), lowBound, upBound);
1193 if(fUsePtWeight) _tempSwap.Sumw2();
1195 else _tempSwap = *fHistSwap; // now _tempSwap holds the desired histo
1196 // non poissonian error when using pt weights
1197 Double_t totalpts(0.), totalptsquares(0.), totalns(0.);
1198 for(Int_t i(0); i < iTracks; i++) {
1199 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1200 if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta) < GetJetContainer()->GetJetRadius()*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - GetJetContainer()->GetJetRadius() - GetJetContainer()->GetJetEtaMax() ) > 0 )) continue;
1201 if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
1203 _tempSwap.Fill(track->Phi(), track->Pt());
1204 if(fUsePtWeightErrorPropagation) {
1205 totalpts += track->Pt();
1206 totalptsquares += track->Pt()*track->Pt();
1208 _tempSwapN.Fill(track->Phi());
1211 else _tempSwap.Fill(track->Phi());
1213 if(fUsePtWeight && fUsePtWeightErrorPropagation) {
1214 // in the case of pt weights overwrite the poissonian error estimate which is assigned by root by a more sophisticated appraoch
1215 // 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
1216 // 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
1217 // event, for the latter use a poissonian estimate. the two contrubitions are assumed to be uncorrelated
1218 if(totalns < 1) return kFALSE; // not one track passes the cuts
1219 for(Int_t l = 0; l < _tempSwap.GetNbinsX(); l++) {
1220 if(_tempSwapN.GetBinContent(l+1) == 0) {
1221 _tempSwap.SetBinContent(l+1,0);
1222 _tempSwap.SetBinError(l+1,0);
1225 Double_t vartimesnsq = totalptsquares*totalns - totalpts*totalpts;
1226 Double_t variance = vartimesnsq/(totalns*(totalns-1.));
1227 Double_t SDOMSq = variance / _tempSwapN.GetBinContent(l+1);
1228 Double_t SDOMSqOverMeanSq = SDOMSq * _tempSwapN.GetBinContent(l+1) * _tempSwapN.GetBinContent(l+1) / (_tempSwapN.GetBinContent(l+1) * _tempSwapN.GetBinContent(l+1));
1229 Double_t poissonfrac = 1./_tempSwapN.GetBinContent(l+1);
1230 Double_t vartotalfrac = SDOMSqOverMeanSq + poissonfrac;
1231 Double_t vartotal = vartotalfrac * _tempSwap.GetBinContent(l+1) * _tempSwap.GetBinContent(l+1);
1232 if(vartotal > 0.0001) _tempSwap.SetBinError(l+1,TMath::Sqrt(vartotal));
1234 _tempSwap.SetBinContent(l+1,0);
1235 _tempSwap.SetBinError(l+1,0);
1241 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1242 switch (fFitModulationType) {
1244 fFitModulation->FixParameter(0, fLocalRho->GetVal() );
1248 fFitModulation->FixParameter(4, psi2);
1252 fFitModulation->FixParameter(4, psi3);
1256 fFitModulation->FixParameter(4, psi2);
1257 fFitModulation->FixParameter(6, psi3);
1260 case kFourierSeries : {
1261 // in this approach, an explicit calculation will be made of vn = sqrt(xn^2+yn^2)
1262 // where x[y] = Integrate[r(phi)cos[sin](n phi)dphi, 0, 2pi]
1263 Double_t cos2(0), sin2(0), cos3(0), sin3(0), sumPt(0);
1264 for(Int_t i(0); i < iTracks; i++) {
1265 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1266 if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
1267 sumPt += track->Pt();
1268 cos2 += track->Pt()*TMath::Cos(2*PhaseShift(track->Phi()-psi2));
1269 sin2 += track->Pt()*TMath::Sin(2*PhaseShift(track->Phi()-psi2));
1270 cos3 += track->Pt()*TMath::Cos(3*PhaseShift(track->Phi()-psi3));
1271 sin3 += track->Pt()*TMath::Sin(3*PhaseShift(track->Phi()-psi3));
1273 fFitModulation->SetParameter(3, TMath::Sqrt(cos2*cos2+sin2*sin2)/fLocalRho->GetVal());
1274 fFitModulation->SetParameter(4, psi2);
1275 fFitModulation->SetParameter(6, psi3);
1276 fFitModulation->SetParameter(7, TMath::Sqrt(cos3*cos3+sin3*sin3)/fLocalRho->GetVal());
1281 // toy mc, just here to check procedure, azimuthal profile is filled from hypothesis so p-value distribution should be flat
1282 Int_t _bins = _tempSwap.GetXaxis()->GetNbins();
1283 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());
1284 _tempFit->SetParameter(0, fFitModulation->GetParameter(0)); // normalization
1285 _tempFit->SetParameter(3, 0.1); // v2
1286 _tempFit->FixParameter(1, 1.); // constant
1287 _tempFit->FixParameter(2, 2.); // constant
1288 _tempFit->FixParameter(5, 3.); // constant
1289 _tempFit->FixParameter(4, fFitModulation->GetParameter(4));
1290 _tempFit->FixParameter(6, fFitModulation->GetParameter(6));
1291 _tempFit->SetParameter(7, 0.1); // v3
1292 _tempSwap.Reset(); // rese bin content
1293 for(int _binsI = 0; _binsI < _bins*_bins; _binsI++) _tempSwap.Fill(_tempFit->GetRandom());
1295 _tempSwap.Fit(fFitModulation, fFitModulationOptions.Data(), "", lowBound, upBound);
1296 // the quality of the fit is evaluated from 1 - the cdf of the chi square distribution
1297 // three methods are available, all with their drawbacks. all are stored, one is selected to do the cut
1298 Int_t NDF(_tempSwap.GetXaxis()->GetNbins()-freeParams);
1299 if(NDF == 0) return kFALSE;
1300 Double_t CDF(1.-ChiSquareCDF(NDF, ChiSquare(_tempSwap, fFitModulation)));
1301 Double_t CDFROOT(1.-ChiSquareCDF(NDF, fFitModulation->GetChisquare()));
1302 Double_t CDFKolmogorov(KolmogorovTest(_tempSwap, fFitModulation));
1303 // fill the values and centrality correlation (redundant but easy on the eyes)
1304 fHistPvalueCDF->Fill(CDF);
1305 fHistPvalueCDFCent->Fill(fCent, CDF);
1306 fHistPvalueCDFROOT->Fill(CDFROOT);
1307 fHistPvalueCDFROOTCent->Fill(fCent, CDFROOT);
1308 fHistKolmogorovTest->Fill(CDFKolmogorov);
1309 fHistChi2ROOTCent->Fill(fCent, fFitModulation->GetChisquare()/((float)NDF));
1310 fHistChi2Cent->Fill(fCent, ChiSquare(_tempSwap, fFitModulation)/((float)NDF));
1311 fHistKolmogorovTestCent->Fill(fCent, CDFKolmogorov);
1312 fHistPChi2Root->Fill(CDFROOT, fFitModulation->GetChisquare()/((float)NDF));
1313 fHistPChi2->Fill(CDF, ChiSquare(_tempSwap, fFitModulation)/((float)NDF));
1314 fHistPKolmogorov->Fill(CDF, CDFKolmogorov);
1316 // variable CDF is used for making cuts, so we fill it with the selected p-value
1317 switch (fFitGoodnessTest) {
1321 case kChi2Poisson : break; // CDF is already CDF
1322 case kKolmogorov : {
1323 CDF = CDFKolmogorov;
1329 // as an additional quality check, see if fitting a control fit has a higher significance
1330 _tempSwap.Fit(fFitControl, fFitModulationOptions.Data(), "", lowBound, upBound);
1331 Double_t CDFControl(-1.);
1332 switch (fFitGoodnessTest) {
1334 CDFControl = 1.-ChiSquareCDF(fFitControl->GetNDF(), fFitModulation->GetChisquare());
1336 case kChi2Poisson : {
1337 CDFControl = 1.-ChiSquareCDF(fFitControl->GetNDF(), ChiSquare(_tempSwap, fFitModulation));
1339 case kKolmogorov : {
1340 CDFControl = KolmogorovTest(_tempSwap, fFitControl);
1344 if(CDFControl > CDF) {
1345 CDF = -1.; // control fit is more significant, so throw out the 'old' fit
1346 fHistRhoStatusCent->Fill(fCent, -1);
1349 if(CDF >= fMinPvalue && CDF <= fMaxPvalue && ( fFitModulation->GetMinimum(0, TMath::TwoPi()) > 0)) { // fit quality. not that although with limited acceptance the fit is performed on just
1350 // part of phase space, the requirement that energy desntiy is larger than zero is applied
1351 // to the FULL spectrum
1352 fHistRhoStatusCent->Fill(fCent, 0.);
1353 // for LOCAL didactic purposes, save the best and the worst fits
1354 // this routine can produce a lot of output histograms (it's not memory 'safe') and will not work on GRID
1355 // since the output will become unmergeable (i.e. different nodes may produce conflicting output)
1356 switch (fRunModeType) {
1358 if(fRandom->Uniform(0, 100) > fPercentageOfFits) break;
1359 static Int_t didacticCounterBest(0);
1360 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
1361 TF1* didacticFit = (TF1*)fFitModulation->Clone(Form("fit_%i_CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
1362 switch(fFitModulationType) {
1364 // to make a nice picture also plot the separate components (v2 and v3) of the fit
1365 // only done for cobined fit where there are actually components to split ...
1366 TF1* v0(new TF1("dfit_kV2", "[0]", 0, TMath::TwoPi()));
1367 v0->SetParameter(0, didacticFit->GetParameter(0)); // normalization
1368 v0->SetLineColor(kMagenta);
1369 v0->SetLineStyle(7);
1370 didacticProfile->GetListOfFunctions()->Add(v0);
1371 TF1* v2(new TF1("dfit_kV2", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
1372 v2->SetParameter(0, didacticFit->GetParameter(0)); // normalization
1373 v2->SetParameter(3, didacticFit->GetParameter(3)); // v2
1374 v2->FixParameter(1, 1.); // constant
1375 v2->FixParameter(2, 2.); // constant
1376 v2->FixParameter(4, didacticFit->GetParameter(4)); // psi2
1377 v2->SetLineColor(kGreen);
1378 didacticProfile->GetListOfFunctions()->Add(v2);
1379 TF1* v3(new TF1("dfit_kV3", "[0]*([1]+[2]*[3]*TMath::Cos([5]*(x-[4])))", 0, TMath::TwoPi()));
1380 v3->SetParameter(0, didacticFit->GetParameter(0)); // normalization
1381 v3->SetParameter(3, didacticFit->GetParameter(7)); // v3
1382 v3->FixParameter(1, 1.); // constant
1383 v3->FixParameter(2, 2.); // constant
1384 v3->FixParameter(4, didacticFit->GetParameter(6)); // psi3
1385 v3->FixParameter(5, 3.); // constant
1386 v3->SetLineColor(kCyan);
1387 didacticProfile->GetListOfFunctions()->Add(v3);
1391 didacticProfile->GetListOfFunctions()->Add(didacticFit);
1392 didacticProfile->GetYaxis()->SetTitle("#frac{d #sum #it{p}_{T}}{d #varphi} [GeV/#it{c}]");
1393 didacticProfile->GetXaxis()->SetTitle("#varphi");
1394 fOutputListGood->Add(didacticProfile);
1395 didacticCounterBest++;
1396 TH2F* didacticSurface = BookTH2F(Form("surface_%s", didacticProfile->GetName()), "#phi", "#eta", 50, 0, TMath::TwoPi(), 50, -1, 1, -1, kFALSE);
1397 for(Int_t i(0); i < iTracks; i++) {
1398 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1399 if(PassesCuts(track)) {
1400 if(fUsePtWeight) didacticSurface->Fill(track->Phi(), track->Eta(), track->Pt());
1401 else didacticSurface->Fill(track->Phi(), track->Eta());
1404 if(fExcludeLeadingJetsFromFit) { // visualize the excluded region
1405 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);
1406 f2->SetParameters(excludeInPt/3.,excludeInPhi,.1,excludeInEta,.1);
1407 didacticSurface->GetListOfFunctions()->Add(f2);
1409 fOutputListGood->Add(didacticSurface);
1413 } else { // if the fit is of poor quality revert to the original rho estimate
1414 switch (fRunModeType) { // again see if we want to save the fit
1416 static Int_t didacticCounterWorst(0);
1417 if(fRandom->Uniform(0, 100) > fPercentageOfFits) break;
1418 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data() ));
1419 TF1* didacticFit = (TF1*)fFitModulation->Clone(Form("fit_%i_p_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data()));
1420 didacticProfile->GetListOfFunctions()->Add(didacticFit);
1421 fOutputListBad->Add(didacticProfile);
1422 didacticCounterWorst++;
1426 switch (fFitModulationType) {
1427 case kNoFit : break; // nothing to do
1428 case kCombined : fFitModulation->SetParameter(7, 0); // no break
1429 case kFourierSeries : fFitModulation->SetParameter(7, 0); // no break
1430 default : { // needs to be done if there was a poor fit
1431 fFitModulation->SetParameter(3, 0);
1432 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1435 if(CDF > -.5) fHistRhoStatusCent->Fill(fCent, 1.);
1436 return kFALSE; // return false if the fit is rejected
1440 //_____________________________________________________________________________
1441 Bool_t AliAnalysisTaskJetV2::PassesCuts(AliVEvent* event)
1444 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1445 if(!event || !AliAnalysisTaskEmcal::IsEventSelected()) return kFALSE;
1446 if(TMath::Abs(InputEvent()->GetPrimaryVertex()->GetZ()) > 10.) return kFALSE;
1447 // aod and esd specific checks
1448 switch (fDataType) {
1450 AliESDEvent* esdEvent = static_cast<AliESDEvent*>(InputEvent());
1451 if( (!esdEvent) || (TMath::Abs(esdEvent->GetPrimaryVertexSPD()->GetZ() - esdEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE;
1454 AliAODEvent* aodEvent = static_cast<AliAODEvent*>(InputEvent());
1455 if( (!aodEvent) || (TMath::Abs(aodEvent->GetPrimaryVertexSPD()->GetZ() - aodEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE;
1459 fCent = InputEvent()->GetCentrality()->GetCentralityPercentile("V0M");
1460 if(fCent <= fCentralityClasses->At(0) || fCent >= fCentralityClasses->At(fCentralityClasses->GetSize()-1) || TMath::Abs(fCent-InputEvent()->GetCentrality()->GetCentralityPercentile("TRK")) > 5.) return kFALSE;
1461 // determine centrality class
1462 fInCentralitySelection = -1;
1463 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
1464 if(fCent >= fCentralityClasses->At(i) && fCent <= fCentralityClasses->At(1+i)) {
1465 fInCentralitySelection = i;
1468 if(fInCentralitySelection<0) return kFALSE; // should be null op
1469 /* if(fExplicitOutlierCut == 2010 || fExplicitOutlierCut == 2011) {
1470 if(!PassesCuts(fExplicitOutlierCut)) return kFALSE;
1472 // see if input containers are filled
1473 if(fTracks->GetEntries() < 1) return kFALSE;
1474 if(fRho->GetVal() <= 0 ) return kFALSE;
1475 if(fAnalysisType == AliAnalysisTaskJetV2::kFull && !fClusterCont) return kFALSE;
1478 //_____________________________________________________________________________
1479 /*Bool_t AliAnalysisTaskJetV2::PassesCuts(Int_t year)
1481 // additional centrality cut based on relation between tpc and global multiplicity
1482 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1483 AliAODEvent* event(dynamic_cast<AliAODEvent*>(InputEvent()));
1484 if(!event) return kFALSE;
1485 Int_t multTPC(0), multGlob(0), nTracks(InputEvent()->GetNumberOfTracks());
1486 for(Int_t iTracks = 0; iTracks < nTracks; iTracks++) {
1487 AliAODTrack* track = event->GetTrack(iTracks);
1488 if(!track) continue;
1489 if (!track || track->Pt() < .2 || track->Pt() > 5.0 || TMath::Abs(track->Eta()) > .8 || track->GetTPCNcls() < 70 || !track->GetDetPid() || track->GetDetPid()->GetTPCsignal() < 10.0) continue; // general quality cut
1490 if (track->TestFilterBit(1) && track->Chi2perNDF() > 0.2) multTPC++;
1491 if (!track->TestFilterBit(16) || track->Chi2perNDF() < 0.1) continue;
1492 Double_t b[2] = {-99., -99.};
1493 Double_t bCov[3] = {-99., -99., -99.};
1494 AliAODTrack copy(*track);
1495 if (copy.PropagateToDCA(event->GetPrimaryVertex(), event->GetMagneticField(), 100., b, bCov) && TMath::Abs(b[0]) < 0.3 && TMath::Abs(b[1]) < 0.3) multGlob++;
1497 if(year == 2010 && multTPC > (-40.3+1.22*multGlob) && multTPC < (32.1+1.59*multGlob)) return kTRUE;
1498 if(year == 2011 && multTPC > (-36.73 + 1.48*multGlob) && multTPC < (62.87 + 1.78*multGlob)) return kTRUE;
1501 //_____________________________________________________________________________
1502 void AliAnalysisTaskJetV2::FillHistogramsAfterSubtraction(Double_t psi2, Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc)
1505 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1506 FillTrackHistograms();
1507 if(fAnalysisType == AliAnalysisTaskJetV2::kFull) FillClusterHistograms();
1508 FillJetHistograms(psi2);
1509 if(fFillQAHistograms) FillEventPlaneHistograms(vzero, vzeroComb, tpc);
1510 FillRhoHistograms();
1511 FillDeltaPtHistograms(psi2);
1513 //_____________________________________________________________________________
1514 void AliAnalysisTaskJetV2::FillTrackHistograms() const
1516 // fill track histograms
1517 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1518 Int_t iTracks(fTracks->GetEntriesFast()), iAcceptedTracks(0);
1519 for(Int_t i(0); i < iTracks; i++) {
1520 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1521 if(!PassesCuts(track)) continue;
1523 fHistPicoTrackPt[fInCentralitySelection]->Fill(track->Pt());
1524 if(fFillQAHistograms) FillQAHistograms(track);
1526 fHistPicoTrackMult[fInCentralitySelection]->Fill(iAcceptedTracks);
1528 //_____________________________________________________________________________
1529 void AliAnalysisTaskJetV2::FillClusterHistograms() const
1531 // fill cluster histograms
1532 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1533 if(!fClusterCont) return;
1534 Int_t iClusters(fClusterCont->GetNClusters());
1535 for(Int_t i(0); i < iClusters; i++) {
1536 AliVCluster* cluster = fClusterCont->GetCluster(i);
1537 if (!PassesCuts(cluster)) continue;
1538 TLorentzVector clusterLorentzVector;
1539 cluster->GetMomentum(clusterLorentzVector, const_cast<Double_t*>(fVertex));
1540 fHistClusterPt[fInCentralitySelection]->Fill(clusterLorentzVector.Pt());
1541 fHistClusterEtaPhi[fInCentralitySelection]->Fill(clusterLorentzVector.Eta(), clusterLorentzVector.Phi());
1542 fHistClusterEtaPhiWeighted[fInCentralitySelection]->Fill(clusterLorentzVector.Eta(), clusterLorentzVector.Phi(), clusterLorentzVector.Pt());
1546 //_____________________________________________________________________________
1547 void AliAnalysisTaskJetV2::FillEventPlaneHistograms(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc) const
1549 // fill event plane histograms
1550 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1551 fHistPsiControl->Fill(0.5, vzero[0][0]); // vzero a psi2
1552 fHistPsiControl->Fill(1.5, vzero[1][0]); // vzero c psi2
1553 fHistPsiControl->Fill(2.5, tpc[0]); // tpc psi 2
1554 fHistPsiControl->Fill(5.5, vzero[0][1]); // vzero a psi3
1555 fHistPsiControl->Fill(6.5, vzero[1][1]); // vzero b psi3
1556 fHistPsiControl->Fill(7.5, tpc[1]); // tpc psi 3
1557 fHistPsiVZEROA->Fill(vzero[0][0]);
1558 fHistPsiVZEROC->Fill(vzero[1][0]);
1559 fHistPsiVZERO->Fill(vzeroComb[0]);
1560 fHistPsiTPC->Fill(tpc[0]);
1561 fHistPsiSpread->Fill(0.5, TMath::Abs(vzero[0][0]-vzero[1][0]));
1562 fHistPsiSpread->Fill(1.5, TMath::Abs(vzero[0][0]-tpc[0]));
1563 fHistPsiSpread->Fill(2.5, TMath::Abs(vzero[1][0]-tpc[0]));
1564 // event plane vs centrality QA histo's to check recentering
1565 Double_t TRK(InputEvent()->GetCentrality()->GetCentralityPercentile("TRK"));
1566 Double_t V0M(InputEvent()->GetCentrality()->GetCentralityPercentile("V0M"));
1567 fHistPsiVZEROAV0M->Fill(V0M, vzero[0][0]);
1568 fHistPsiVZEROCV0M->Fill(V0M, vzero[1][0]);
1569 fHistPsiVZEROVV0M->Fill(V0M, vzeroComb[0]);
1570 fHistPsiTPCiV0M->Fill(V0M, tpc[0]);
1571 fHistPsiVZEROATRK->Fill(TRK, vzero[0][0]);
1572 fHistPsiVZEROCTRK->Fill(TRK, vzero[1][0]);
1573 fHistPsiVZEROTRK->Fill(TRK, vzeroComb[0]);
1574 fHistPsiTPCTRK->Fill(TRK, tpc[0]);
1576 //_____________________________________________________________________________
1577 void AliAnalysisTaskJetV2::FillRhoHistograms()
1579 // fill rho histograms
1580 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1581 fHistRhoPackage[fInCentralitySelection]->Fill(fLocalRho->GetVal()); // save the rho estimate from the emcal jet package
1582 // get multiplicity FIXME inefficient
1583 Int_t iJets(fJets->GetEntriesFast());
1584 Double_t rho(fLocalRho->GetLocalVal(TMath::Pi(), TMath::Pi(), fLocalRho->GetVal()));
1585 fHistRho[fInCentralitySelection]->Fill(rho);
1586 fHistRhoVsMult->Fill(fTracks->GetEntries(), rho);
1587 fHistRhoVsCent->Fill(fCent, rho);
1588 for(Int_t i(0); i < iJets; i++) {
1589 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
1590 if(!PassesCuts(jet)) continue;
1591 fHistRhoAVsMult->Fill(fTracks->GetEntries(), rho * jet->Area());
1592 fHistRhoAVsCent->Fill(fCent, rho * jet->Area());
1595 //_____________________________________________________________________________
1596 void AliAnalysisTaskJetV2::FillDeltaPtHistograms(Double_t psi2) const
1598 // fill delta pt histograms
1599 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1601 const Float_t areaRC = GetJetRadius()*GetJetRadius()*TMath::Pi();
1602 // we're retrieved the leading jet, now get a random cone
1603 for(i = 0; i < fMaxCones; i++) {
1604 Float_t pt(0), eta(0), phi(0);
1605 // get a random cone without constraints on leading jet position
1606 CalculateRandomCone(pt, eta, phi, fTracksCont, fClusterCont, 0x0);
1608 if(fFillQAHistograms) fHistRCPhiEta[fInCentralitySelection]->Fill(phi, eta);
1609 fHistRhoVsRCPt[fInCentralitySelection]->Fill(pt, fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal())*areaRC);
1610 fHistRCPt[fInCentralitySelection]->Fill(pt);
1611 fHistDeltaPtDeltaPhi2[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1612 fHistDeltaPtDeltaPhi2Rho0[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetVal());
1615 // get a random cone excluding leading jet area
1616 CalculateRandomCone(pt, eta, phi, fTracksCont, fClusterCont, fLeadingJet);
1618 if(fFillQAHistograms) fHistRCPhiEtaExLJ[fInCentralitySelection]->Fill(phi, eta);
1619 fHistRhoVsRCPtExLJ[fInCentralitySelection]->Fill(pt, fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal())*areaRC);
1620 fHistRCPtExLJ[fInCentralitySelection]->Fill(pt);
1621 fHistDeltaPtDeltaPhi2ExLJ[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1622 fHistDeltaPtDeltaPhi2ExLJRho0[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetVal());
1626 //_____________________________________________________________________________
1627 void AliAnalysisTaskJetV2::FillJetHistograms(Double_t psi2)
1629 // fill jet histograms
1630 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1631 Int_t iJets(fJets->GetEntriesFast());
1632 for(Int_t i(0); i < iJets; i++) {
1633 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
1634 if(PassesCuts(jet)) {
1635 Double_t pt(jet->Pt()), area(jet->Area()), eta(jet->Eta()), phi(jet->Phi());
1636 Double_t rho(fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1637 fHistJetPtRaw[fInCentralitySelection]->Fill(pt);
1638 fHistJetPt[fInCentralitySelection]->Fill(pt-area*rho);
1639 if(fFillQAHistograms) fHistJetEtaPhi[fInCentralitySelection]->Fill(eta, phi);
1640 fHistJetPtArea[fInCentralitySelection]->Fill(pt-area*rho, area);
1641 fHistJetPtEta[fInCentralitySelection]->Fill(pt-area*rho, eta);
1642 fHistJetPsi2Pt[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt-area*rho);
1643 fHistJetPsi2PtRho0[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt-area*fLocalRho->GetVal());
1644 fHistJetPtConstituents[fInCentralitySelection]->Fill(pt-area*rho, jet->Nch());
1645 fHistJetEtaRho[fInCentralitySelection]->Fill(eta, pt/area);
1649 //_____________________________________________________________________________
1650 void AliAnalysisTaskJetV2::FillQAHistograms(AliVTrack* vtrack) const
1652 // fill qa histograms for pico tracks
1654 AliPicoTrack* track = static_cast<AliPicoTrack*>(vtrack);
1655 fHistRunnumbersPhi->Fill(fMappedRunNumber, track->Phi());
1656 fHistRunnumbersEta->Fill(fMappedRunNumber, track->Eta());
1657 Int_t type((int)(track->GetTrackType()));
1660 fHistPicoCat1[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1663 fHistPicoCat2[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1666 fHistPicoCat3[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1671 //_____________________________________________________________________________
1672 void AliAnalysisTaskJetV2::FillQAHistograms(AliVEvent* vevent)
1674 // fill qa histograms for events
1676 fHistVertexz->Fill(vevent->GetPrimaryVertex()->GetZ());
1677 fHistCentrality->Fill(fCent);
1678 Int_t runNumber(InputEvent()->GetRunNumber());
1679 for(fMappedRunNumber = 0; fExpectedRuns->GetSize()+1; fMappedRunNumber++) {
1680 if(fExpectedRuns->At(fMappedRunNumber) == runNumber) break;
1683 //_____________________________________________________________________________
1684 void AliAnalysisTaskJetV2::FillAnalysisSummaryHistogram() const
1686 // fill the analysis summary histrogram, saves all relevant analysis settigns
1687 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1688 fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fJetRadius");
1689 fHistAnalysisSummary->SetBinContent(2, GetJetContainer()->GetJetRadius());
1690 fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fJetEtaMin");
1691 fHistAnalysisSummary->SetBinContent(3, GetJetContainer()->GetJetEtaMin());
1692 fHistAnalysisSummary->GetXaxis()->SetBinLabel(4, "fJetEtaMax");
1693 fHistAnalysisSummary->SetBinContent(4, GetJetContainer()->GetJetEtaMax());
1694 fHistAnalysisSummary->GetXaxis()->SetBinLabel(5, "fJetPhiMin");
1695 fHistAnalysisSummary->SetBinContent(5, GetJetContainer()->GetJetPhiMin());
1696 fHistAnalysisSummary->GetXaxis()->SetBinLabel(6, "fJetPhiMax");
1697 fHistAnalysisSummary->SetBinContent(6, GetJetContainer()->GetJetPhiMin());
1698 fHistAnalysisSummary->GetXaxis()->SetBinLabel(16, "fForceBeamType");
1699 fHistAnalysisSummary->SetBinContent(16, fForceBeamType);
1700 fHistAnalysisSummary->GetXaxis()->SetBinLabel(17, "fMinCent");
1701 fHistAnalysisSummary->SetBinContent(17, fMinCent);
1702 fHistAnalysisSummary->GetXaxis()->SetBinLabel(18, "fMaxCent");
1703 fHistAnalysisSummary->SetBinContent(18, fMaxCent);
1704 fHistAnalysisSummary->GetXaxis()->SetBinLabel(19, "fMinVz");
1705 fHistAnalysisSummary->SetBinContent(19, fMinVz);
1706 fHistAnalysisSummary->GetXaxis()->SetBinLabel(20, "fMaxVz");
1707 fHistAnalysisSummary->SetBinContent(20, fMaxVz);
1708 fHistAnalysisSummary->GetXaxis()->SetBinLabel(21, "fOffTrigger");
1709 fHistAnalysisSummary->SetBinContent(21, fOffTrigger);
1710 fHistAnalysisSummary->GetXaxis()->SetBinLabel(34, "fitModulationType");
1711 fHistAnalysisSummary->SetBinContent(34, (int)fFitModulationType);
1712 fHistAnalysisSummary->GetXaxis()->SetBinLabel(35, "runModeType");
1713 fHistAnalysisSummary->SetBinContent(35, (int)fRunModeType);
1714 fHistAnalysisSummary->GetXaxis()->SetBinLabel(36, "data type");
1715 fHistAnalysisSummary->SetBinContent(36, (int)fDataType);
1716 fHistAnalysisSummary->GetXaxis()->SetBinLabel(37, "iterator");
1717 fHistAnalysisSummary->SetBinContent(37, 1.);
1718 fHistAnalysisSummary->GetXaxis()->SetBinLabel(38, "fMinPvalue");
1719 fHistAnalysisSummary->SetBinContent(38, fMinPvalue);
1720 fHistAnalysisSummary->GetXaxis()->SetBinLabel(39, "fMaxPvalue");
1721 fHistAnalysisSummary->SetBinContent(39, fMaxPvalue);
1722 fHistAnalysisSummary->GetXaxis()->SetBinLabel(40, "fExcludeLeadingJetsFromFit");
1723 fHistAnalysisSummary->SetBinContent(40, fExcludeLeadingJetsFromFit);
1724 fHistAnalysisSummary->GetXaxis()->SetBinLabel(41, "fRebinSwapHistoOnTheFly");
1725 fHistAnalysisSummary->SetBinContent(41, (int)fRebinSwapHistoOnTheFly);
1726 fHistAnalysisSummary->GetXaxis()->SetBinLabel(42, "fUsePtWeight");
1727 fHistAnalysisSummary->SetBinContent(42, (int)fUsePtWeight);
1728 // fHistAnalysisSummary->GetXaxis()->SetBinLabel(43, "fExplicitOutlierCut");
1729 // fHistAnalysisSummary->SetBinContent(43, fExplicitOutlierCut);
1730 fHistAnalysisSummary->GetXaxis()->SetBinLabel(44, "fSoftTrackMinPt");
1731 fHistAnalysisSummary->SetBinContent(44, fSoftTrackMinPt);
1732 fHistAnalysisSummary->GetXaxis()->SetBinLabel(45, "fSoftTrackMaxPt");
1733 fHistAnalysisSummary->SetBinContent(45, fSoftTrackMaxPt);
1734 fHistAnalysisSummary->GetXaxis()->SetBinLabel(46, "fMaxCones");
1735 fHistAnalysisSummary->SetBinContent(46, fMaxCones);
1736 fHistAnalysisSummary->GetXaxis()->SetBinLabel(47, "used rho");
1737 fHistAnalysisSummary->GetXaxis()->SetBinLabel(48, "used small rho");
1739 //_____________________________________________________________________________
1740 void AliAnalysisTaskJetV2::Terminate(Option_t *)
1743 switch (fRunModeType) {
1745 printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1746 AliAnalysisTaskJetV2::Dump();
1747 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));
1752 //_____________________________________________________________________________
1753 void AliAnalysisTaskJetV2::SetModulationFit(TF1* fit)
1755 // set modulation fit
1756 if (fFitModulation) delete fFitModulation;
1757 fFitModulation = fit;
1759 //_____________________________________________________________________________
1760 void AliAnalysisTaskJetV2::SetUseControlFit(Bool_t c)
1763 if (fFitControl) delete fFitControl;
1765 fFitControl = new TF1("controlFit", "pol0", 0, TMath::TwoPi());
1766 } else fFitControl = 0x0;
1768 //_____________________________________________________________________________
1769 TH1F* AliAnalysisTaskJetV2::GetResolutionFromOuptutFile(detectorType det, Int_t h, TArrayD* cen)
1771 // INTERFACE METHOD FOR OUTPUTFILE
1772 // get the detector resolution, user has ownership of the returned histogram
1774 printf(" > Please add fOutputList first < \n");
1778 (cen) ? r = new TH1F("R", "R", cen->GetSize()-1, cen->GetArray()) : r = new TH1F("R", "R", 10, 0, 10);
1779 if(!cen) r->GetXaxis()->SetTitle("number of centrality bin");
1780 r->GetYaxis()->SetTitle(Form("Resolution #Psi_{%i}", h));
1781 for(Int_t i(0); i < 10; i++) {
1782 TProfile* temp((TProfile*)fOutputList->FindObject(Form("fProfV%iResolution_%i", h, i)));
1784 Double_t a(temp->GetBinContent(3)), b(temp->GetBinContent(5)), c(temp->GetBinContent(7));
1785 Double_t d(temp->GetBinContent(9)), e(temp->GetBinContent(10)), f(temp->GetBinContent(11));
1786 Double_t _a(temp->GetBinError(3)), _b(temp->GetBinError(5)), _c(temp->GetBinError(7));
1787 Double_t _d(temp->GetBinError(9)), _e(temp->GetBinError(10)), _f(temp->GetBinError(11));
1788 if(a <= 0 || b <= 0 || c <= 0 || d <= 0 || e <= 0 || f <= 0) continue;
1791 r->SetBinContent(1+i, TMath::Sqrt((a*b)/c));
1792 if(i==0) r->SetNameTitle("VZEROA resolution", "VZEROA resolution");
1793 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
1796 r->SetBinContent(1+i, TMath::Sqrt((a*c)/b));
1797 if(i==0) r->SetNameTitle("VZEROC resolution", "VZEROC resolution");
1798 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
1801 r->SetBinContent(1+i, TMath::Sqrt((b*c)/a));
1802 if(i==0) r->SetNameTitle("TPC resolution", "TPC resolution");
1803 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
1806 r->SetBinContent(1+i, TMath::Sqrt((d*e)/f));
1807 if(i==0) r->SetNameTitle("VZEROComb resolution", "VZEROComb resolution");
1808 r->SetBinError(1+i, TMath::Sqrt(_d*_d+_e*_e+_f*_f));
1815 //_____________________________________________________________________________
1816 TH1F* AliAnalysisTaskJetV2::CorrectForResolutionDiff(TH1F* v, detectorType det, TArrayD* cen, Int_t c, Int_t h)
1818 // INTERFACE METHOD FOR OUTPUT FILE
1819 // correct the supplied differential vn histogram v for detector resolution
1820 TH1F* r(GetResolutionFromOuptutFile(det, h, cen));
1822 printf(" > Couldn't find resolution < \n");
1825 Double_t res(1./r->GetBinContent(1+r->FindBin(c)));
1826 TF1* line = new TF1("line", "pol0", 0, 200);
1827 line->SetParameter(0, res);
1831 //_____________________________________________________________________________
1832 TH1F* AliAnalysisTaskJetV2::CorrectForResolutionInt(TH1F* v, detectorType det, TArrayD* cen, Int_t h)
1834 // INTERFACE METHOD FOR OUTPUT FILE
1835 // correct the supplied intetrated vn histogram v for detector resolution
1836 // integrated vn must have the same centrality binning as the resolotion correction
1837 TH1F* r(GetResolutionFromOuptutFile(det, h, cen));
1841 //_____________________________________________________________________________
1842 TH1F* AliAnalysisTaskJetV2::GetDifferentialQC(TProfile* refCumulants, TProfile* diffCumlants, TArrayD* ptBins, Int_t h)
1844 // get differential QC
1845 Double_t r(refCumulants->GetBinContent(h-1)); // v2 reference flow
1846 if(r > 0) r = TMath::Sqrt(r);
1847 TH1F* qc = new TH1F(Form("QC2v%i", h), Form("QC2v%i", h), ptBins->GetSize()-1, ptBins->GetArray());
1848 Double_t a(0), b(0), c(0); // dummy variables
1849 for(Int_t i(0); i < ptBins->GetSize(); i++) {
1851 a = diffCumlants->GetBinContent(1+i);
1852 b = diffCumlants->GetBinError(1+i);
1854 qc->SetBinContent(1+i, c);
1855 (a <= 0 || b <= 0) ? qc->SetBinError(1+i, b) : qc->SetBinError(1+i, TMath::Sqrt(c*c*b*b/(a*a)));
1861 //_____________________________________________________________________________