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;
161 if(fOutputListGood) delete fOutputListGood;
162 if(fOutputListBad) delete fOutputListBad;
163 if(fFitModulation) delete fFitModulation;
164 if(fHistSwap) delete fHistSwap;
165 if(fCentralityClasses) delete fCentralityClasses;
166 if(fExpectedRuns) delete fExpectedRuns;
167 if(fExpectedSemiGoodRuns) delete fExpectedSemiGoodRuns;
168 if(fFitControl) delete fFitControl;
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::InitializeAnalysis()
189 // initialize the anaysis
190 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
191 // if not set, estimate the number of cones that would fit into the selected acceptance
192 if(fMaxCones <= 0) fMaxCones = TMath::CeilNint((TMath::Abs(GetJetContainer()->GetJetEtaMax()-GetJetContainer()->GetJetEtaMin())*TMath::Abs(GetJetContainer()->GetJetPhiMax()-GetJetContainer()->GetJetPhiMin()))/(TMath::Pi()*GetJetRadius()*GetJetRadius()));
193 // manually 'override' the default acceptance cuts of the emcal framework (use with caution)
194 if(fMinDisanceRCtoLJ==0) fMinDisanceRCtoLJ = GetJetRadius();
195 if(dynamic_cast<AliAODEvent*>(InputEvent())) fDataType = kAOD; // determine the datatype
196 else if(dynamic_cast<AliESDEvent*>(InputEvent())) fDataType = kESD;
197 fHistAnalysisSummary->SetBinContent(36, (int)fDataType);
198 if(!fRandom) fRandom = new TRandom3(0); // set randomizer and random seed
199 switch (fFitModulationType) {
200 case kNoFit : { SetModulationFit(new TF1("fix_kNoFit", "[0]", 0, TMath::TwoPi())); } break;
202 SetModulationFit(new TF1("fit_kV2", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
203 fFitModulation->SetParameter(0, 0.); // normalization
204 fFitModulation->SetParameter(3, 0.2); // v2
205 fFitModulation->FixParameter(1, 1.); // constant
206 fFitModulation->FixParameter(2, 2.); // constant
209 SetModulationFit(new TF1("fit_kV3", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
210 fFitModulation->SetParameter(0, 0.); // normalization
211 fFitModulation->SetParameter(3, 0.2); // v3
212 fFitModulation->FixParameter(1, 1.); // constant
213 fFitModulation->FixParameter(2, 3.); // constant
215 default : { // for the combined fit, the 'direct fourier series' or the user supplied vn values we use v2 and v3
216 SetModulationFit(new TF1("fit_kCombined", "[0]*([1]+[2]*([3]*TMath::Cos([2]*(x-[4]))+[7]*TMath::Cos([5]*(x-[6]))))", 0, TMath::TwoPi()));
217 fFitModulation->SetParameter(0, 0.); // normalization
218 fFitModulation->SetParameter(3, 0.2); // v2
219 fFitModulation->FixParameter(1, 1.); // constant
220 fFitModulation->FixParameter(2, 2.); // constant
221 fFitModulation->FixParameter(5, 3.); // constant
222 fFitModulation->SetParameter(7, 0.2); // v3
225 switch (fRunModeType) {
226 case kGrid : { fFitModulationOptions += "N0"; } break;
229 FillAnalysisSummaryHistogram();
232 //_____________________________________________________________________________
233 TH1F* AliAnalysisTaskJetV2::BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Int_t c, Bool_t append)
235 // book a TH1F and connect it to the output container
236 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
237 if(fReduceBinsXByFactor > 0 ) bins = TMath::Nint(bins/fReduceBinsXByFactor);
238 if(!fOutputList) return 0x0;
240 if(c!=-1) { // format centrality dependent histograms accordingly
241 name = Form("%s_%i", name, c);
242 title += Form("_%i-%i", (int)(fCentralityClasses->At(c)), (int)(fCentralityClasses->At((1+c))));
244 title += Form(";%s;[counts]", x);
245 TH1F* histogram = new TH1F(name, title.Data(), bins, min, max);
247 if(append) fOutputList->Add(histogram);
250 //_____________________________________________________________________________
251 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)
253 // book a TH2F and connect it to the output container
254 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
255 if(fReduceBinsXByFactor > 0 ) binsx = TMath::Nint(binsx/fReduceBinsXByFactor);
256 if(fReduceBinsYByFactor > 0 ) binsy = TMath::Nint(binsy/fReduceBinsYByFactor);
257 if(!fOutputList) return 0x0;
259 if(c!=-1) { // format centrality dependent histograms accordingly
260 name = Form("%s_%i", name, c);
261 title += Form("_%i-%i", (int)fCentralityClasses->At(c), (int)(fCentralityClasses->At((1+c))));
263 title += Form(";%s;%s", x, y);
264 TH2F* histogram = new TH2F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy);
266 if(append) fOutputList->Add(histogram);
269 //_____________________________________________________________________________
270 void AliAnalysisTaskJetV2::UserCreateOutputObjects()
272 // create output objects. also initializes some default values in case they aren't
273 // loaded via the AddTask macro
274 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
275 fOutputList = new TList();
276 fOutputList->SetOwner(kTRUE);
277 if(!fCentralityClasses) { // classes must be defined at this point
278 Double_t c[] = {0., 20., 40., 60., 80., 100.};
279 fCentralityClasses = new TArrayD(sizeof(c)/sizeof(c[0]), c);
281 if(!fExpectedRuns) { // expected runs must be defined at this point
282 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 */};
283 fExpectedRuns = new TArrayI(sizeof(r)/sizeof(r[0]), r);
285 if(!fExpectedSemiGoodRuns) {
286 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};
287 fExpectedSemiGoodRuns = new TArrayI(sizeof(r)/sizeof(r[0]), r);
290 fHistCentrality = BookTH1F("fHistCentrality", "centrality", 102, -2, 100);
291 fHistVertexz = BookTH1F("fHistVertexz", "vertex z (cm)", 100, -12, 12);
293 // pico track and emcal cluster kinematics
294 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
295 fHistPicoTrackPt[i] = BookTH1F("fHistPicoTrackPt", "p_{t} [GeV/c]", 100, 0, 100, i);
296 fHistPicoTrackMult[i] = BookTH1F("fHistPicoTrackMult", "multiplicity", 100, 0, 5000, i);
297 if(fFillQAHistograms) {
298 fHistPicoCat1[i] = BookTH2F("fHistPicoCat1", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
299 fHistPicoCat2[i] = BookTH2F("fHistPicoCat2", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
300 fHistPicoCat3[i] = BookTH2F("fHistPicoCat3", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
301 if(fAnalysisType == AliAnalysisTaskJetV2::kFull) {
302 fHistClusterPt[i] = BookTH1F("fHistClusterPt", "p_{t} [GeV/c]", 100, 0, 100, i);
303 fHistClusterEtaPhi[i] = BookTH2F("fHistClusterEtaPhi", "#eta", "#phi", 100, -1., 1., 100, 0, TMath::TwoPi(), i);
304 fHistClusterEtaPhiWeighted[i] = BookTH2F("fHistClusterEtaPhiWeighted", "#eta", "#phi", 100, -1., 1., 100, 0, TMath::TwoPi(), i);
309 if(fFillQAHistograms) {
310 // event plane estimates and quality
311 fHistPsiControl = new TProfile("fHistPsiControl", "fHistPsiControl", 10, 0, 10);
312 fHistPsiControl->Sumw2();
313 fHistPsiSpread = new TProfile("fHistPsiSpread", "fHistPsiSpread", 4, 0, 4);
314 fHistPsiSpread->Sumw2();
315 fHistPsiControl->GetXaxis()->SetBinLabel(1, "<#Psi_{2, VZEROA}>");
316 fHistPsiControl->GetXaxis()->SetBinLabel(2, "<#Psi_{2, VZEROC}>");
317 fHistPsiControl->GetXaxis()->SetBinLabel(3, "<#Psi_{2, TPC}>");
318 fHistPsiControl->GetXaxis()->SetBinLabel(4, "<#Psi_{2, TPC, #eta < 0}>");
319 fHistPsiControl->GetXaxis()->SetBinLabel(5, "<#Psi_{2, TPC, #eta > 0}>");
320 fHistPsiControl->GetXaxis()->SetBinLabel(6, "<#Psi_{3, VZEROA}>");
321 fHistPsiControl->GetXaxis()->SetBinLabel(7, "<#Psi_{3, VZEROC}>");
322 fHistPsiControl->GetXaxis()->SetBinLabel(8, "<#Psi_{3, TPC}>");
323 fHistPsiControl->GetXaxis()->SetBinLabel(9, "<#Psi_{3, TPC, #eta < 0}>");
324 fHistPsiControl->GetXaxis()->SetBinLabel(10, "<#Psi_{3, TPC, #eta > 0}>");
325 fHistPsiSpread->GetXaxis()->SetBinLabel(1, "<#Psi_{2, VZEROA} - #Psi_{2, VZEROC}>");
326 fHistPsiSpread->GetXaxis()->SetBinLabel(2, "<#Psi_{2, VZEROC} - #Psi_{2, TPC}>");
327 fHistPsiSpread->GetXaxis()->SetBinLabel(3, "<#Psi_{2, VZEROC} - #Psi_{2, TPC}>");
328 fHistPsiSpread->GetXaxis()->SetBinLabel(4, "<#Psi_{2, TPC, #eta < 0} - #Psi_{2, TPC, #eta > 0}>");
329 fOutputList->Add(fHistPsiControl);
330 fOutputList->Add(fHistPsiSpread);
331 fHistPsiVZEROA = BookTH1F("fHistPsiVZEROA", "#Psi_{VZEROA}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
332 fHistPsiVZEROC = BookTH1F("fHistPsiVZEROC", "#Psi_{VZEROC}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
333 fHistPsiVZERO = BookTH1F("fHistPsiVZERO", "#Psi_{VZERO}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
334 fHistPsiTPC = BookTH1F("fHistPsiTPC", "#Psi_{TPC}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
335 fHistPsiVZEROAV0M = BookTH2F("fHistPsiVZEROAV0M", "V0M", "#Psi_{2, VZEROA}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
336 fHistPsiVZEROCV0M = BookTH2F("fHistPsiVZEROCV0M", "V0M", "#Psi_{2, VZEROC}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
337 fHistPsiVZEROVV0M = BookTH2F("fHistPsiVZEROV0M", "V0M", "#Psi_{2, VZERO}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
338 fHistPsiTPCiV0M = BookTH2F("fHistPsiTPCV0M", "V0M", "#Psi_{2, TRK}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
339 fHistPsiVZEROATRK = BookTH2F("fHistPsiVZEROATRK", "TRK", "#Psi_{2, VZEROA}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
340 fHistPsiVZEROCTRK = BookTH2F("fHistPsiVZEROCTRK", "TRK", "#Psi_{2, VZEROC}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
341 fHistPsiVZEROTRK = BookTH2F("fHistPsiVZEROTRK", "TRK", "#Psi_{2, VZERO}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
342 fHistPsiTPCTRK = BookTH2F("fHistPsiTPCTRK", "TRK", "#Psi_{2, TRK}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
345 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i ++) {
346 fHistRhoPackage[i] = BookTH1F("fHistRhoPackage", "#rho [GeV/c]", 100, 0, 150, i);
347 fHistRho[i] = BookTH1F("fHistRho", "#rho [GeV/c]", 100, 0, 150, i);
349 fHistRhoVsMult = BookTH2F("fHistRhoVsMult", "multiplicity", "#rho [GeV/c]", 100, 0, 4000, 100, 0, 250);
350 fHistRhoVsCent = BookTH2F("fHistRhoVsCent", "centrality", "#rho [GeV/c]", 100, 0, 100, 100, 0, 250);
351 fHistRhoAVsMult = BookTH2F("fHistRhoAVsMult", "multiplicity", "#rho * A (jet) [GeV/c]", 100, 0, 4000, 100, 0, 50);
352 fHistRhoAVsCent = BookTH2F("fHistRhoAVsCent", "centrality", "#rho * A (jet) [GeV/c]", 100, 0, 100, 100, 0, 50);
354 TString detector("");
355 switch (fDetectorType) {
356 case kTPC : detector+="TPC";
358 case kVZEROA : detector+="VZEROA";
360 case kVZEROC : detector+="VZEROC";
362 case kVZEROComb : detector+="VZEROComb";
366 // delta pt distributions
367 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i ++) {
368 if(fFillQAHistograms) fHistRCPhiEta[i] = BookTH2F("fHistRCPhiEta", "#phi (RC)", "#eta (RC)", 40, 0, TMath::TwoPi(), 40, -1, 1, i);
369 fHistRhoVsRCPt[i] = BookTH2F("fHistRhoVsRCPt", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i);
370 fHistRCPt[i] = BookTH1F("fHistRCPt", "p_{t} (RC) [GeV/c]", 130, -20, 150, i);
371 if(fFillQAHistograms) fHistRCPhiEtaExLJ[i] = BookTH2F("fHistRCPhiEtaExLJ", "#phi (RC)", "#eta (RC)", 40, 0, TMath::TwoPi(), 40, -1, 1, i);
372 fHistDeltaPtDeltaPhi2[i] = BookTH2F("fHistDeltaPtDeltaPhi2", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
373 fHistDeltaPtDeltaPhi2Rho0[i] = BookTH2F("fHistDeltaPtDeltaPhi2Rho0", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
374 fHistRhoVsRCPtExLJ[i] = BookTH2F("fHistRhoVsRCPtExLJ", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i);
375 fHistRCPtExLJ[i] = BookTH1F("fHistRCPtExLJ", "p_{t} (RC) [GeV/c]", 130, -20, 150, i);
376 fHistDeltaPtDeltaPhi2ExLJ[i] = BookTH2F("fHistDeltaPtDeltaPhi2ExLJ", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
377 fHistDeltaPtDeltaPhi2ExLJRho0[i] = BookTH2F("fHistDeltaPtDeltaPhi2ExLJRho0", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
378 // jet histograms (after kinematic cuts)
379 fHistJetPtRaw[i] = BookTH1F("fHistJetPtRaw", "p_{t, jet} RAW [GeV/c]", 200, -50, 150, i);
380 fHistJetPt[i] = BookTH1F("fHistJetPt", "p_{t, jet} [GeV/c]", 350, -100, 250, i);
381 if(fFillQAHistograms) fHistJetEtaPhi[i] = BookTH2F("fHistJetEtaPhi", "#eta", "#phi", 100, -1, 1, 100, 0, TMath::TwoPi(), i);
382 fHistJetPtArea[i] = BookTH2F("fHistJetPtArea", "p_{t, jet} [GeV/c]", "Area", 175, -100, 250, 30, 0, 0.9, i);
383 fHistJetPtEta[i] = BookTH2F("fHistJetPtEta", "p_{t, jet} [GeV/c]", "Eta", 175, -100, 250, 30, -0.9, 0.9, i);
384 fHistJetPtConstituents[i] = BookTH2F("fHistJetPtConstituents", "p_{t, jet} [GeV/c]", "Area", 350, -100, 250, 60, 0, 150, i);
385 fHistJetEtaRho[i] = BookTH2F("fHistJetEtaRho", "#eta", "#rho", 100, -1, 1, 100, 0, 300, i);
386 // in plane and out of plane spectra
387 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);
388 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);
389 // profiles for all correlator permutations which are necessary to calculate each second and third order event plane resolution
390 fProfV2Resolution[i] = new TProfile(Form("fProfV2Resolution_%i", i), Form("fProfV2Resolution_%i", i), 11, -0.5, 10.5);
391 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(3, "<cos(2(#Psi_{VZEROA} - #Psi_{VZEROC}))>");
392 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(4, "<cos(2(#Psi_{VZEROC} - #Psi_{VZEROA}))>");
393 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(5, "<cos(2(#Psi_{VZEROA} - #Psi_{TPC}))>");
394 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(6, "<cos(2(#Psi_{TPC} - #Psi_{VZEROA}))>");
395 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(7, "<cos(2(#Psi_{VZEROC} - #Psi_{TPC}))>");
396 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(8, "<cos(2(#Psi_{TPC} - #Psi_{VZEROC}))>");
397 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(9, "<cos(2(#Psi_{VZERO} - #Psi_{TPC_A}))>");
398 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(10, "<cos(2(#Psi_{VZERO} - #Psi_{TPC_B}))>");
399 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(11, "<cos(2(#Psi_{TPC_A} - #Psi_{TPC_B}))>");
400 fOutputList->Add(fProfV2Resolution[i]);
401 fProfV3Resolution[i] = new TProfile(Form("fProfV3Resolution_%i", i), Form("fProfV3Resolution_%i", i), 11, -0.5, 10.5);
402 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(3, "<cos(3(#Psi_{VZEROA} - #Psi_{VZEROC}))>");
403 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(4, "<cos(3(#Psi_{VZEROC} - #Psi_{VZEROA}))>");
404 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(5, "<cos(3(#Psi_{VZEROA} - #Psi_{TPC}))>");
405 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(6, "<cos(3(#Psi_{TPC} - #Psi_{VZEROA}))>");
406 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(7, "<cos(3(#Psi_{VZEROC} - #Psi_{TPC}))>");
407 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(8, "<cos(3(#Psi_{TPC} - #Psi_{VZEROC}))>");
408 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(9, "<cos(3(#Psi_{VZERO} - #Psi_{TPC_A}))>");
409 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(10, "<cos(3(#Psi_{VZERO} - #Psi_{TPC_B}))>");
410 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(11, "<cos(3(#Psi_{TPC_A} - #Psi_{TPC_B}))>");
411 fOutputList->Add(fProfV3Resolution[i]);
414 Float_t temp[fCentralityClasses->GetSize()];
415 for(Int_t i(0); i < fCentralityClasses->GetSize(); i++) temp[i] = fCentralityClasses->At(i);
416 fProfV2 = new TProfile("fProfV2", "fProfV2", fCentralityClasses->GetSize()-1, temp);
417 fProfV3 = new TProfile("fProfV3", "fProfV3", fCentralityClasses->GetSize()-1, temp);
418 fOutputList->Add(fProfV2);
419 fOutputList->Add(fProfV3);
420 switch (fFitModulationType) {
422 fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
423 fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
424 fOutputList->Add(fProfV2Cumulant);
425 fOutputList->Add(fProfV3Cumulant);
428 fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
429 fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
430 fOutputList->Add(fProfV2Cumulant);
431 fOutputList->Add(fProfV3Cumulant);
435 // for the histograms initialized below, binning is fixed to runnumbers or flags
436 fReduceBinsXByFactor = 1;
437 fReduceBinsYByFactor = 1;
438 if(fFillQAHistograms) {
439 fHistRunnumbersEta = new TH2F("fHistRunnumbersEta", "fHistRunnumbersEta", fExpectedRuns->GetSize()+1, -.5, fExpectedRuns->GetSize()+.5, 100, -1.1, 1.1);
440 fHistRunnumbersEta->Sumw2();
441 fOutputList->Add(fHistRunnumbersEta);
442 fHistRunnumbersPhi = new TH2F("fHistRunnumbersPhi", "fHistRunnumbersPhi", fExpectedRuns->GetSize()+1, -.5, fExpectedRuns->GetSize()+.5, 100, -0.2, TMath::TwoPi()+0.2);
443 fHistRunnumbersPhi->Sumw2();
444 fOutputList->Add(fHistRunnumbersPhi);
445 for(Int_t i(0); i < fExpectedRuns->GetSize(); i++) {
446 fHistRunnumbersPhi->GetXaxis()->SetBinLabel(i+1, Form("%i", fExpectedRuns->At(i)));
447 fHistRunnumbersEta->GetXaxis()->SetBinLabel(i+1, Form("%i", fExpectedRuns->At(i)));
449 fHistRunnumbersPhi->GetXaxis()->SetBinLabel(fExpectedRuns->GetSize()+1, "undetermined");
450 fHistRunnumbersEta->GetXaxis()->SetBinLabel(fExpectedRuns->GetSize()+1, "undetermined");
452 fHistAnalysisSummary = BookTH1F("fHistAnalysisSummary", "flag", 54, -0.5, 54.5);
453 fHistSwap = new TH1F("fHistSwap", "fHistSwap", 20, 0, TMath::TwoPi());
454 if(fUsePtWeight) fHistSwap->Sumw2();
456 if(fUserSuppliedV2) fOutputList->Add(fUserSuppliedV2);
457 if(fUserSuppliedV3) fOutputList->Add(fUserSuppliedV3);
458 if(fUserSuppliedR2) fOutputList->Add(fUserSuppliedR2);
459 if(fUserSuppliedR3) fOutputList->Add(fUserSuppliedR3);
460 // increase readability of output list
462 // cdf and pdf of chisquare distribution
463 fHistPvalueCDF = BookTH1F("fHistPvalueCDF", "CDF #chi^{2}", 50, 0, 1);
464 fHistPvalueCDFCent = BookTH2F("fHistPvalueCDFCent", "centrality", "p-value", 40, 0, 100, 40, 0, 1);
465 fHistChi2Cent = BookTH2F("fHistChi2Cent", "centrality", "#tilde{#chi^{2}}", 100, 0, 100, 100, 0, 5);
466 fHistPChi2 = BookTH2F("fHistPChi2", "p-value", "#tilde{#chi^{2}}", 1000, 0, 1, 100, 0, 5);
467 fHistKolmogorovTest = BookTH1F("fHistKolmogorovTest", "KolmogorovTest", 50, 0, 1);
468 fHistKolmogorovTestCent = BookTH2F("fHistKolmogorovTestCent", "centrality", "Kolmogorov p", 40, 0, 100, 45, 0, 1);
469 fHistPvalueCDFROOT = BookTH1F("fHistPvalueCDFROOT", "CDF #chi^{2} ROOT", 50, 0, 1);
470 fHistPvalueCDFROOTCent = BookTH2F("fHistPvalueCDFROOTCent", "centrality", "p-value ROOT", 40, 0, 100, 45, 0, 1);
471 fHistChi2ROOTCent = BookTH2F("fHistChi2ROOTCent", "centrality", "#tilde{#chi^{2}}", 40, 0, 100, 45, 0, 5);
472 fHistPChi2Root = BookTH2F("fHistPChi2Root", "p-value", "#tilde{#chi^{2}} ROOT", 1000, 0, 1, 100, 0, 5);
473 fHistPKolmogorov = BookTH2F("fHistPKolmogorov", "p-value", "kolmogorov p",40, 0, 1, 40, 0, 1);
474 fHistRhoStatusCent = BookTH2F("fHistRhoStatusCent", "centrality", "status [-1=lin was better, 0=ok, 1 = failed]", 101, -1, 100, 3, -1.5, 1.5);
475 fHistUndeterminedRunQA = BookTH1F("fHistUndeterminedRunQA", "runnumber", 10, 0, 10);
477 PostData(1, fOutputList);
479 switch (fRunModeType) {
481 fOutputListGood = new TList();
482 fOutputListGood->SetOwner(kTRUE);
483 fOutputListBad = new TList();
484 fOutputListBad->SetOwner(kTRUE);
485 PostData(2, fOutputListGood);
486 PostData(3, fOutputListBad);
491 // get the containers
492 fTracksCont = GetParticleContainer("Tracks");
493 fClusterCont = GetClusterContainer(0); // get the default cluster container
494 fJetsCont = GetJetContainer("Jets");
496 //_____________________________________________________________________________
497 Bool_t AliAnalysisTaskJetV2::Run()
499 // user exec: execute once for each event
500 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
501 if(!fTracks||!fJets||!fRho) return kFALSE;
502 if(!fLocalInit) fLocalInit = InitializeAnalysis();
503 // reject the event if expected data is missing
504 if(!PassesCuts(InputEvent())) return kFALSE;
505 fLeadingJet = GetLeadingJet(); // store the leading jet
507 fLocalRho->SetVal(fRho->GetVal());
508 // [0][0] psi2a [1,0] psi2c
509 // [0][1] psi3a [1,1] psi3c
510 Double_t vzero[2][2];
511 CalculateEventPlaneVZERO(vzero);
512 /* for the combined vzero event plane
514 * not fully implmemented yet, use with caution ! */
515 Double_t vzeroComb[2];
516 CalculateEventPlaneCombinedVZERO(vzeroComb);
519 CalculateEventPlaneTPC(tpc);
520 Double_t psi2(-1), psi3(-1);
521 // arrays which will hold the fit parameters
522 switch (fDetectorType) { // determine the detector type for the rho fit
523 case kTPC : { psi2 = tpc[0]; psi3 = tpc[1]; } break;
524 case kVZEROA : { psi2 = vzero[0][0]; psi3 = vzero[0][1]; } break;
525 case kVZEROC : { psi2 = vzero[1][0]; psi3 = vzero[1][1]; } break;
526 case kVZEROComb : { psi2 = vzeroComb[0]; psi3 = vzeroComb[1];} break;
529 switch (fFitModulationType) { // do the fits
531 switch (fCollisionType) {
532 case kPythia : { // background is zero for pp jets
533 fFitModulation->FixParameter(0, 0);
534 fLocalRho->SetVal(0);
537 fFitModulation->FixParameter(0, fLocalRho->GetVal());
541 case kV2 : { // only v2
542 if(CorrectRho(psi2, psi3)) {
543 fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
544 if(fUserSuppliedR2) {
545 Double_t r(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
546 if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
548 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
551 case kV3 : { // only v3
552 if(CorrectRho(psi2, psi3)) {
553 if(fUserSuppliedR3) {
554 Double_t r(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
555 if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
557 fProfV3->Fill(fCent, fFitModulation->GetParameter(3));
558 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
561 case kQC2 : { // qc2 analysis
562 if(CorrectRho(psi2, psi3)) {
563 if(fUserSuppliedR2 && fUserSuppliedR3) {
564 // note for the qc method, resolution is REVERSED to go back to v2obs
565 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
566 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
567 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
568 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)*r3);
570 if (fUsePtWeight) { // use weighted weights
571 Double_t dQCnM11 = (fNoEventWeightsForQC) ? 1. : QCnM11();
572 fProfV2->Fill(fCent, fFitModulation->GetParameter(3), dQCnM11);
573 fProfV3->Fill(fCent, fFitModulation->GetParameter(7), dQCnM11);
575 Double_t dQCnM = (fNoEventWeightsForQC) ? 2. : QCnM();
576 fProfV2->Fill(fCent, fFitModulation->GetParameter(3), dQCnM*(dQCnM-1));
577 fProfV3->Fill(fCent, fFitModulation->GetParameter(7), dQCnM*(dQCnM-1));
579 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
583 if(CorrectRho(psi2, psi3)) {
584 if(fUserSuppliedR2 && fUserSuppliedR3) {
585 // note for the qc method, resolution is REVERSED to go back to v2obs
586 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
587 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
588 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
589 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)*r3);
591 if (fUsePtWeight) { // use weighted weights
592 fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM1111()*/);
593 fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM1111()*/);
595 fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
596 fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
599 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
602 if(CorrectRho(psi2, psi3)) {
603 if(fUserSuppliedR2 && fUserSuppliedR3) {
604 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
605 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
606 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r2);
607 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)/r3);
609 fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
610 fProfV3->Fill(fCent, fFitModulation->GetParameter(7));
611 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
615 // if all went well, update the local rho parameter
616 fLocalRho->SetLocalRho(fFitModulation);
617 // fill a number of histograms. event qa needs to be filled first as it also determines the runnumber for the track qa
618 if(fFillQAHistograms) FillQAHistograms(InputEvent());
619 if(fFillHistograms) FillHistogramsAfterSubtraction(psi2, vzero, vzeroComb, tpc);
620 // send the output to the connected output container
621 PostData(1, fOutputList);
622 switch (fRunModeType) {
624 PostData(2, fOutputListGood);
625 PostData(3, fOutputListBad);
632 //_____________________________________________________________________________
633 void AliAnalysisTaskJetV2::CalculateEventPlaneVZERO(Double_t vzero[2][2]) const
635 // get the vzero event plane
636 if(fUseV0EventPlaneFromHeader) { // use the vzero from the header
637 Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0), h(0);
638 vzero[0][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 2, a, b);
639 vzero[1][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 2, c, d);
640 vzero[0][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 3, e, f);
641 vzero[1][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 3, g, h);
644 // grab the vzero event plane without recentering
645 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
646 Double_t qxa2(0), qya2(0), qxc2(0), qyc2(0); // for psi2
647 Double_t qxa3(0), qya3(0), qxc3(0), qyc3(0); // for psi3
648 for(Int_t iVZERO(0); iVZERO < 64; iVZERO++) {
649 Double_t phi(TMath::PiOver4()*(.5+iVZERO%8)), /* eta(0), */ weight(InputEvent()->GetVZEROEqMultiplicity(iVZERO));
651 qxa2 += weight*TMath::Cos(2.*phi);
652 qya2 += weight*TMath::Sin(2.*phi);
653 qxa3 += weight*TMath::Cos(3.*phi);
654 qya3 += weight*TMath::Sin(3.*phi);
657 qxc2 += weight*TMath::Cos(2.*phi);
658 qyc2 += weight*TMath::Sin(2.*phi);
659 qxc3 += weight*TMath::Cos(3.*phi);
660 qyc3 += weight*TMath::Sin(3.*phi);
663 vzero[0][0] = .5*TMath::ATan2(qya2, qxa2);
664 vzero[1][0] = .5*TMath::ATan2(qyc2, qxc2);
665 vzero[0][1] = (1./3.)*TMath::ATan2(qya3, qxa3);
666 vzero[1][1] = (1./3.)*TMath::ATan2(qyc3, qxc3);
668 //_____________________________________________________________________________
669 void AliAnalysisTaskJetV2::CalculateEventPlaneTPC(Double_t* tpc)
671 // grab the TPC event plane
672 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
673 fNAcceptedTracks = 0; // reset the track counter
674 Double_t qx2(0), qy2(0); // for psi2
675 Double_t qx3(0), qy3(0); // for psi3
677 Float_t excludeInEta = -999;
678 if(fExcludeLeadingJetsFromFit > 0 ) { // remove the leading jet from ep estimate
679 if(fLeadingJet) excludeInEta = fLeadingJet->Eta();
681 for(Int_t iTPC(0); iTPC < fTracksCont->GetNEntries(); iTPC++) {
682 AliVParticle* track = fTracksCont->GetParticle(iTPC);
683 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
684 if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta) < GetJetContainer()->GetJetRadius()*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - GetJetContainer()->GetJetRadius() - GetJetContainer()->GetJetEtaMax() ) > 0 )) continue;
686 qx2+= TMath::Cos(2.*track->Phi());
687 qy2+= TMath::Sin(2.*track->Phi());
688 qx3+= TMath::Cos(3.*track->Phi());
689 qy3+= TMath::Sin(3.*track->Phi());
692 tpc[0] = .5*TMath::ATan2(qy2, qx2);
693 tpc[1] = (1./3.)*TMath::ATan2(qy3, qx3);
695 //_____________________________________________________________________________
696 void AliAnalysisTaskJetV2::CalculateEventPlaneCombinedVZERO(Double_t* comb) const
698 // grab the combined vzero event plane
699 Double_t a(0), b(0), c(0), d(0);
700 comb[0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 10, 2, a, b);
701 comb[1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 10, 3, c, d);
703 //_____________________________________________________________________________
704 void AliAnalysisTaskJetV2::CalculateEventPlaneResolution(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc)
706 // fill the profiles for the resolution parameters
707 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
708 fProfV2Resolution[fInCentralitySelection]->Fill(2., TMath::Cos(2.*(vzero[0][0] - vzero[1][0])));
709 fProfV2Resolution[fInCentralitySelection]->Fill(3., TMath::Cos(2.*(vzero[1][0] - vzero[0][0])));
710 fProfV2Resolution[fInCentralitySelection]->Fill(4., TMath::Cos(2.*(vzero[0][0] - tpc[0])));
711 fProfV2Resolution[fInCentralitySelection]->Fill(5., TMath::Cos(2.*(tpc[0] - vzero[0][0])));
712 fProfV2Resolution[fInCentralitySelection]->Fill(6., TMath::Cos(2.*(vzero[1][0] - tpc[0])));
713 fProfV2Resolution[fInCentralitySelection]->Fill(7., TMath::Cos(2.*(tpc[0] - vzero[1][0])));
714 fProfV3Resolution[fInCentralitySelection]->Fill(2., TMath::Cos(3.*(vzero[0][0] - vzero[1][0])));
715 fProfV3Resolution[fInCentralitySelection]->Fill(3., TMath::Cos(3.*(vzero[1][0] - vzero[0][0])));
716 fProfV3Resolution[fInCentralitySelection]->Fill(4., TMath::Cos(3.*(vzero[0][0] - tpc[0])));
717 fProfV3Resolution[fInCentralitySelection]->Fill(5., TMath::Cos(3.*(tpc[0] - vzero[0][0])));
718 fProfV3Resolution[fInCentralitySelection]->Fill(6., TMath::Cos(3.*(vzero[1][0] - tpc[0])));
719 fProfV3Resolution[fInCentralitySelection]->Fill(7., TMath::Cos(3.*(tpc[0] - vzero[1][0])));
720 // for the resolution of the combined vzero event plane, use two tpc halves as uncorrelated subdetectors
721 Double_t qx2a(0), qy2a(0); // for psi2a, negative eta
722 Double_t qx3a(0), qy3a(0); // for psi3a, negative eta
723 Double_t qx2b(0), qy2b(0); // for psi2a, positive eta
724 Double_t qx3b(0), qy3b(0); // for psi3a, positive eta
726 Int_t iTracks(fTracks->GetEntriesFast());
727 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
728 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
729 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
730 if(track->Eta() < 0 ) {
731 qx2a+= TMath::Cos(2.*track->Phi());
732 qy2a+= TMath::Sin(2.*track->Phi());
733 qx3a+= TMath::Cos(3.*track->Phi());
734 qy3a+= TMath::Sin(3.*track->Phi());
735 } else if (track->Eta() > 0) {
736 qx2b+= TMath::Cos(2.*track->Phi());
737 qy2b+= TMath::Sin(2.*track->Phi());
738 qx3b+= TMath::Cos(3.*track->Phi());
739 qy3b+= TMath::Sin(3.*track->Phi());
743 Double_t tpca2(.5*TMath::ATan2(qy2a, qx2a));
744 Double_t tpca3((1./3.)*TMath::ATan2(qy3a, qx3a));
745 Double_t tpcb2(.5*TMath::ATan2(qy2b, qx2b));
746 Double_t tpcb3((1./3.)*TMath::ATan2(qy3b, qx3b));
747 fProfV2Resolution[fInCentralitySelection]->Fill(8., TMath::Cos(2.*(vzeroComb[0] - tpca2)));
748 fProfV2Resolution[fInCentralitySelection]->Fill(9., TMath::Cos(2.*(vzeroComb[0] - tpcb2)));
749 fProfV2Resolution[fInCentralitySelection]->Fill(10., TMath::Cos(2.*(tpca2 - tpcb2)));
750 fProfV3Resolution[fInCentralitySelection]->Fill(8., TMath::Cos(3.*(vzeroComb[1] - tpca3)));
751 fProfV3Resolution[fInCentralitySelection]->Fill(9., TMath::Cos(3.*(vzeroComb[1] - tpcb3)));
752 fProfV3Resolution[fInCentralitySelection]->Fill(10., TMath::Cos(3.*(tpca3 - tpcb3)));
754 //_____________________________________________________________________________
755 void AliAnalysisTaskJetV2::CalculateRandomCone(Float_t &pt, Float_t &eta, Float_t &phi,
756 AliParticleContainer* tracksCont, AliClusterContainer* clusterCont, AliEmcalJet* jet) const
759 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
760 pt = 0; eta = 0; phi = 0;
761 Float_t etaJet(999), phiJet(999), dJet(999); // no jet: same as jet very far away
762 if(jet) { // if a leading jet is given, use its kinematic properties to exclude it
766 // the random cone acceptance has to equal the jet acceptance
767 // this also insures safety when runnnig on the semi-good tpc runs for 11h data,
768 // where jet acceptance is adjusted to reduced acceptance - hence random cone acceptance as well
769 Float_t minPhi(GetJetContainer()->GetJetPhiMin()), maxPhi(GetJetContainer()->GetJetPhiMax());
770 if(maxPhi > TMath::TwoPi()) maxPhi = TMath::TwoPi();
771 if(minPhi < 0 ) minPhi = 0.;
772 // construct a random cone and see if it's far away enough from the leading jet
773 Int_t attempts(1000);
776 eta = gRandom->Uniform(GetJetContainer()->GetJetEtaMin(), GetJetContainer()->GetJetEtaMax());
777 phi = gRandom->Uniform(minPhi, maxPhi);
779 dJet = TMath::Sqrt((etaJet-eta)*(etaJet-eta)+(phiJet-phi)*(phiJet-phi));
780 if(dJet > fMinDisanceRCtoLJ) break;
781 else if (attempts == 0) {
782 printf(" > No random cone after 1000 tries, giving up ... !\n");
786 // get the charged energy (if tracks are provided)
788 AliVParticle* track = tracksCont->GetNextAcceptParticle(0);
790 Float_t etaTrack(track->Eta()), phiTrack(track->Phi());
791 // get distance from cone
792 if(TMath::Abs(phiTrack-phi) > TMath::Abs(phiTrack - phi + TMath::TwoPi())) phiTrack+=TMath::TwoPi();
793 if(TMath::Abs(phiTrack-phi) > TMath::Abs(phiTrack - phi - TMath::TwoPi())) phiTrack-=TMath::TwoPi();
794 if(TMath::Sqrt(TMath::Abs((etaTrack-eta)*(etaTrack-eta)+(phiTrack-phi)*(phiTrack-phi))) <= GetJetRadius()) pt += track->Pt();
795 track = tracksCont->GetNextAcceptParticle();
798 // get the neutral energy (if clusters are provided)
800 AliVCluster* cluster = clusterCont->GetNextAcceptCluster(0);
802 TLorentzVector momentum;
803 cluster->GetMomentum(momentum, const_cast<Double_t*>(fVertex));
804 Float_t etaClus(momentum.Eta()), phiClus(momentum.Phi());
805 // get distance from cone
806 if(TMath::Abs(phiClus-phi) > TMath::Abs(phiClus - phi + TMath::TwoPi())) phiClus+=TMath::TwoPi();
807 if(TMath::Abs(phiClus-phi) > TMath::Abs(phiClus - phi - TMath::TwoPi())) phiClus-=TMath::TwoPi();
808 if(TMath::Sqrt(TMath::Abs((etaClus-eta)*(etaClus-eta)+(phiClus-phi)*(phiClus-phi))) <= GetJetRadius()) pt += momentum.Pt();
809 cluster = clusterCont->GetNextAcceptCluster();
813 //_____________________________________________________________________________
814 Double_t AliAnalysisTaskJetV2::CalculateQC2(Int_t harm) {
815 // get the second order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho
816 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
817 Double_t reQ(0), imQ(0), modQ(0), M11(0), M(0);
818 if(fUsePtWeight) { // for the weighted 2-nd order q-cumulant
819 QCnQnk(harm, 1, reQ, imQ); // get the weighted 2-nd order q-vectors
820 modQ = reQ*reQ+imQ*imQ; // get abs Q-squared
821 M11 = QCnM11(); // equals S2,1 - S1,2
822 return (M11 > 0) ? ((modQ - QCnS(1,2))/M11) : -999;
823 } // else return the non-weighted 2-nd order q-cumulant
824 QCnQnk(harm, 0, reQ, imQ); // get the non-weighted 2-nd order q-vectors
825 modQ = reQ*reQ+imQ*imQ; // get abs Q-squared
827 return (M > 1) ? (modQ - M)/(M*(M-1)) : -999;
829 //_____________________________________________________________________________
830 Double_t AliAnalysisTaskJetV2::CalculateQC4(Int_t harm) {
831 // get the fourth order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho
832 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
833 Double_t reQn1(0), imQn1(0), reQ2n2(0), imQ2n2(0), reQn3(0), imQn3(0), M1111(0), M(0);
834 Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0); // terms of the calculation
835 if(fUsePtWeight) { // for the weighted 4-th order q-cumulant
836 QCnQnk(harm, 1, reQn1, imQn1);
837 QCnQnk(harm*2, 2, reQ2n2, imQ2n2);
838 QCnQnk(harm, 3, reQn3, imQn3);
839 // fill in the terms ...
840 a = (reQn1*reQn1+imQn1*imQn1)*(reQn1*reQn1+imQn1*imQn1);
841 b = reQ2n2*reQ2n2 + imQ2n2*imQ2n2;
842 c = -2.*(reQ2n2*reQn1*reQn1-reQ2n2*imQn1*imQn1+2.*imQ2n2*reQn1*imQn1);
843 d = 8.*(reQn3*reQn1+imQn3*imQn1);
844 e = -4.*QCnS(1,2)*(reQn1*reQn1+imQn1*imQn1);
848 return (M1111 > 0) ? (a+b+c+d+e+f+g)/M1111 : -999;
849 } // else return the unweighted case
850 Double_t reQn(0), imQn(0), reQ2n(0), imQ2n(0);
851 QCnQnk(harm, 0, reQn, imQn);
852 QCnQnk(harm*2, 0, reQ2n, imQ2n);
853 // fill in the terms ...
855 if(M < 4) return -999;
856 a = (reQn*reQn+imQn*imQn)*(reQn*reQn+imQn*imQn);
857 b = reQ2n*reQ2n + imQ2n*imQ2n;
858 c = -2.*(reQ2n*reQn*reQn-reQ2n*imQn*imQn+2.*imQ2n*reQn*imQn);
859 e = -4.*(M-2)*(reQn*reQn+imQn*imQn);
861 return (a+b+c+e+f)/(M*(M-1)*(M-2)*(M-3));
863 //_____________________________________________________________________________
864 void AliAnalysisTaskJetV2::QCnQnk(Int_t n, Int_t k, Double_t &reQ, Double_t &imQ) {
865 // get the weighted n-th order q-vector, pass real and imaginary part as reference
866 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
868 fNAcceptedTracksQCn = 0;
869 Int_t iTracks(fTracks->GetEntriesFast());
870 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
871 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
872 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
873 fNAcceptedTracksQCn++;
874 // for the unweighted case, k equals zero and the weight doesn't contribute to the equation below
875 reQ += TMath::Power(track->Pt(), k) * TMath::Cos(((double)n)*track->Phi());
876 imQ += TMath::Power(track->Pt(), k) * TMath::Sin(((double)n)*track->Phi());
879 //_____________________________________________________________________________
880 void AliAnalysisTaskJetV2::QCnDiffentialFlowVectors(
881 TClonesArray* pois, TArrayD* ptBins, Bool_t vpart, Double_t* repn, Double_t* impn,
882 Double_t *mp, Double_t *reqn, Double_t *imqn, Double_t* mq, Int_t n)
884 // get unweighted differential flow vectors
885 Int_t iPois(pois->GetEntriesFast());
887 for(Int_t i(0); i < iPois; i++) {
888 for(Int_t ptBin(0); ptBin < ptBins->GetSize()-1; ptBin++) {
889 AliVTrack* poi = static_cast<AliVTrack*>(pois->At(i));
890 if(PassesCuts(poi)) {
891 if(poi->Pt() >= ptBins->At(ptBin) && poi->Pt() < ptBins->At(ptBin+1)) {
892 // fill the flow vectors assuming that all poi's are in the rp selection (true by design)
893 repn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
894 impn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
896 reqn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
897 imqn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
904 for(Int_t i(0); i < iPois; i++) {
905 for(Int_t ptBin(0); ptBin < ptBins->GetSize()-1; ptBin++) {
906 AliEmcalJet* poi = static_cast<AliEmcalJet*>(pois->At(i));
907 if(PassesCuts(poi)) {
908 Double_t pt(poi->Pt()-poi->Area()*fLocalRho->GetLocalVal(poi->Phi(), GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
909 if(pt >= ptBins->At(ptBin) && pt < ptBins->At(ptBin+1)) {
910 repn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
911 impn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
912 mp[ptBin]++; // qn isn't filled, no overlap between poi's and rp's
919 //_____________________________________________________________________________
920 Double_t AliAnalysisTaskJetV2::QCnS(Int_t i, Int_t j) {
921 // get the weighted ij-th order autocorrelation correction
922 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
923 if(!fTracks || i <= 0 || j <= 0) return -999;
924 Int_t iTracks(fTracks->GetEntriesFast());
926 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
927 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
928 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
929 Sij+=TMath::Power(track->Pt(), j);
931 return TMath::Power(Sij, i);
933 //_____________________________________________________________________________
934 Double_t AliAnalysisTaskJetV2::QCnM() {
935 // get multiplicity for unweighted q-cumulants. function QCnQnk should be called first
936 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
937 return (Double_t) fNAcceptedTracksQCn;
939 //_____________________________________________________________________________
940 Double_t AliAnalysisTaskJetV2::QCnM11() {
941 // get multiplicity weights for the weighted two particle cumulant
942 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
943 return (QCnS(2,1) - QCnS(1,2));
945 //_____________________________________________________________________________
946 Double_t AliAnalysisTaskJetV2::QCnM1111() {
947 // get multiplicity weights for the weighted four particle cumulant
948 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
949 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));
951 //_____________________________________________________________________________
952 Bool_t AliAnalysisTaskJetV2::QCnRecovery(Double_t psi2, Double_t psi3) {
953 // decides how to deal with the situation where c2 or c3 is negative
954 // returns kTRUE depending on whether or not a modulated rho is used for the jet background
955 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
956 if(TMath::AreEqualAbs(fFitModulation->GetParameter(3), .0, 1e-10) && TMath::AreEqualAbs(fFitModulation->GetParameter(7), .0,1e-10)) {
957 fFitModulation->SetParameter(7, 0);
958 fFitModulation->SetParameter(3, 0);
959 fFitModulation->SetParameter(0, fLocalRho->GetVal());
960 return kTRUE; // v2 and v3 have physical null values
962 switch (fQCRecovery) {
963 case kFixedRho : { // roll back to the original rho
964 fFitModulation->SetParameter(7, 0);
965 fFitModulation->SetParameter(3, 0);
966 fFitModulation->SetParameter(0, fLocalRho->GetVal());
967 return kFALSE; // rho is forced to be fixed
970 Double_t c2(fFitModulation->GetParameter(3));
971 Double_t c3(fFitModulation->GetParameter(7));
972 if( c2 < 0 ) c2 = -1.*TMath::Sqrt(-1.*c2);
973 if( c3 < 0 ) c3 = -1.*TMath::Sqrt(-1.*c3);
974 fFitModulation->SetParameter(3, c2);
975 fFitModulation->SetParameter(7, c3);
976 return kTRUE; // is this a physical quantity ?
979 fitModulationType tempType(fFitModulationType); // store temporarily
980 fFitModulationType = kCombined;
981 fFitModulation->SetParameter(7, 0);
982 fFitModulation->SetParameter(3, 0);
983 Bool_t pass(CorrectRho(psi2, psi3)); // do the fit and all quality checks
984 fFitModulationType = tempType; // roll back for next event
987 default : return kFALSE;
991 //_____________________________________________________________________________
992 Bool_t AliAnalysisTaskJetV2::CorrectRho(Double_t psi2, Double_t psi3)
994 // get rho' -> rho(phi)
995 // two routines are available, both can be used with or without pt weights
996 // [1] get vn from q-cumulants or as an integrated value from a user supplied histogram
997 // in case of cumulants, both cumulants and vn values are stored. in both cases, v2 and v3
998 // are expected. a check is performed to see if rho has no negative local minimum
999 // for full description, see Phys. Rev. C 83, 044913
1000 // since the cn distribution has negative values, vn = sqrt(cn) can be imaginary sometimes
1001 // in this case one can either roll back to the 'original' rixed rho, do a fit for vn or take use
1002 // vn = - sqrt(|cn|)
1003 // [2] fitting a fourier expansion to the de/dphi distribution
1004 // the fit can be done with either v2, v3 or a combination.
1005 // in all cases, a cut can be made on the p-value of the chi-squared value of the fit
1006 // and a check can be performed to see if rho has no negative local minimum
1007 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1008 Int_t freeParams(2); // free parameters of the fit (for NDF)
1009 switch (fFitModulationType) { // for approaches where no fitting is required
1011 fFitModulation->FixParameter(4, psi2);
1012 fFitModulation->FixParameter(6, psi3);
1013 fFitModulation->FixParameter(3, CalculateQC2(2)); // set here with cn, vn = sqrt(cn)
1014 fFitModulation->FixParameter(7, CalculateQC2(3));
1015 // first fill the histos of the raw cumulant distribution
1016 if (fUsePtWeight) { // use weighted weights
1017 Double_t dQCnM11 = (fNoEventWeightsForQC) ? 1. : QCnM11();
1018 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), dQCnM11);
1019 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), dQCnM11);
1021 Double_t dQCnM = (fNoEventWeightsForQC) ? 2. : QCnM();
1022 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), dQCnM*(dQCnM-1));
1023 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), dQCnM*(dQCnM-1));
1025 // then see if one of the cn value is larger than zero and vn is readily available
1026 if(fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) {
1027 fFitModulation->FixParameter(3, TMath::Sqrt(fFitModulation->GetParameter(3)));
1028 fFitModulation->FixParameter(7, TMath::Sqrt(fFitModulation->GetParameter(7)));
1029 } else if (!QCnRecovery(psi2, psi3)) return kFALSE; // try to recover the cumulant, this will set v2 and v3
1030 if(fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { // general check
1031 fFitModulation->SetParameter(7, 0);
1032 fFitModulation->SetParameter(3, 0);
1033 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1039 fFitModulation->FixParameter(4, psi2);
1040 fFitModulation->FixParameter(6, psi3);
1041 fFitModulation->FixParameter(3, CalculateQC4(2)); // set here with cn, vn = sqrt(cn)
1042 fFitModulation->FixParameter(7, CalculateQC4(3));
1043 // first fill the histos of the raw cumulant distribution
1044 if (fUsePtWeight) { // use weighted weights
1045 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
1046 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
1048 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
1049 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
1051 // then see if one of the cn value is larger than zero and vn is readily available
1052 if(fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) {
1053 fFitModulation->FixParameter(3, TMath::Sqrt(fFitModulation->GetParameter(3)));
1054 fFitModulation->FixParameter(7, TMath::Sqrt(fFitModulation->GetParameter(7)));
1055 } else if (!QCnRecovery(psi2, psi3)) return kFALSE; // try to recover the cumulant, this will set v2 and v3
1056 if(fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { // general check
1057 fFitModulation->SetParameter(7, 0);
1058 fFitModulation->SetParameter(3, 0);
1059 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1063 case kIntegratedFlow : {
1064 // use v2 and v3 values from an earlier iteration over the data
1065 fFitModulation->FixParameter(3, fUserSuppliedV2->GetBinContent(fUserSuppliedV2->GetXaxis()->FindBin(fCent)));
1066 fFitModulation->FixParameter(4, psi2);
1067 fFitModulation->FixParameter(6, psi3);
1068 fFitModulation->FixParameter(7, fUserSuppliedV3->GetBinContent(fUserSuppliedV3->GetXaxis()->FindBin(fCent)));
1069 if(fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) {
1070 fFitModulation->SetParameter(7, 0);
1071 fFitModulation->SetParameter(3, 0);
1072 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1079 TString detector("");
1080 switch (fDetectorType) {
1081 case kTPC : detector+="TPC";
1083 case kVZEROA : detector+="VZEROA";
1085 case kVZEROC : detector+="VZEROC";
1087 case kVZEROComb : detector+="VZEROComb";
1091 Int_t iTracks(fTracks->GetEntriesFast());
1092 Double_t excludeInEta = -999;
1093 Double_t excludeInPhi = -999;
1094 Double_t excludeInPt = -999;
1095 if(iTracks <= 0 || fLocalRho->GetVal() <= 0 ) return kFALSE; // no use fitting an empty event ...
1096 if(fExcludeLeadingJetsFromFit > 0 ) {
1098 excludeInEta = fLeadingJet->Eta();
1099 excludeInPhi = fLeadingJet->Phi();
1100 excludeInPt = fLeadingJet->Pt();
1103 // check the acceptance of the track selection that will be used
1104 // if one uses e.g. semi-good tpc tracks, accepance in phi is reduced to 0 < phi < 4
1105 // the defaults (-10 < phi < 10) which accept all, are then overwritten
1106 Double_t lowBound(0.), upBound(TMath::TwoPi()); // bounds for fit
1107 if(GetParticleContainer()->GetParticlePhiMin() > lowBound) lowBound = GetParticleContainer()->GetParticlePhiMin();
1108 if(GetParticleContainer()->GetParticlePhiMax() < upBound) upBound = GetParticleContainer()->GetParticlePhiMax();
1110 fHistSwap->Reset(); // clear the histogram
1111 TH1F _tempSwap; // on stack for quick access
1112 TH1F _tempSwapN; // on stack for quick access, bookkeeping histogram
1113 if(fRebinSwapHistoOnTheFly) {
1114 if(fNAcceptedTracks < 49) fNAcceptedTracks = 49; // avoid aliasing effects
1115 _tempSwap = TH1F("_tempSwap", "_tempSwap", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), lowBound, upBound);
1116 if(fUsePtWeightErrorPropagation) _tempSwapN = TH1F("_tempSwapN", "_tempSwapN", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), lowBound, upBound);
1117 if(fUsePtWeight) _tempSwap.Sumw2();
1119 else _tempSwap = *fHistSwap; // now _tempSwap holds the desired histo
1120 // non poissonian error when using pt weights
1121 Double_t totalpts(0.), totalptsquares(0.), totalns(0.);
1122 for(Int_t i(0); i < iTracks; i++) {
1123 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1124 if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta) < GetJetContainer()->GetJetRadius()*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - GetJetContainer()->GetJetRadius() - GetJetContainer()->GetJetEtaMax() ) > 0 )) continue;
1125 if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
1127 _tempSwap.Fill(track->Phi(), track->Pt());
1128 if(fUsePtWeightErrorPropagation) {
1129 totalpts += track->Pt();
1130 totalptsquares += track->Pt()*track->Pt();
1132 _tempSwapN.Fill(track->Phi());
1135 else _tempSwap.Fill(track->Phi());
1137 if(fUsePtWeight && fUsePtWeightErrorPropagation) {
1138 // in the case of pt weights overwrite the poissonian error estimate which is assigned by root by a more sophisticated appraoch
1139 // 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
1140 // 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
1141 // event, for the latter use a poissonian estimate. the two contrubitions are assumed to be uncorrelated
1142 if(totalns < 1) return kFALSE; // not one track passes the cuts
1143 for(Int_t l = 0; l < _tempSwap.GetNbinsX(); l++) {
1144 if(_tempSwapN.GetBinContent(l+1) == 0) {
1145 _tempSwap.SetBinContent(l+1,0);
1146 _tempSwap.SetBinError(l+1,0);
1149 Double_t vartimesnsq = totalptsquares*totalns - totalpts*totalpts;
1150 Double_t variance = vartimesnsq/(totalns*(totalns-1.));
1151 Double_t SDOMSq = variance / _tempSwapN.GetBinContent(l+1);
1152 Double_t SDOMSqOverMeanSq = SDOMSq * _tempSwapN.GetBinContent(l+1) * _tempSwapN.GetBinContent(l+1) / (_tempSwapN.GetBinContent(l+1) * _tempSwapN.GetBinContent(l+1));
1153 Double_t poissonfrac = 1./_tempSwapN.GetBinContent(l+1);
1154 Double_t vartotalfrac = SDOMSqOverMeanSq + poissonfrac;
1155 Double_t vartotal = vartotalfrac * _tempSwap.GetBinContent(l+1) * _tempSwap.GetBinContent(l+1);
1156 if(vartotal > 0.0001) _tempSwap.SetBinError(l+1,TMath::Sqrt(vartotal));
1158 _tempSwap.SetBinContent(l+1,0);
1159 _tempSwap.SetBinError(l+1,0);
1165 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1166 switch (fFitModulationType) {
1168 fFitModulation->FixParameter(0, fLocalRho->GetVal() );
1172 fFitModulation->FixParameter(4, psi2);
1176 fFitModulation->FixParameter(4, psi3);
1180 fFitModulation->FixParameter(4, psi2);
1181 fFitModulation->FixParameter(6, psi3);
1184 case kFourierSeries : {
1185 // in this approach, an explicit calculation will be made of vn = sqrt(xn^2+yn^2)
1186 // where x[y] = Integrate[r(phi)cos[sin](n phi)dphi, 0, 2pi]
1187 Double_t cos2(0), sin2(0), cos3(0), sin3(0), sumPt(0);
1188 for(Int_t i(0); i < iTracks; i++) {
1189 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1190 if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
1191 sumPt += track->Pt();
1192 cos2 += track->Pt()*TMath::Cos(2*PhaseShift(track->Phi()-psi2));
1193 sin2 += track->Pt()*TMath::Sin(2*PhaseShift(track->Phi()-psi2));
1194 cos3 += track->Pt()*TMath::Cos(3*PhaseShift(track->Phi()-psi3));
1195 sin3 += track->Pt()*TMath::Sin(3*PhaseShift(track->Phi()-psi3));
1197 fFitModulation->SetParameter(3, TMath::Sqrt(cos2*cos2+sin2*sin2)/fLocalRho->GetVal());
1198 fFitModulation->SetParameter(4, psi2);
1199 fFitModulation->SetParameter(6, psi3);
1200 fFitModulation->SetParameter(7, TMath::Sqrt(cos3*cos3+sin3*sin3)/fLocalRho->GetVal());
1205 // toy mc, just here to check procedure, azimuthal profile is filled from hypothesis so p-value distribution should be flat
1206 Int_t _bins = _tempSwap.GetXaxis()->GetNbins();
1207 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());
1208 _tempFit->SetParameter(0, fFitModulation->GetParameter(0)); // normalization
1209 _tempFit->SetParameter(3, 0.1); // v2
1210 _tempFit->FixParameter(1, 1.); // constant
1211 _tempFit->FixParameter(2, 2.); // constant
1212 _tempFit->FixParameter(5, 3.); // constant
1213 _tempFit->FixParameter(4, fFitModulation->GetParameter(4));
1214 _tempFit->FixParameter(6, fFitModulation->GetParameter(6));
1215 _tempFit->SetParameter(7, 0.1); // v3
1216 _tempSwap.Reset(); // rese bin content
1217 for(int _binsI = 0; _binsI < _bins*_bins; _binsI++) _tempSwap.Fill(_tempFit->GetRandom());
1219 _tempSwap.Fit(fFitModulation, fFitModulationOptions.Data(), "", lowBound, upBound);
1220 // the quality of the fit is evaluated from 1 - the cdf of the chi square distribution
1221 // three methods are available, all with their drawbacks. all are stored, one is selected to do the cut
1222 Int_t NDF(_tempSwap.GetXaxis()->GetNbins()-freeParams);
1223 if(NDF == 0) return kFALSE;
1224 Double_t CDF(1.-ChiSquareCDF(NDF, ChiSquare(_tempSwap, fFitModulation)));
1225 Double_t CDFROOT(1.-ChiSquareCDF(NDF, fFitModulation->GetChisquare()));
1226 Double_t CDFKolmogorov(KolmogorovTest(_tempSwap, fFitModulation));
1227 // fill the values and centrality correlation (redundant but easy on the eyes)
1228 fHistPvalueCDF->Fill(CDF);
1229 fHistPvalueCDFCent->Fill(fCent, CDF);
1230 fHistPvalueCDFROOT->Fill(CDFROOT);
1231 fHistPvalueCDFROOTCent->Fill(fCent, CDFROOT);
1232 fHistKolmogorovTest->Fill(CDFKolmogorov);
1233 fHistChi2ROOTCent->Fill(fCent, fFitModulation->GetChisquare()/((float)NDF));
1234 fHistChi2Cent->Fill(fCent, ChiSquare(_tempSwap, fFitModulation)/((float)NDF));
1235 fHistKolmogorovTestCent->Fill(fCent, CDFKolmogorov);
1236 fHistPChi2Root->Fill(CDFROOT, fFitModulation->GetChisquare()/((float)NDF));
1237 fHistPChi2->Fill(CDF, ChiSquare(_tempSwap, fFitModulation)/((float)NDF));
1238 fHistPKolmogorov->Fill(CDF, CDFKolmogorov);
1240 // variable CDF is used for making cuts, so we fill it with the selected p-value
1241 switch (fFitGoodnessTest) {
1245 case kChi2Poisson : break; // CDF is already CDF
1246 case kKolmogorov : {
1247 CDF = CDFKolmogorov;
1253 // as an additional quality check, see if fitting a control fit has a higher significance
1254 _tempSwap.Fit(fFitControl, fFitModulationOptions.Data(), "", lowBound, upBound);
1255 Double_t CDFControl(-1.);
1256 switch (fFitGoodnessTest) {
1258 CDFControl = 1.-ChiSquareCDF(fFitControl->GetNDF(), fFitModulation->GetChisquare());
1260 case kChi2Poisson : {
1261 CDFControl = 1.-ChiSquareCDF(fFitControl->GetNDF(), ChiSquare(_tempSwap, fFitModulation));
1263 case kKolmogorov : {
1264 CDFControl = KolmogorovTest(_tempSwap, fFitControl);
1268 if(CDFControl > CDF) {
1269 CDF = -1.; // control fit is more significant, so throw out the 'old' fit
1270 fHistRhoStatusCent->Fill(fCent, -1);
1273 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
1274 // part of phase space, the requirement that energy desntiy is larger than zero is applied
1275 // to the FULL spectrum
1276 fHistRhoStatusCent->Fill(fCent, 0.);
1277 // for LOCAL didactic purposes, save the best and the worst fits
1278 // this routine can produce a lot of output histograms (it's not memory 'safe') and will not work on GRID
1279 // since the output will become unmergeable (i.e. different nodes may produce conflicting output)
1280 switch (fRunModeType) {
1282 if(fRandom->Uniform(0, 100) > fPercentageOfFits) break;
1283 static Int_t didacticCounterBest(0);
1284 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
1285 TF1* didacticFit = (TF1*)fFitModulation->Clone(Form("fit_%i_CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
1286 switch(fFitModulationType) {
1288 // to make a nice picture also plot the separate components (v2 and v3) of the fit
1289 // only done for cobined fit where there are actually components to split ...
1290 TF1* v0(new TF1("dfit_kV2", "[0]", 0, TMath::TwoPi()));
1291 v0->SetParameter(0, didacticFit->GetParameter(0)); // normalization
1292 v0->SetLineColor(kMagenta);
1293 v0->SetLineStyle(7);
1294 didacticProfile->GetListOfFunctions()->Add(v0);
1295 TF1* v2(new TF1("dfit_kV2", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
1296 v2->SetParameter(0, didacticFit->GetParameter(0)); // normalization
1297 v2->SetParameter(3, didacticFit->GetParameter(3)); // v2
1298 v2->FixParameter(1, 1.); // constant
1299 v2->FixParameter(2, 2.); // constant
1300 v2->FixParameter(4, didacticFit->GetParameter(4)); // psi2
1301 v2->SetLineColor(kGreen);
1302 didacticProfile->GetListOfFunctions()->Add(v2);
1303 TF1* v3(new TF1("dfit_kV3", "[0]*([1]+[2]*[3]*TMath::Cos([5]*(x-[4])))", 0, TMath::TwoPi()));
1304 v3->SetParameter(0, didacticFit->GetParameter(0)); // normalization
1305 v3->SetParameter(3, didacticFit->GetParameter(7)); // v3
1306 v3->FixParameter(1, 1.); // constant
1307 v3->FixParameter(2, 2.); // constant
1308 v3->FixParameter(4, didacticFit->GetParameter(6)); // psi3
1309 v3->FixParameter(5, 3.); // constant
1310 v3->SetLineColor(kCyan);
1311 didacticProfile->GetListOfFunctions()->Add(v3);
1315 didacticProfile->GetListOfFunctions()->Add(didacticFit);
1316 didacticProfile->GetYaxis()->SetTitle("#frac{d #sum #it{p}_{T}}{d #varphi} [GeV/#it{c}]");
1317 didacticProfile->GetXaxis()->SetTitle("#varphi");
1318 fOutputListGood->Add(didacticProfile);
1319 didacticCounterBest++;
1320 TH2F* didacticSurface = BookTH2F(Form("surface_%s", didacticProfile->GetName()), "#phi", "#eta", 50, 0, TMath::TwoPi(), 50, -1, 1, -1, kFALSE);
1321 for(Int_t i(0); i < iTracks; i++) {
1322 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1323 if(PassesCuts(track)) {
1324 if(fUsePtWeight) didacticSurface->Fill(track->Phi(), track->Eta(), track->Pt());
1325 else didacticSurface->Fill(track->Phi(), track->Eta());
1328 if(fExcludeLeadingJetsFromFit) { // visualize the excluded region
1329 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);
1330 f2->SetParameters(excludeInPt/3.,excludeInPhi,.1,excludeInEta,.1);
1331 didacticSurface->GetListOfFunctions()->Add(f2);
1333 fOutputListGood->Add(didacticSurface);
1337 } else { // if the fit is of poor quality revert to the original rho estimate
1338 switch (fRunModeType) { // again see if we want to save the fit
1340 static Int_t didacticCounterWorst(0);
1341 if(fRandom->Uniform(0, 100) > fPercentageOfFits) break;
1342 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data() ));
1343 TF1* didacticFit = (TF1*)fFitModulation->Clone(Form("fit_%i_p_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data()));
1344 didacticProfile->GetListOfFunctions()->Add(didacticFit);
1345 fOutputListBad->Add(didacticProfile);
1346 didacticCounterWorst++;
1350 switch (fFitModulationType) {
1351 case kNoFit : break; // nothing to do
1352 case kCombined : fFitModulation->SetParameter(7, 0); // no break
1353 case kFourierSeries : fFitModulation->SetParameter(7, 0); // no break
1354 default : { // needs to be done if there was a poor fit
1355 fFitModulation->SetParameter(3, 0);
1356 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1359 if(CDF > -.5) fHistRhoStatusCent->Fill(fCent, 1.);
1360 return kFALSE; // return false if the fit is rejected
1364 //_____________________________________________________________________________
1365 Bool_t AliAnalysisTaskJetV2::PassesCuts(AliVEvent* event)
1368 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1369 // determine the run number to see if the track and jet cuts should be refreshed for semi-good TPC runs
1370 // only done if the runnumber changes, could be moved to a call to AliAnalysisTaskSE::Notify()
1371 if(fRunNumber != InputEvent()->GetRunNumber()) {
1372 fRunNumber = InputEvent()->GetRunNumber(); // set the current run number
1373 if(fDebug > 0) printf("__FUNC__ %s > NEW RUNNUMBER DETECTED \n ", __func__);
1374 // reset the cuts. should be a pointless operation except for the case where the run number changes
1375 // from semi-good back to good on one node, which is not a likely scenario (unless trains will
1376 // run as one masterjob)
1377 AliAnalysisTaskEmcal::SetTrackPhiLimits(-10., 10.);
1378 switch (fAnalysisType) {
1380 AliAnalysisTaskEmcalJet::SetJetPhiLimits(-10., 10.);
1383 AliAnalysisTaskEmcalJet::SetJetPhiLimits(1.405 + GetJetRadius(), 3.135 - GetJetRadius());
1387 if(fCachedRho) { // if there's a cached rho, it's the default, so switch back
1388 if(fDebug > 0) printf("__FUNC__ %s > replacing rho with cached rho \n ", __func__);
1389 fRho = fCachedRho; // reset rho back to cached value. again, should be pointless
1391 Bool_t flaggedAsSemiGood(kFALSE); // not flagged as anything
1392 for(Int_t i(0); i < fExpectedSemiGoodRuns->GetSize(); i++) {
1393 if(fExpectedSemiGoodRuns->At(i) == fRunNumber) { // run is semi-good
1394 if(fDebug > 0) printf("__FUNC__ %s > semi-good tpc run detected, adjusting acceptance \n ", __func__);
1395 flaggedAsSemiGood = kTRUE;
1396 switch (fAnalysisType) {
1397 // for full jets the jet acceptance does not have to be changed as emcal does not
1398 // cover the tpc low voltage readout strips
1400 AliAnalysisTaskEmcalJet::SetJetPhiLimits(fSemiGoodJetMinPhi, fSemiGoodJetMaxPhi); // just an acceptance cut, jets are obtained from full azimuth, so no edge effects
1404 AliAnalysisTaskEmcal::SetTrackPhiLimits(fSemiGoodTrackMinPhi, fSemiGoodTrackMaxPhi); // only affects vn extraction, NOT jet finding
1405 // for semi-good runs, also try to get the 'small rho' estimate, if it is available
1406 AliRhoParameter* tempRho(dynamic_cast<AliRhoParameter*>(InputEvent()->FindListObject(fNameSmallRho.Data())));
1408 if(fDebug > 0) printf("__FUNC__ %s > switching to small rho, caching normal rho \n ", __func__);
1409 fHistAnalysisSummary->SetBinContent(54, 1.); // bookkeep the fact that small rho is used
1410 fCachedRho = fRho; // cache the original rho ...
1411 fRho = tempRho; // ... and use the small rho
1415 if(!flaggedAsSemiGood) {
1416 // in case the run is not a semi-good run, check if it is recognized as another run
1417 // only done to catch unexpected runs
1418 for(Int_t i(0); i < fExpectedRuns->GetSize(); i++) {
1419 if(fExpectedRuns->At(i) == fRunNumber) break; // run is known, break the loop else store the number in a random bin
1420 fHistUndeterminedRunQA->SetBinContent(TMath::Nint(10.*gRandom->Uniform(0.,.9))+1, fRunNumber);
1422 fHistAnalysisSummary->SetBinContent(53, 1.); // bookkeep which rho estimate is used
1425 // continue with event selection
1426 if(!event || !AliAnalysisTaskEmcal::IsEventSelected()) return kFALSE;
1427 if(TMath::Abs(InputEvent()->GetPrimaryVertex()->GetZ()) > 10.) return kFALSE;
1428 // aod and esd specific checks
1429 switch (fDataType) {
1431 AliESDEvent* esdEvent = static_cast<AliESDEvent*>(InputEvent());
1432 if( (!esdEvent) || (TMath::Abs(esdEvent->GetPrimaryVertexSPD()->GetZ() - esdEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE;
1435 AliAODEvent* aodEvent = static_cast<AliAODEvent*>(InputEvent());
1436 if( (!aodEvent) || (TMath::Abs(aodEvent->GetPrimaryVertexSPD()->GetZ() - aodEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE;
1440 fCent = InputEvent()->GetCentrality()->GetCentralityPercentile("V0M");
1441 if(fCent <= fCentralityClasses->At(0) || fCent >= fCentralityClasses->At(fCentralityClasses->GetSize()-1) || TMath::Abs(fCent-InputEvent()->GetCentrality()->GetCentralityPercentile("TRK")) > 5.) return kFALSE;
1442 // determine centrality class
1443 fInCentralitySelection = -1;
1444 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
1445 if(fCent >= fCentralityClasses->At(i) && fCent <= fCentralityClasses->At(1+i)) {
1446 fInCentralitySelection = i;
1449 if(fInCentralitySelection<0) return kFALSE; // should be null op
1450 if(fExplicitOutlierCut == 2010 || fExplicitOutlierCut == 2011) {
1451 if(!PassesCuts(fExplicitOutlierCut)) return kFALSE;
1453 // see if input containers are filled
1454 if(fTracks->GetEntries() < 1) return kFALSE;
1455 if(fRho->GetVal() <= 0 ) return kFALSE;
1456 if(fAnalysisType == AliAnalysisTaskJetV2::kFull && !fClusterCont) return kFALSE;
1459 //_____________________________________________________________________________
1460 Bool_t AliAnalysisTaskJetV2::PassesCuts(Int_t year)
1462 // additional centrality cut based on relation between tpc and global multiplicity
1463 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1464 AliAODEvent* event(dynamic_cast<AliAODEvent*>(InputEvent()));
1465 if(!event) return kFALSE;
1466 Int_t multTPC(0), multGlob(0), nTracks(InputEvent()->GetNumberOfTracks());
1467 for(Int_t iTracks = 0; iTracks < nTracks; iTracks++) {
1468 AliAODTrack* track = event->GetTrack(iTracks);
1469 if(!track) continue;
1470 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
1471 if (track->TestFilterBit(1) && track->Chi2perNDF() > 0.2) multTPC++;
1472 if (!track->TestFilterBit(16) || track->Chi2perNDF() < 0.1) continue;
1473 Double_t b[2] = {-99., -99.};
1474 Double_t bCov[3] = {-99., -99., -99.};
1475 AliAODTrack copy(*track);
1476 if (copy.PropagateToDCA(event->GetPrimaryVertex(), event->GetMagneticField(), 100., b, bCov) && TMath::Abs(b[0]) < 0.3 && TMath::Abs(b[1]) < 0.3) multGlob++;
1478 if(year == 2010 && multTPC > (-40.3+1.22*multGlob) && multTPC < (32.1+1.59*multGlob)) return kTRUE;
1479 if(year == 2011 && multTPC > (-36.73 + 1.48*multGlob) && multTPC < (62.87 + 1.78*multGlob)) return kTRUE;
1482 //_____________________________________________________________________________
1483 void AliAnalysisTaskJetV2::FillHistogramsAfterSubtraction(Double_t psi2, Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc)
1486 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1487 FillTrackHistograms();
1488 if(fAnalysisType == AliAnalysisTaskJetV2::kFull) FillClusterHistograms();
1489 FillJetHistograms(psi2);
1490 if(fFillQAHistograms) FillEventPlaneHistograms(vzero, vzeroComb, tpc);
1491 FillRhoHistograms();
1492 FillDeltaPtHistograms(psi2);
1494 //_____________________________________________________________________________
1495 void AliAnalysisTaskJetV2::FillTrackHistograms() const
1497 // fill track histograms
1498 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1499 Int_t iTracks(fTracks->GetEntriesFast()), iAcceptedTracks(0);
1500 for(Int_t i(0); i < iTracks; i++) {
1501 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1502 if(!PassesCuts(track)) continue;
1504 fHistPicoTrackPt[fInCentralitySelection]->Fill(track->Pt());
1505 if(fFillQAHistograms) FillQAHistograms(track);
1507 fHistPicoTrackMult[fInCentralitySelection]->Fill(iAcceptedTracks);
1509 //_____________________________________________________________________________
1510 void AliAnalysisTaskJetV2::FillClusterHistograms() const
1512 // fill cluster histograms
1513 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1514 if(!fClusterCont) return;
1515 Int_t iClusters(fClusterCont->GetNClusters());
1516 for(Int_t i(0); i < iClusters; i++) {
1517 AliVCluster* cluster = fClusterCont->GetCluster(i);
1518 if (!PassesCuts(cluster)) continue;
1519 TLorentzVector clusterLorentzVector;
1520 cluster->GetMomentum(clusterLorentzVector, const_cast<Double_t*>(fVertex));
1521 fHistClusterPt[fInCentralitySelection]->Fill(clusterLorentzVector.Pt());
1522 fHistClusterEtaPhi[fInCentralitySelection]->Fill(clusterLorentzVector.Eta(), clusterLorentzVector.Phi());
1523 fHistClusterEtaPhiWeighted[fInCentralitySelection]->Fill(clusterLorentzVector.Eta(), clusterLorentzVector.Phi(), clusterLorentzVector.Pt());
1527 //_____________________________________________________________________________
1528 void AliAnalysisTaskJetV2::FillEventPlaneHistograms(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc) const
1530 // fill event plane histograms
1531 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1532 fHistPsiControl->Fill(0.5, vzero[0][0]); // vzero a psi2
1533 fHistPsiControl->Fill(1.5, vzero[1][0]); // vzero c psi2
1534 fHistPsiControl->Fill(2.5, tpc[0]); // tpc psi 2
1535 fHistPsiControl->Fill(5.5, vzero[0][1]); // vzero a psi3
1536 fHistPsiControl->Fill(6.5, vzero[1][1]); // vzero b psi3
1537 fHistPsiControl->Fill(7.5, tpc[1]); // tpc psi 3
1538 fHistPsiVZEROA->Fill(vzero[0][0]);
1539 fHistPsiVZEROC->Fill(vzero[1][0]);
1540 fHistPsiVZERO->Fill(vzeroComb[0]);
1541 fHistPsiTPC->Fill(tpc[0]);
1542 fHistPsiSpread->Fill(0.5, TMath::Abs(vzero[0][0]-vzero[1][0]));
1543 fHistPsiSpread->Fill(1.5, TMath::Abs(vzero[0][0]-tpc[0]));
1544 fHistPsiSpread->Fill(2.5, TMath::Abs(vzero[1][0]-tpc[0]));
1545 // event plane vs centrality QA histo's to check recentering
1546 Double_t TRK(InputEvent()->GetCentrality()->GetCentralityPercentile("TRK"));
1547 Double_t V0M(InputEvent()->GetCentrality()->GetCentralityPercentile("V0M"));
1548 fHistPsiVZEROAV0M->Fill(V0M, vzero[0][0]);
1549 fHistPsiVZEROCV0M->Fill(V0M, vzero[1][0]);
1550 fHistPsiVZEROVV0M->Fill(V0M, vzeroComb[0]);
1551 fHistPsiTPCiV0M->Fill(V0M, tpc[0]);
1552 fHistPsiVZEROATRK->Fill(TRK, vzero[0][0]);
1553 fHistPsiVZEROCTRK->Fill(TRK, vzero[1][0]);
1554 fHistPsiVZEROTRK->Fill(TRK, vzeroComb[0]);
1555 fHistPsiTPCTRK->Fill(TRK, tpc[0]);
1557 //_____________________________________________________________________________
1558 void AliAnalysisTaskJetV2::FillRhoHistograms()
1560 // fill rho histograms
1561 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1562 fHistRhoPackage[fInCentralitySelection]->Fill(fLocalRho->GetVal()); // save the rho estimate from the emcal jet package
1563 // get multiplicity FIXME inefficient
1564 Int_t iJets(fJets->GetEntriesFast());
1565 Double_t rho(fLocalRho->GetLocalVal(TMath::Pi(), TMath::Pi(), fLocalRho->GetVal()));
1566 fHistRho[fInCentralitySelection]->Fill(rho);
1567 fHistRhoVsMult->Fill(fTracks->GetEntries(), rho);
1568 fHistRhoVsCent->Fill(fCent, rho);
1569 for(Int_t i(0); i < iJets; i++) {
1570 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
1571 if(!PassesCuts(jet)) continue;
1572 fHistRhoAVsMult->Fill(fTracks->GetEntries(), rho * jet->Area());
1573 fHistRhoAVsCent->Fill(fCent, rho * jet->Area());
1576 //_____________________________________________________________________________
1577 void AliAnalysisTaskJetV2::FillDeltaPtHistograms(Double_t psi2) const
1579 // fill delta pt histograms
1580 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1582 const Float_t areaRC = GetJetRadius()*GetJetRadius()*TMath::Pi();
1583 // we're retrieved the leading jet, now get a random cone
1584 for(i = 0; i < fMaxCones; i++) {
1585 Float_t pt(0), eta(0), phi(0);
1586 // get a random cone without constraints on leading jet position
1587 CalculateRandomCone(pt, eta, phi, fTracksCont, fClusterCont, 0x0);
1589 if(fFillQAHistograms) fHistRCPhiEta[fInCentralitySelection]->Fill(phi, eta);
1590 fHistRhoVsRCPt[fInCentralitySelection]->Fill(pt, fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal())*areaRC);
1591 fHistRCPt[fInCentralitySelection]->Fill(pt);
1592 fHistDeltaPtDeltaPhi2[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1593 fHistDeltaPtDeltaPhi2Rho0[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetVal());
1596 // get a random cone excluding leading jet area
1597 CalculateRandomCone(pt, eta, phi, fTracksCont, fClusterCont, fLeadingJet);
1599 if(fFillQAHistograms) fHistRCPhiEtaExLJ[fInCentralitySelection]->Fill(phi, eta);
1600 fHistRhoVsRCPtExLJ[fInCentralitySelection]->Fill(pt, fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal())*areaRC);
1601 fHistRCPtExLJ[fInCentralitySelection]->Fill(pt);
1602 fHistDeltaPtDeltaPhi2ExLJ[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1603 fHistDeltaPtDeltaPhi2ExLJRho0[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetVal());
1607 //_____________________________________________________________________________
1608 void AliAnalysisTaskJetV2::FillJetHistograms(Double_t psi2)
1610 // fill jet histograms
1611 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1612 Int_t iJets(fJets->GetEntriesFast());
1613 for(Int_t i(0); i < iJets; i++) {
1614 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
1615 if(PassesCuts(jet)) {
1616 Double_t pt(jet->Pt()), area(jet->Area()), eta(jet->Eta()), phi(jet->Phi());
1617 Double_t rho(fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1618 fHistJetPtRaw[fInCentralitySelection]->Fill(pt);
1619 fHistJetPt[fInCentralitySelection]->Fill(pt-area*rho);
1620 if(fFillQAHistograms) fHistJetEtaPhi[fInCentralitySelection]->Fill(eta, phi);
1621 fHistJetPtArea[fInCentralitySelection]->Fill(pt-area*rho, area);
1622 fHistJetPtEta[fInCentralitySelection]->Fill(pt-area*rho, eta);
1623 fHistJetPsi2Pt[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt-area*rho);
1624 fHistJetPsi2PtRho0[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt-area*fLocalRho->GetVal());
1625 fHistJetPtConstituents[fInCentralitySelection]->Fill(pt-area*rho, jet->Nch());
1626 fHistJetEtaRho[fInCentralitySelection]->Fill(eta, pt/area);
1630 //_____________________________________________________________________________
1631 void AliAnalysisTaskJetV2::FillQAHistograms(AliVTrack* vtrack) const
1633 // fill qa histograms for pico tracks
1635 AliPicoTrack* track = static_cast<AliPicoTrack*>(vtrack);
1636 fHistRunnumbersPhi->Fill(fMappedRunNumber, track->Phi());
1637 fHistRunnumbersEta->Fill(fMappedRunNumber, track->Eta());
1638 Int_t type((int)(track->GetTrackType()));
1641 fHistPicoCat1[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1644 fHistPicoCat2[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1647 fHistPicoCat3[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1652 //_____________________________________________________________________________
1653 void AliAnalysisTaskJetV2::FillQAHistograms(AliVEvent* vevent)
1655 // fill qa histograms for events
1657 fHistVertexz->Fill(vevent->GetPrimaryVertex()->GetZ());
1658 fHistCentrality->Fill(fCent);
1659 Int_t runNumber(InputEvent()->GetRunNumber());
1660 for(fMappedRunNumber = 0; fExpectedRuns->GetSize()+1; fMappedRunNumber++) {
1661 if(fExpectedRuns->At(fMappedRunNumber) == runNumber) break;
1664 //_____________________________________________________________________________
1665 void AliAnalysisTaskJetV2::FillAnalysisSummaryHistogram() const
1667 // fill the analysis summary histrogram, saves all relevant analysis settigns
1668 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1669 fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fJetRadius");
1670 fHistAnalysisSummary->SetBinContent(2, GetJetContainer()->GetJetRadius());
1671 fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fJetEtaMin");
1672 fHistAnalysisSummary->SetBinContent(3, GetJetContainer()->GetJetEtaMin());
1673 fHistAnalysisSummary->GetXaxis()->SetBinLabel(4, "fJetEtaMax");
1674 fHistAnalysisSummary->SetBinContent(4, GetJetContainer()->GetJetEtaMax());
1675 fHistAnalysisSummary->GetXaxis()->SetBinLabel(5, "fJetPhiMin");
1676 fHistAnalysisSummary->SetBinContent(5, GetJetContainer()->GetJetPhiMin());
1677 fHistAnalysisSummary->GetXaxis()->SetBinLabel(6, "fJetPhiMax");
1678 fHistAnalysisSummary->SetBinContent(6, GetJetContainer()->GetJetPhiMin());
1679 fHistAnalysisSummary->GetXaxis()->SetBinLabel(16, "fForceBeamType");
1680 fHistAnalysisSummary->SetBinContent(16, fForceBeamType);
1681 fHistAnalysisSummary->GetXaxis()->SetBinLabel(17, "fMinCent");
1682 fHistAnalysisSummary->SetBinContent(17, fMinCent);
1683 fHistAnalysisSummary->GetXaxis()->SetBinLabel(18, "fMaxCent");
1684 fHistAnalysisSummary->SetBinContent(18, fMaxCent);
1685 fHistAnalysisSummary->GetXaxis()->SetBinLabel(19, "fMinVz");
1686 fHistAnalysisSummary->SetBinContent(19, fMinVz);
1687 fHistAnalysisSummary->GetXaxis()->SetBinLabel(20, "fMaxVz");
1688 fHistAnalysisSummary->SetBinContent(20, fMaxVz);
1689 fHistAnalysisSummary->GetXaxis()->SetBinLabel(21, "fOffTrigger");
1690 fHistAnalysisSummary->SetBinContent(21, fOffTrigger);
1691 fHistAnalysisSummary->GetXaxis()->SetBinLabel(34, "fitModulationType");
1692 fHistAnalysisSummary->SetBinContent(34, (int)fFitModulationType);
1693 fHistAnalysisSummary->GetXaxis()->SetBinLabel(35, "runModeType");
1694 fHistAnalysisSummary->SetBinContent(35, (int)fRunModeType);
1695 fHistAnalysisSummary->GetXaxis()->SetBinLabel(36, "data type");
1696 fHistAnalysisSummary->SetBinContent(36, (int)fDataType);
1697 fHistAnalysisSummary->GetXaxis()->SetBinLabel(37, "iterator");
1698 fHistAnalysisSummary->SetBinContent(37, 1.);
1699 fHistAnalysisSummary->GetXaxis()->SetBinLabel(38, "fMinPvalue");
1700 fHistAnalysisSummary->SetBinContent(38, fMinPvalue);
1701 fHistAnalysisSummary->GetXaxis()->SetBinLabel(39, "fMaxPvalue");
1702 fHistAnalysisSummary->SetBinContent(39, fMaxPvalue);
1703 fHistAnalysisSummary->GetXaxis()->SetBinLabel(40, "fExcludeLeadingJetsFromFit");
1704 fHistAnalysisSummary->SetBinContent(40, fExcludeLeadingJetsFromFit);
1705 fHistAnalysisSummary->GetXaxis()->SetBinLabel(41, "fRebinSwapHistoOnTheFly");
1706 fHistAnalysisSummary->SetBinContent(41, (int)fRebinSwapHistoOnTheFly);
1707 fHistAnalysisSummary->GetXaxis()->SetBinLabel(42, "fUsePtWeight");
1708 fHistAnalysisSummary->SetBinContent(42, (int)fUsePtWeight);
1709 fHistAnalysisSummary->GetXaxis()->SetBinLabel(43, "fExplicitOutlierCut");
1710 fHistAnalysisSummary->SetBinContent(43, fExplicitOutlierCut);
1711 fHistAnalysisSummary->GetXaxis()->SetBinLabel(44, "fSoftTrackMinPt");
1712 fHistAnalysisSummary->SetBinContent(44, fSoftTrackMinPt);
1713 fHistAnalysisSummary->GetXaxis()->SetBinLabel(45, "fSoftTrackMaxPt");
1714 fHistAnalysisSummary->SetBinContent(45, fSoftTrackMaxPt);
1715 fHistAnalysisSummary->GetXaxis()->SetBinLabel(46, "fMaxCones");
1716 fHistAnalysisSummary->SetBinContent(46, fMaxCones);
1717 fHistAnalysisSummary->GetXaxis()->SetBinLabel(47, "used rho");
1718 fHistAnalysisSummary->GetXaxis()->SetBinLabel(48, "used small rho");
1720 //_____________________________________________________________________________
1721 void AliAnalysisTaskJetV2::Terminate(Option_t *)
1724 switch (fRunModeType) {
1726 printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1727 AliAnalysisTaskJetV2::Dump();
1728 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));
1733 //_____________________________________________________________________________
1734 void AliAnalysisTaskJetV2::SetModulationFit(TF1* fit)
1736 // set modulation fit
1737 if (fFitModulation) delete fFitModulation;
1738 fFitModulation = fit;
1740 //_____________________________________________________________________________
1741 void AliAnalysisTaskJetV2::SetUseControlFit(Bool_t c)
1744 if (fFitControl) delete fFitControl;
1746 fFitControl = new TF1("controlFit", "pol0", 0, TMath::TwoPi());
1747 } else fFitControl = 0x0;
1749 //_____________________________________________________________________________
1750 TH1F* AliAnalysisTaskJetV2::GetResolutionFromOuptutFile(detectorType det, Int_t h, TArrayD* cen)
1752 // INTERFACE METHOD FOR OUTPUTFILE
1753 // get the detector resolution, user has ownership of the returned histogram
1755 printf(" > Please add fOutputList first < \n");
1759 (cen) ? r = new TH1F("R", "R", cen->GetSize()-1, cen->GetArray()) : r = new TH1F("R", "R", 10, 0, 10);
1760 if(!cen) r->GetXaxis()->SetTitle("number of centrality bin");
1761 r->GetYaxis()->SetTitle(Form("Resolution #Psi_{%i}", h));
1762 for(Int_t i(0); i < 10; i++) {
1763 TProfile* temp((TProfile*)fOutputList->FindObject(Form("fProfV%iResolution_%i", h, i)));
1765 Double_t a(temp->GetBinContent(3)), b(temp->GetBinContent(5)), c(temp->GetBinContent(7));
1766 Double_t d(temp->GetBinContent(9)), e(temp->GetBinContent(10)), f(temp->GetBinContent(11));
1767 Double_t _a(temp->GetBinError(3)), _b(temp->GetBinError(5)), _c(temp->GetBinError(7));
1768 Double_t _d(temp->GetBinError(9)), _e(temp->GetBinError(10)), _f(temp->GetBinError(11));
1769 if(a <= 0 || b <= 0 || c <= 0 || d <= 0 || e <= 0 || f <= 0) continue;
1772 r->SetBinContent(1+i, TMath::Sqrt((a*b)/c));
1773 if(i==0) r->SetNameTitle("VZEROA resolution", "VZEROA resolution");
1774 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
1777 r->SetBinContent(1+i, TMath::Sqrt((a*c)/b));
1778 if(i==0) r->SetNameTitle("VZEROC resolution", "VZEROC resolution");
1779 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
1782 r->SetBinContent(1+i, TMath::Sqrt((b*c)/a));
1783 if(i==0) r->SetNameTitle("TPC resolution", "TPC resolution");
1784 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
1787 r->SetBinContent(1+i, TMath::Sqrt((d*e)/f));
1788 if(i==0) r->SetNameTitle("VZEROComb resolution", "VZEROComb resolution");
1789 r->SetBinError(1+i, TMath::Sqrt(_d*_d+_e*_e+_f*_f));
1796 //_____________________________________________________________________________
1797 TH1F* AliAnalysisTaskJetV2::CorrectForResolutionDiff(TH1F* v, detectorType det, TArrayD* cen, Int_t c, Int_t h)
1799 // INTERFACE METHOD FOR OUTPUT FILE
1800 // correct the supplied differential vn histogram v for detector resolution
1801 TH1F* r(GetResolutionFromOuptutFile(det, h, cen));
1803 printf(" > Couldn't find resolution < \n");
1806 Double_t res(1./r->GetBinContent(1+r->FindBin(c)));
1807 TF1* line = new TF1("line", "pol0", 0, 200);
1808 line->SetParameter(0, res);
1812 //_____________________________________________________________________________
1813 TH1F* AliAnalysisTaskJetV2::CorrectForResolutionInt(TH1F* v, detectorType det, TArrayD* cen, Int_t h)
1815 // INTERFACE METHOD FOR OUTPUT FILE
1816 // correct the supplied intetrated vn histogram v for detector resolution
1817 // integrated vn must have the same centrality binning as the resolotion correction
1818 TH1F* r(GetResolutionFromOuptutFile(det, h, cen));
1822 //_____________________________________________________________________________
1823 TH1F* AliAnalysisTaskJetV2::GetDifferentialQC(TProfile* refCumulants, TProfile* diffCumlants, TArrayD* ptBins, Int_t h)
1825 // get differential QC
1826 Double_t r(refCumulants->GetBinContent(h-1)); // v2 reference flow
1827 if(r > 0) r = TMath::Sqrt(r);
1828 TH1F* qc = new TH1F(Form("QC2v%i", h), Form("QC2v%i", h), ptBins->GetSize()-1, ptBins->GetArray());
1829 Double_t a(0), b(0), c(0); // dummy variables
1830 for(Int_t i(0); i < ptBins->GetSize(); i++) {
1832 a = diffCumlants->GetBinContent(1+i);
1833 b = diffCumlants->GetBinError(1+i);
1835 qc->SetBinContent(1+i, c);
1836 (a <= 0 || b <= 0) ? qc->SetBinError(1+i, b) : qc->SetBinError(1+i, TMath::Sqrt(c*c*b*b/(a*a)));
1842 //_____________________________________________________________________________