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 **************************************************************************/
17 * analysis task for jet flow preparation
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
42 #include <AliAnalysisTask.h>
43 #include <AliAnalysisManager.h>
44 #include <AliCentrality.h>
45 #include <AliVVertex.h>
46 #include <AliESDEvent.h>
47 #include <AliAODEvent.h>
48 #include <AliAODTrack.h>
50 #include <AliPicoTrack.h>
51 #include <AliEmcalJet.h>
52 #include <AliRhoParameter.h>
54 #include "AliAnalysisTaskRhoVnModulation.h"
57 class AliAnalysisTaskRhoVnModulation;
60 ClassImp(AliAnalysisTaskRhoVnModulation)
62 AliAnalysisTaskRhoVnModulation::AliAnalysisTaskRhoVnModulation() : AliAnalysisTaskEmcalJet("AliAnalysisTaskRhoVnModulation", kTRUE),
63 fDebug(0), fInitialized(0), fFillQAHistograms(kTRUE), fCentralityClasses(0), fUserSuppliedV2(0), fUserSuppliedV3(0), fNAcceptedTracks(0), fFitModulationType(kNoFit), fUsePtWeight(kTRUE), fDetectorType(kTPC), fFitModulationOptions("Q"), fRunModeType(kGrid), fDataType(kESD), fRandom(0), fMappedRunNumber(0), fInCentralitySelection(-1), fFitModulation(0), fMinPvalue(0), fMaxPvalue(1), fNameJetClones(0), fNamePicoTrackClones(0), fNameRho(0), fAbsVertexZ(10), fHistCentrality(0), fHistVertexz(0), fHistRunnumbersPhi(0), fHistRunnumbersEta(0), fHistPvaluePDF(0), fHistPvalueCDF(0), fMinDisanceRCtoLJ(0), fRandomConeRadius(0.4), fAbsVnHarmonics(kTRUE), fExcludeLeadingJetsFromFit(1.), fRebinSwapHistoOnTheFly(kTRUE), fPercentageOfFits(10.), fUseV0EventPlaneFromHeader(kFALSE), fSetPtSub(kFALSE), fExplicitOutlierCut(-1), fMinLeadingHadronPt(0), fOutputList(0), fOutputListGood(0), fOutputListBad(0), fHistAnalysisSummary(0), fHistSwap(0), fProfV2(0), fProfV3(0), fHistPsiControl(0), fHistPsiSpread(0), fHistPsiVZEROA(0), fHistPsiVZEROC(0), fHistPsiTPC(0), fHistPsiTPCSUBA(0), fHistPsiTPCSUBB(0),
64 fHistRhoVsMult(0), fHistRhoVsCent(0), fHistRhoAVsMult(0), fHistRhoAVsCent(0) {
65 for(Int_t i(0); i < 10; i++) {
66 fProfV2Resolution[i] = 0;
67 fProfV3Resolution[i] = 0;
68 fHistPicoTrackPt[i] = 0;
72 /* fHistClusterPt[i] = 0; */
73 /* fHistClusterPhi[i] = 0; */
74 /* fHistClusterEta[i] = 0; */
75 /* fHistClusterCorrPt[i] = 0; */
76 /* fHistClusterCorrPhi[i] = 0; */
77 /* fHistClusterCorrEta[i] = 0; */
78 fHistRhoPackage[i] = 0;
81 fHistRhoVsRCPt[i] = 0;
83 fHistDeltaPtDeltaPhi2[i] = 0;
84 fHistDeltaPtDeltaPhi3[i] = 0;
85 fHistRCPhiEtaExLJ[i] = 0;
86 fHistRhoVsRCPtExLJ[i] = 0;
88 fHistDeltaPtDeltaPhi2ExLJ[i] = 0;
89 fHistDeltaPtDeltaPhi3ExLJ[i] = 0;
90 /* fHistRCPhiEtaRand[i] = 0; */
91 /* fHistRhoVsRCPtRand[i] = 0; */
92 /* fHistRCPtRand[i] = 0; */
93 /* fHistDeltaPtDeltaPhi2Rand[i] = 0; */
94 /* fHistDeltaPtDeltaPhi3Rand[i] = 0; */
97 fHistJetEtaPhi[i] = 0;
98 fHistJetPtArea[i] = 0;
99 fHistJetPtConstituents[i] = 0;
100 fHistJetEtaRho[i] = 0;
101 fHistJetPsiTPCPt[i] = 0;
102 fHistJetPsiVZEROAPt[i] = 0;
103 fHistJetPsiVZEROCPt[i] = 0;
104 fHistDeltaPhi2VZEROA[i] = 0;
105 fHistDeltaPhi2VZEROC[i] = 0;
106 fHistDeltaPhi2TPC[i] = 0;
107 fHistDeltaPhi3VZEROA[i] = 0;
108 fHistDeltaPhi3VZEROC[i] = 0;
109 fHistDeltaPhi3TPC[i] = 0;
111 // default constructor
113 //_____________________________________________________________________________
114 AliAnalysisTaskRhoVnModulation::AliAnalysisTaskRhoVnModulation(const char* name, runModeType type) : AliAnalysisTaskEmcalJet(name, kTRUE),
115 fDebug(0), fInitialized(0), fFillQAHistograms(kTRUE), fCentralityClasses(0), fUserSuppliedV2(0), fUserSuppliedV3(0), fNAcceptedTracks(0), fFitModulationType(kNoFit), fUsePtWeight(kTRUE), fDetectorType(kTPC), fFitModulationOptions("Q"), fRunModeType(type), fDataType(kESD), fRandom(0), fMappedRunNumber(0), fInCentralitySelection(-1), fFitModulation(0), fMinPvalue(0), fMaxPvalue(1), fNameJetClones(0), fNamePicoTrackClones(0), fNameRho(0), fAbsVertexZ(10), fHistCentrality(0), fHistVertexz(0), fHistRunnumbersPhi(0), fHistRunnumbersEta(0), fHistPvaluePDF(0), fHistPvalueCDF(0), fMinDisanceRCtoLJ(0), fRandomConeRadius(0.4), fAbsVnHarmonics(kTRUE), fExcludeLeadingJetsFromFit(1.), fRebinSwapHistoOnTheFly(kTRUE), fPercentageOfFits(10.), fUseV0EventPlaneFromHeader(kFALSE), fSetPtSub(kFALSE), fExplicitOutlierCut(-1), fMinLeadingHadronPt(0), fOutputList(0), fOutputListGood(0), fOutputListBad(0), fHistAnalysisSummary(0), fHistSwap(0), fProfV2(0), fProfV3(0), fHistPsiControl(0), fHistPsiSpread(0), fHistPsiVZEROA(0), fHistPsiVZEROC(0), fHistPsiTPC(0), fHistPsiTPCSUBA(0), fHistPsiTPCSUBB(0),
116 fHistRhoVsMult(0), fHistRhoVsCent(0), fHistRhoAVsMult(0), fHistRhoAVsCent(0) {
117 for(Int_t i(0); i < 10; i++) {
118 fProfV2Resolution[i] = 0;
119 fProfV3Resolution[i] = 0;
120 fHistPicoTrackPt[i] = 0;
121 fHistPicoCat1[i] = 0;
122 fHistPicoCat2[i] = 0;
123 fHistPicoCat3[i] = 0;
124 /* fHistClusterPt[i] = 0; */
125 /* fHistClusterPhi[i] = 0; */
126 /* fHistClusterEta[i] = 0; */
127 /* fHistClusterCorrPt[i] = 0; */
128 /* fHistClusterCorrPhi[i] = 0; */
129 /* fHistClusterCorrEta[i] = 0; */
130 fHistRhoPackage[i] = 0;
132 fHistRCPhiEta[i] = 0;
133 fHistRhoVsRCPt[i] = 0;
135 fHistDeltaPtDeltaPhi2[i] = 0;
136 fHistDeltaPtDeltaPhi3[i] = 0;
137 fHistRCPhiEtaExLJ[i] = 0;
138 fHistRhoVsRCPtExLJ[i] = 0;
139 fHistRCPtExLJ[i] = 0;
140 fHistDeltaPtDeltaPhi2ExLJ[i] = 0;
141 fHistDeltaPtDeltaPhi3ExLJ[i] = 0;
142 /* fHistRCPhiEtaRand[i] = 0; */
143 /* fHistRhoVsRCPtRand[i] = 0; */
144 /* fHistRCPtRand[i] = 0; */
145 /* fHistDeltaPtDeltaPhi2Rand[i] = 0; */
146 /* fHistDeltaPtDeltaPhi3Rand[i] = 0; */
147 fHistJetPtRaw[i] = 0;
149 fHistJetEtaPhi[i] = 0;
150 fHistJetPtArea[i] = 0;
151 fHistJetPtConstituents[i] = 0;
152 fHistJetEtaRho[i] = 0;
153 fHistJetPsiTPCPt[i] = 0;
154 fHistJetPsiVZEROAPt[i] = 0;
155 fHistJetPsiVZEROCPt[i] = 0;
156 fHistDeltaPhi2VZEROA[i] = 0;
157 fHistDeltaPhi2VZEROC[i] = 0;
158 fHistDeltaPhi2TPC[i] = 0;
159 fHistDeltaPhi3VZEROA[i] = 0;
160 fHistDeltaPhi3VZEROC[i] = 0;
161 fHistDeltaPhi3TPC[i] = 0;
164 DefineInput(0, TChain::Class());
165 DefineOutput(1, TList::Class());
166 switch (fRunModeType) {
168 gStyle->SetOptFit(1);
169 DefineOutput(2, TList::Class());
170 DefineOutput(3, TList::Class());
172 default: fDebug = -1; // suppress debug info explicitely when not running locally
175 //_____________________________________________________________________________
176 AliAnalysisTaskRhoVnModulation::~AliAnalysisTaskRhoVnModulation()
179 if(fOutputList) delete fOutputList;
180 if(fOutputListGood) delete fOutputListGood;
181 if(fOutputListBad) delete fOutputListBad;
182 if(fFitModulation) delete fFitModulation;
183 if(fHistSwap) delete fHistSwap;
184 if(fCentralityClasses) delete fCentralityClasses;
186 //_____________________________________________________________________________
187 Bool_t AliAnalysisTaskRhoVnModulation::InitializeAnalysis()
189 // initialize the anaysis
190 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
191 if(fMinDisanceRCtoLJ==0) fMinDisanceRCtoLJ = .5*fJetRadius;
192 if(dynamic_cast<AliAODEvent*>(InputEvent())) fDataType = kAOD; // determine the datatype
193 else if(dynamic_cast<AliESDEvent*>(InputEvent())) fDataType = kESD;
194 fHistAnalysisSummary->SetBinContent(36, (int)fDataType);
195 if(!fRandom) fRandom = new TRandom3(0); // get a randomized if one hasn't been user-supplied
196 switch (fFitModulationType) {
197 case kNoFit : { SetModulationFit(new TF1("fix_kNoFit", "[0]", 0, TMath::TwoPi())); } break;
199 SetModulationFit(new TF1("fit_kV2", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
200 fFitModulation->SetParameter(0, 0.); // normalization
201 fFitModulation->SetParameter(3, 0.2); // v2
202 fFitModulation->FixParameter(1, 1.); // constant
203 fFitModulation->FixParameter(2, 2.); // constant
206 SetModulationFit(new TF1("fit_kV3", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
207 fFitModulation->SetParameter(0, 0.); // normalization
208 fFitModulation->SetParameter(3, 0.2); // v3
209 fFitModulation->FixParameter(1, 1.); // constant
210 fFitModulation->FixParameter(2, 3.); // constant
212 default : { // for the combined fit, the 'direct fourier series' or the user supplied vn values we use v2 and v3
213 SetModulationFit(new TF1("fit_kCombined", "[0]*([1]+[2]*([3]*TMath::Cos([2]*(x-[4]))+[7]*TMath::Cos([5]*(x-[6]))))", 0, TMath::TwoPi()));
214 fFitModulation->SetParameter(0, 0.); // normalization
215 fFitModulation->SetParameter(3, 0.2); // v2
216 fFitModulation->FixParameter(1, 1.); // constant
217 fFitModulation->FixParameter(2, 2.); // constant
218 fFitModulation->FixParameter(5, 3.); // constant
219 fFitModulation->SetParameter(7, 0.2); // v3
222 switch (fRunModeType) {
223 case kGrid : { fFitModulationOptions += "N0"; } break;
228 //_____________________________________________________________________________
229 TH1F* AliAnalysisTaskRhoVnModulation::BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Int_t c, Bool_t append)
231 // book a TH1F and connect it to the output container
232 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
233 if(!fOutputList) return 0x0;
235 if(c!=-1) { // format centrality dependent histograms accordingly
236 name = Form("%s_%i", name, c);
237 title += Form("_%i-%i", fCentralityClasses->At(c), fCentralityClasses->At(1+c));
239 title += Form(";%s;[counts]", x);
240 TH1F* histogram = new TH1F(name, title.Data(), bins, min, max);
242 if(append) fOutputList->Add(histogram);
245 //_____________________________________________________________________________
246 TH2F* AliAnalysisTaskRhoVnModulation::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)
248 // book a TH2F and connect it to the output container
249 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
250 if(!fOutputList) return 0x0;
252 if(c!=-1) { // format centrality dependent histograms accordingly
253 name = Form("%s_%i", name, c);
254 title += Form("_%i-%i", fCentralityClasses->At(c), fCentralityClasses->At(1+c));
256 title += Form(";%s;%s", x, y);
257 TH2F* histogram = new TH2F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy);
259 if(append) fOutputList->Add(histogram);
262 //_____________________________________________________________________________
263 void AliAnalysisTaskRhoVnModulation::UserCreateOutputObjects()
265 // create output objects
266 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
267 fOutputList = new TList();
268 fOutputList->SetOwner(kTRUE);
269 if(!fCentralityClasses) { // classes must be defined at this point
270 Int_t c[] = {0, 20, 40, 60, 80, 100};
271 fCentralityClasses = new TArrayI(sizeof(c)/sizeof(c[0]), c);
274 fHistCentrality = BookTH1F("fHistCentrality", "centrality", 102, -2, 100);
275 fHistVertexz = BookTH1F("fHistVertexz", "vertex z (cm)", 100, -12, 12);
277 // pico track kinematics
278 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
279 fHistPicoTrackPt[i] = BookTH1F("fHistPicoTrackPt", "p_{t} [GeV/c]", 100, 0, 50, i);
280 if(fFillQAHistograms) {
281 fHistPicoCat1[i] = BookTH2F("fHistPicoCat1", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
282 fHistPicoCat2[i] = BookTH2F("fHistPicoCat2", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
283 fHistPicoCat3[i] = BookTH2F("fHistPicoCat3", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
286 /* fHistClusterPt[i] = BookTH1F("fHistClusterPt", "p_{t} [GeV/c]", 100, 0, 100, i); */
287 /* fHistClusterPhi[i] = BookTH1F("fHistClusterPhi", "#phi", 100, 0, TMath::TwoPi(), i); */
288 /* fHistClusterEta[i] = BookTH1F("fHistClusterEta", "#eta", 100, -5, 5); */
290 // emcal kinematics after hadronic correction
291 /* fHistClusterCorrPt[i] = BookTH1F("fHistClusterCorrPt", "p_{t} [GeV/c]", 100, 0, 100, i); */
292 /* fHistClusterCorrPhi[i] = BookTH1F("fHistClusterCorrPhi", "#phi", 100, 0, TMath::TwoPi(), i); */
293 /* fHistClusterCorrEta[i] = BookTH1F("fHistClusterCorrEta", "#eta", 100, -5, 5, i); */
296 // event plane estimates and quality
297 fHistPsiControl = new TProfile("fHistPsiControl", "fHistPsiControl", 10, 0, 10);
298 fHistPsiControl->Sumw2();
299 fHistPsiSpread = new TProfile("fHistPsiSpread", "fHistPsiSpread", 4, 0, 4);
300 fHistPsiSpread->Sumw2();
301 fHistPsiControl->GetXaxis()->SetBinLabel(1, "<#Psi_{2, VZEROA}>");
302 fHistPsiControl->GetXaxis()->SetBinLabel(2, "<#Psi_{2, VZEROC}>");
303 fHistPsiControl->GetXaxis()->SetBinLabel(3, "<#Psi_{2, TPC}>");
304 fHistPsiControl->GetXaxis()->SetBinLabel(4, "<#Psi_{2, TPC, #eta < 0}>");
305 fHistPsiControl->GetXaxis()->SetBinLabel(5, "<#Psi_{2, TPC, #eta > 0}>");
306 fHistPsiControl->GetXaxis()->SetBinLabel(1, "<#Psi_{3, VZEROA}>");
307 fHistPsiControl->GetXaxis()->SetBinLabel(2, "<#Psi_{3, VZEROC}>");
308 fHistPsiControl->GetXaxis()->SetBinLabel(3, "<#Psi_{3, TPC}>");
309 fHistPsiControl->GetXaxis()->SetBinLabel(4, "<#Psi_{3, TPC, #eta < 0}>");
310 fHistPsiControl->GetXaxis()->SetBinLabel(5, "<#Psi_{3, TPC, #eta > 0}>");
311 fHistPsiSpread->GetXaxis()->SetBinLabel(1, "<#Psi_{2, VZEROA} - #Psi_{2, VZEROC}>");
312 fHistPsiSpread->GetXaxis()->SetBinLabel(2, "<#Psi_{2, VZEROC} - #Psi_{2, TPC}>");
313 fHistPsiSpread->GetXaxis()->SetBinLabel(3, "<#Psi_{2, VZEROC} - #Psi_{2, TPC}>");
314 fHistPsiSpread->GetXaxis()->SetBinLabel(4, "<#Psi_{2, TPC, #eta < 0} - #Psi_{2, TPC, #eta > 0}>");
315 fOutputList->Add(fHistPsiControl);
316 fOutputList->Add(fHistPsiSpread);
317 fHistPsiVZEROA = BookTH1F("fHistPsiVZEROA", "#Psi_{VZEROA}", 100, -.5*TMath::Pi(), .5*TMath::Pi());
318 fHistPsiVZEROC = BookTH1F("fHistPsiVZEROC", "#Psi_{VZEROC}", 100, -.5*TMath::Pi(), .5*TMath::Pi());
319 fHistPsiTPC = BookTH1F("fHistPsiTPC", "#Psi_{TPC}", 100, -.5*TMath::Pi(), .5*TMath::Pi());
320 fHistPsiTPCSUBA = BookTH1F("fHistPsiTPCSUBA", "#Psi_{TPC, #eta < 0}", 100, -.5*TMath::Pi(), .5*TMath::Pi());
321 fHistPsiTPCSUBB = BookTH1F("fHistPsiTPCSUBB", "#Psi_{TPC, #eta > 0}", 100, -.5*TMath::Pi(), .5*TMath::Pi());
324 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i ++) {
325 fHistRhoPackage[i] = BookTH1F("fHistRhoPackage", "#rho [GeV/c]", 100, 0, 150, i);
326 fHistRho[i] = BookTH1F("fHistRho", "#rho [GeV/c]", 100, 0, 150, i);
328 fHistRhoVsMult = BookTH2F("fHistRhoVsMult", "multiplicity", "#rho [GeV/c]", 100, 0, 4000, 100, 0, 250);
329 fHistRhoVsCent = BookTH2F("fHistRhoVsCent", "centrality", "#rho [GeV/c]", 100, 0, 100, 100, 0, 250);
330 fHistRhoAVsMult = BookTH2F("fHistRhoAVsMult", "multiplicity", "#rho * A (jet) [GeV/c]", 100, 0, 4000, 100, 0, 50);
331 fHistRhoAVsCent = BookTH2F("fHistRhoAVsCent", "centrality", "#rho * A (jet) [GeV/c]", 100, 0, 100, 100, 0, 50);
333 // delta pt distributions
334 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i ++) {
335 fHistRCPhiEta[i] = BookTH2F("fHistRCPhiEta", "#phi (RC)", "#eta (RC)", 100, 0, TMath::TwoPi(), 100, -1, 1, i);
336 fHistRhoVsRCPt[i] = BookTH2F("fHistRhoVsRCPt", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i);
337 fHistRCPt[i] = BookTH1F("fHistRCPt", "p_{t} (RC) [GeV/c]", 130, -20, 150, i);
338 fHistRCPhiEtaExLJ[i] = BookTH2F("fHistRCPhiEtaExLJ", "#phi (RC)", "#eta (RC)", 100, 0, TMath::TwoPi(), 100, -1, 1, i);
339 fHistDeltaPtDeltaPhi2[i] = BookTH2F("fHistDeltaPtDeltaPhi2", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::Pi(), 100, -50, 100, i);
340 fHistDeltaPtDeltaPhi3[i] = BookTH2F("fHistDeltaPtDeltaPhi3", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::TwoPi()/3., 100, -50, 100, i);
341 fHistRhoVsRCPtExLJ[i] = BookTH2F("fHistRhoVsRCPtExLJ", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i);
342 fHistRCPtExLJ[i] = BookTH1F("fHistRCPtExLJ", "p_{t} (RC) [GeV/c]", 130, -20, 150, i);
343 /* fHistRCPhiEtaRand[i] = BookTH2F("fHistRCPhiEtaRand", "#phi (RC)", "#eta (RC)", 100, 0, TMath::TwoPi(), 100, -1, 1, i); */
344 fHistDeltaPtDeltaPhi2ExLJ[i] = BookTH2F("fHistDeltaPtDeltaPhi2ExLJ", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::Pi(), 100, -50, 100, i);
345 fHistDeltaPtDeltaPhi3ExLJ[i] = BookTH2F("fHistDeltaPtDeltaPhi3ExLJ", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::TwoPi()/3., 100, -50, 100, i);
346 /* fHistRhoVsRCPtRand[i] = BookTH2F("fHistRhoVsRCPtRand", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i); */
347 /* fHistRCPtRand[i] = BookTH1F("fHistRCPtRand", "p_{t} (RC) [GeV/c]", 130, -20, 150, i); */
348 /* fHistDeltaPtDeltaPhi2Rand[i] = BookTH2F("fHistDeltaPtDeltaPhi2Rand", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::Pi(), 100, -50, 100, i); */
349 /* fHistDeltaPtDeltaPhi3Rand[i] = BookTH2F("fHistDeltaPtDeltaPhi3Rand", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::TwoPi()/3., 100, -50, 100, i); */
350 // jet histograms (after kinematic cuts)
351 fHistJetPtRaw[i] = BookTH1F("fHistJetPtRaw", "p_{t} RAW [GeV/c]", 200, -50, 150, i);
352 fHistJetPt[i] = BookTH1F("fHistJetPt", "p_{t} [GeV/c]", 350, -100, 250, i);
353 fHistJetEtaPhi[i] = BookTH2F("fHistJetEtaPhi", "#eta", "#phi", 100, -1, 1, 100, 0, TMath::TwoPi(), i);
354 fHistJetPtArea[i] = BookTH2F("fHistJetPtArea", "p_{t} [GeV/c]", "Area", 175, -100, 250, 30, 0, 0.9, i);
355 fHistJetPtConstituents[i] = BookTH2F("fHistJetPtConstituents", "p_{t} [GeV/c]", "Area", 350, -100, 250, 60, 0, 150, i);
356 fHistJetEtaRho[i] = BookTH2F("fHistJetEtaRho", "#eta", "#rho", 100, -1, 1, 100, 0, 300, i);
357 // in plane and out of plane spectra
358 fHistJetPsiTPCPt[i] = BookTH2F("fHistJetPsiTPCPt", "#phi_{jet} - #Psi_{2, TPC}", "p_{t} [GeV/c]", 50, 0., TMath::Pi(), 700, -100, 250, i);
359 fHistJetPsiVZEROAPt[i] = BookTH2F("fHistJetPsiVZEROAPt", "#phi_{jet} - #Psi_{2, VZEROA}", "p_{t} [GeV/c]", 50, 0., TMath::Pi(), 700, -100, 250, i);
360 fHistJetPsiVZEROCPt[i] = BookTH2F("fHistJetPsiVZEROCPt", "#phi_{jet} - #Psi_{V2, ZEROC}", "p_{t} [GeV/c]", 50, 0., TMath::Pi(), 700, -100, 250, i);
362 fHistDeltaPhi2VZEROA[i] = BookTH1F("fHistDeltaPhi2VZEROA", "#phi_{jet} - #Psi_{2, VZEROA}", 50, 0, TMath::Pi(), i);
363 fHistDeltaPhi2VZEROC[i] = BookTH1F("fHistDeltaPhi2VZEROC", "#phi_{jet} - #Psi_{2, VZEROC}", 50, 0, TMath::Pi(), i);
364 fHistDeltaPhi2TPC[i] = BookTH1F("fHistDeltaPhi2TPC", "#phi_{jet} - #Psi_{2, TPC}", 50, 0, TMath::Pi(), i);
365 fHistDeltaPhi3VZEROA[i] = BookTH1F("fHistDeltaPhi3VZEROA", "#phi_{jet} - #Psi_{2, VZEROA}", 50, 0, TMath::TwoPi()/3., i);
366 fHistDeltaPhi3VZEROC[i] = BookTH1F("fHistDeltaPhi3VZEROC", "#phi_{jet} - #Psi_{2, VZEROC}", 50, 0, TMath::TwoPi()/3., i);
367 fHistDeltaPhi3TPC[i] = BookTH1F("fHistDeltaPhi3TPC", "#phi_{jet} - #Psi_{2, TPC}", 50, 0, TMath::TwoPi()/3., i);
369 fProfV2Resolution[i] = new TProfile(Form("fProfV2Resolution_%i", i), Form("fProfV2Resolution_%i", i), 8, -0.5, 7.5);
370 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(1, "<cos(2(#Psi_{a} - #Psi_{b}))>");
371 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(2, "<cos(2(#Psi_{b} - #Psi_{a}))>");
372 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(3, "<cos(2(#Psi_{VZEROA} - #Psi_{VZEROC}))>");
373 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(4, "<cos(2(#Psi_{VZEROC} - #Psi_{VZEROA}))>");
374 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(5, "<cos(2(#Psi_{VZEROA} - #Psi_{TPC}))>");
375 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(6, "<cos(2(#Psi_{TPC} - #Psi_{VZEROA}))>");
376 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(7, "<cos(2(#Psi_{VZEROC} - #Psi_{TPC}))>");
377 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(8, "<cos(2(#Psi_{TPC} - #Psi_{VZEROC}))>");
378 fOutputList->Add(fProfV2Resolution[i]);
379 fProfV3Resolution[i] = new TProfile(Form("fProfV3Resolution_%i", i), Form("fProfV3Resolution_%i", i), 8, -0.5, 7.5);
380 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(1, "<cos(3(#Psi_{a} - #Psi_{b}))>");
381 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(2, "<cos(3(#Psi_{b} - #Psi_{a}))>");
382 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(3, "<cos(3(#Psi_{VZEROA} - #Psi_{VZEROC}))>");
383 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(4, "<cos(3(#Psi_{VZEROC} - #Psi_{VZEROA}))>");
384 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(5, "<cos(3(#Psi_{VZEROA} - #Psi_{TPC}))>");
385 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(6, "<cos(3(#Psi_{TPC} - #Psi_{VZEROA}))>");
386 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(7, "<cos(3(#Psi_{VZEROC} - #Psi_{TPC}))>");
387 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(8, "<cos(3(#Psi_{TPC} - #Psi_{VZEROC}))>");
388 fOutputList->Add(fProfV3Resolution[i]);
390 // cdf and pdf of chisquare distribution
391 fHistPvaluePDF = BookTH1F("fHistPvaluePDF", "PDF #chi^{2}", 500, 0, 1);
392 fHistPvalueCDF = BookTH1F("fHistPvalueCDF", "CDF #chi^{2}", 500, 0, 1);
394 Float_t temp[fCentralityClasses->GetSize()];
395 for(Int_t i(0); i < fCentralityClasses->GetSize(); i++) temp[i] = fCentralityClasses->At(i);
396 fProfV2 = new TProfile("fProfV2", "fProfV2", fCentralityClasses->GetSize()-1, temp);
397 fProfV3 = new TProfile("fProfV3", "fProfV3", fCentralityClasses->GetSize()-1, temp);
398 fOutputList->Add(fProfV2);
399 fOutputList->Add(fProfV3);
401 // analysis summary histrogram, saves all relevant analysis settigns
402 fHistAnalysisSummary = BookTH1F("fHistAnalysisSummary", "flag", 44, -0.5, 44.5);
403 fHistAnalysisSummary->GetXaxis()->SetBinLabel(1, "fJetRadius");
404 fHistAnalysisSummary->SetBinContent(1, fJetRadius);
405 fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fPtBiasJetTrack");
406 fHistAnalysisSummary->SetBinContent(2, fPtBiasJetTrack);
407 fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fPtBiasJetClus");
408 fHistAnalysisSummary->SetBinContent(3, fPtBiasJetClus);
409 fHistAnalysisSummary->GetXaxis()->SetBinLabel(4, "fJetPtCut");
410 fHistAnalysisSummary->SetBinContent(4, fJetPtCut);
411 fHistAnalysisSummary->GetXaxis()->SetBinLabel(5, "fJetAreaCut");
412 fHistAnalysisSummary->SetBinContent(5, fJetAreaCut);
413 fHistAnalysisSummary->GetXaxis()->SetBinLabel(6, "fPercAreaCut");
414 fHistAnalysisSummary->SetBinContent(6, fPercAreaCut);
415 fHistAnalysisSummary->GetXaxis()->SetBinLabel(7, "fAreaEmcCut");
416 fHistAnalysisSummary->SetBinContent(7, fAreaEmcCut);
417 fHistAnalysisSummary->GetXaxis()->SetBinLabel(8, "fJetMinEta");
418 fHistAnalysisSummary->SetBinContent(8, fJetMinEta);
419 fHistAnalysisSummary->GetXaxis()->SetBinLabel(9, "fJetMaxEta");
420 fHistAnalysisSummary->SetBinContent(9, fJetMaxEta);
421 fHistAnalysisSummary->GetXaxis()->SetBinLabel(10, "fJetMinPhi");
422 fHistAnalysisSummary->SetBinContent(10, fJetMinPhi);
423 fHistAnalysisSummary->GetXaxis()->SetBinLabel(11, "fJetMaxPhi");
424 fHistAnalysisSummary->SetBinContent(11, fJetMaxPhi);
425 fHistAnalysisSummary->GetXaxis()->SetBinLabel(12, "fMaxClusterPt");
426 fHistAnalysisSummary->SetBinContent(12, fMaxClusterPt);
427 fHistAnalysisSummary->GetXaxis()->SetBinLabel(13, "fMaxTrackPt");
428 fHistAnalysisSummary->SetBinContent(13, fMaxTrackPt);
429 fHistAnalysisSummary->GetXaxis()->SetBinLabel(14, "fLeadingHadronType");
430 fHistAnalysisSummary->SetBinContent(14, fLeadingHadronType);
431 fHistAnalysisSummary->GetXaxis()->SetBinLabel(15, "fAnaType");
432 fHistAnalysisSummary->SetBinContent(15, fAnaType);
433 fHistAnalysisSummary->GetXaxis()->SetBinLabel(16, "fForceBeamType");
434 fHistAnalysisSummary->SetBinContent(16, fForceBeamType);
435 fHistAnalysisSummary->GetXaxis()->SetBinLabel(17, "fMinCent");
436 fHistAnalysisSummary->SetBinContent(17, fMinCent);
437 fHistAnalysisSummary->GetXaxis()->SetBinLabel(18, "fMaxCent");
438 fHistAnalysisSummary->SetBinContent(18, fMaxCent);
439 fHistAnalysisSummary->GetXaxis()->SetBinLabel(19, "fMinVz");
440 fHistAnalysisSummary->SetBinContent(19, fMinVz);
441 fHistAnalysisSummary->GetXaxis()->SetBinLabel(20, "fMaxVz");
442 fHistAnalysisSummary->SetBinContent(20, fMaxVz);
443 fHistAnalysisSummary->GetXaxis()->SetBinLabel(21, "fOffTrigger");
444 fHistAnalysisSummary->SetBinContent(21, fOffTrigger);
445 fHistAnalysisSummary->GetXaxis()->SetBinLabel(22, "fClusPtCut");
446 fHistAnalysisSummary->SetBinContent(22, fClusPtCut);
447 fHistAnalysisSummary->GetXaxis()->SetBinLabel(23, "fTrackPtCut");
448 fHistAnalysisSummary->SetBinContent(23, fTrackPtCut);
449 fHistAnalysisSummary->GetXaxis()->SetBinLabel(24, "fTrackMinEta");
450 fHistAnalysisSummary->SetBinContent(24, fTrackMinEta);
451 fHistAnalysisSummary->GetXaxis()->SetBinLabel(25, "fTrackMaxEta");
452 fHistAnalysisSummary->SetBinContent(25, fTrackMaxEta);
453 fHistAnalysisSummary->GetXaxis()->SetBinLabel(26, "fTrackMinPhi");
454 fHistAnalysisSummary->SetBinContent(26, fTrackMinPhi);
455 fHistAnalysisSummary->GetXaxis()->SetBinLabel(27, "fTrackMaxPhi");
456 fHistAnalysisSummary->SetBinContent(27, fTrackMaxPhi);
457 fHistAnalysisSummary->GetXaxis()->SetBinLabel(28, "fClusTimeCutLow");
458 fHistAnalysisSummary->SetBinContent(28, fClusTimeCutLow);
459 fHistAnalysisSummary->GetXaxis()->SetBinLabel(29, "fClusTimeCutUp");
460 fHistAnalysisSummary->SetBinContent(29, fClusTimeCutUp);
461 fHistAnalysisSummary->GetXaxis()->SetBinLabel(30, "fMinPtTrackInEmcal");
462 fHistAnalysisSummary->SetBinContent(30, fMinPtTrackInEmcal);
463 fHistAnalysisSummary->GetXaxis()->SetBinLabel(31, "fEventPlaneVsEmcal");
464 fHistAnalysisSummary->SetBinContent(31, fEventPlaneVsEmcal);
465 fHistAnalysisSummary->GetXaxis()->SetBinLabel(32, "fMinEventPlane");
466 fHistAnalysisSummary->SetBinContent(32, fMaxEventPlane);
467 fHistAnalysisSummary->GetXaxis()->SetBinLabel(33, "fRandomConeRadius");
468 fHistAnalysisSummary->SetBinContent(33, fRandomConeRadius);
469 fHistAnalysisSummary->GetXaxis()->SetBinLabel(34, "fitModulationType");
470 fHistAnalysisSummary->SetBinContent(34, (int)fFitModulationType);
471 fHistAnalysisSummary->GetXaxis()->SetBinLabel(35, "runModeType");
472 fHistAnalysisSummary->SetBinContent(35, (int)fRunModeType);
473 fHistAnalysisSummary->GetXaxis()->SetBinLabel(36, "data type");
474 fHistAnalysisSummary->SetBinContent(36, (int)fDataType);
475 fHistAnalysisSummary->GetXaxis()->SetBinLabel(37, "iterator");
476 fHistAnalysisSummary->SetBinContent(37, 1.);
477 fHistAnalysisSummary->GetXaxis()->SetBinLabel(38, "fMinPvalue");
478 fHistAnalysisSummary->SetBinContent(38, fMinPvalue);
479 fHistAnalysisSummary->GetXaxis()->SetBinLabel(39, "fMaxPvalue");
480 fHistAnalysisSummary->SetBinContent(39, fMaxPvalue);
481 fHistAnalysisSummary->GetXaxis()->SetBinLabel(40, "fExcludeLeadingJetsFromFit");
482 fHistAnalysisSummary->SetBinContent(40, fExcludeLeadingJetsFromFit);
483 fHistAnalysisSummary->GetXaxis()->SetBinLabel(41, "fRebinSwapHistoOnTheFly");
484 fHistAnalysisSummary->SetBinContent(41, (int)fRebinSwapHistoOnTheFly);
485 fHistAnalysisSummary->GetXaxis()->SetBinLabel(42, "fUsePtWeight");
486 fHistAnalysisSummary->SetBinContent(42, (int)fUsePtWeight);
487 fHistAnalysisSummary->GetXaxis()->SetBinLabel(43, "fMinLeadingHadronPt");
488 fHistAnalysisSummary->SetBinContent(43, fMinLeadingHadronPt);
489 fHistAnalysisSummary->GetXaxis()->SetBinLabel(44, "fExplicitOutlierCut");
490 fHistAnalysisSummary->SetBinContent(44, fExplicitOutlierCut);
492 if(fFillQAHistograms) {
493 fHistRunnumbersEta = new TH2F("fHistRunnumbersEta", "fHistRunnumbersEta", 100, -.5, 99.5, 100, -1.1, 1.1);
494 fHistRunnumbersEta->Sumw2();
495 fOutputList->Add(fHistRunnumbersEta);
496 fHistRunnumbersPhi = new TH2F("fHistRunnumbersPhi", "fHistRunnumbersPhi", 100, -.5, 99.5, 100, -0.2, TMath::TwoPi()+0.2);
497 fHistRunnumbersPhi->Sumw2();
498 fOutputList->Add(fHistRunnumbersPhi);
501 fHistSwap = new TH1F("fHistSwap", "fHistSwap", 20, 0, TMath::TwoPi());
504 if(fUserSuppliedV2) fOutputList->Add(fUserSuppliedV2);
505 if(fUserSuppliedV3) fOutputList->Add(fUserSuppliedV3);
506 // increase readability of output list
508 PostData(1, fOutputList);
510 switch (fRunModeType) {
512 fOutputListGood = new TList();
513 fOutputListGood->SetOwner(kTRUE);
514 fOutputListBad = new TList();
515 fOutputListBad->SetOwner(kTRUE);
516 PostData(2, fOutputListGood);
517 PostData(3, fOutputListBad);
522 //_____________________________________________________________________________
523 Bool_t AliAnalysisTaskRhoVnModulation::Run()
525 // user exec: execute once for each event
526 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
527 if(!fInitialized) fInitialized = InitializeAnalysis();
528 // reject the event if expected data is missing
529 if(!PassesCuts(InputEvent())) return kFALSE;
530 if(!(fTracks||fJets||fRho)) return kFALSE;
531 if(!fCaloClusters && fDebug > 0) printf(" > Warning: couldn't retreive calo clusters! < \n");
532 // [0][0] psi2a [1,0] psi2c
533 // [0][1] psi3a [1,1] psi3c
534 Double_t vzero[2][2];
535 CalculateEventPlaneVZERO(vzero);
537 // [2] psi2 a [3] psi2 b
538 // [4] psi3 a [5] psi3 b
540 CalculateEventPlaneTPC(tpc);
541 // arrays which will hold the fit parameters
542 Double_t fitParameters[] = {0,0,0,0,0,0,0,0,0};
543 Double_t psi2(-1), psi3(-1), psi2b(-1), psi3b(-1);
544 switch (fDetectorType) { // determine the detector type for the rho fit
545 case kTPC : { psi2 = tpc[0]; psi3 = tpc[1]; } break;
546 case kTPCSUB : { psi2 = tpc[2]; psi3 = tpc[4];
547 psi2b = tpc[3]; psi3b = tpc[5]; } break;
548 case kVZEROA : { psi2 = vzero[0][0]; psi3 = vzero[0][1]; } break;
549 case kVZEROC : { psi2 = vzero[1][0]; psi3 = vzero[1][1]; } break;
553 switch (fFitModulationType) { // do the fits
554 case kNoFit : { fFitModulation->FixParameter(0, RhoVal()); } break;
556 if(CorrectRho(fitParameters, psi2, psi3, psi2b, psi3b)) {
557 fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
558 CalculateEventPlaneResolution(vzero, tpc);
562 if(CorrectRho(fitParameters, psi2, psi3, psi2b, psi3b)) {
563 fProfV3->Fill(fCent, fFitModulation->GetParameter(3));
564 CalculateEventPlaneResolution(vzero, tpc);
568 CorrectRho(fitParameters, psi2, psi3, psi2b, psi3b);
571 if(CorrectRho(fitParameters, psi2, psi3, psi2b, psi3b)) {
572 fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
573 fProfV3->Fill(fCent, fFitModulation->GetParameter(7));
574 CalculateEventPlaneResolution(vzero, tpc);
578 // fill a number of histograms
579 FillHistogramsAfterSubtraction(vzero, tpc);
581 // send the output to the connected output container
582 PostData(1, fOutputList);
583 switch (fRunModeType) {
585 PostData(2, fOutputListGood);
586 PostData(3, fOutputListBad);
592 //_____________________________________________________________________________
593 void AliAnalysisTaskRhoVnModulation::CalculateEventPlaneVZERO(Double_t vzero[2][2]) const
595 // get the vzero event plane
596 if(fUseV0EventPlaneFromHeader) { // use the vzero from the header
597 Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0), h(0);
598 vzero[0][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 2, a, b);
599 vzero[1][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 2, c, d);
600 vzero[0][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 3, e, f);
601 vzero[1][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 3, g, h);
604 // grab the vzero event plane without recentering
605 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
606 Double_t qxa2(0), qya2(0), qxc2(0), qyc2(0); // for psi2
607 Double_t qxa3(0), qya3(0), qxc3(0), qyc3(0); // for psi3
608 for(Int_t iVZERO(0); iVZERO < 64; iVZERO++) {
609 Double_t phi(TMath::PiOver4()*(.5+iVZERO%8)), /* eta(0), */ weight(InputEvent()->GetVZEROEqMultiplicity(iVZERO));
610 // (iVZERO<32) ? eta = -3.45+.5*(iVZERO/8) : eta = 4.8-.6*((iVZERO/8)-4);
612 qxa2 += weight*TMath::Cos(2.*phi);
613 qya2 += weight*TMath::Sin(2.*phi);
614 qxa3 += weight*TMath::Cos(3.*phi);
615 qya3 += weight*TMath::Sin(3.*phi);
618 qxc2 += weight*TMath::Cos(2.*phi);
619 qyc2 += weight*TMath::Sin(2.*phi);
620 qxc3 += weight*TMath::Cos(3.*phi);
621 qyc3 += weight*TMath::Sin(3.*phi);
624 vzero[0][0] = .5*TMath::ATan2(qya2, qxa2);
625 vzero[1][0] = .5*TMath::ATan2(qyc2, qxc2);
626 vzero[0][1] = (1./3.)*TMath::ATan2(qya3, qxa3);
627 vzero[1][1] = (1./3.)*TMath::ATan2(qyc3, qxc3);
629 //_____________________________________________________________________________
630 void AliAnalysisTaskRhoVnModulation::CalculateEventPlaneTPC(Double_t* tpc)
632 // grab the TPC event plane
633 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
634 fNAcceptedTracks = 0; // reset the track counter
635 Double_t qx2(0), qy2(0); // for psi2
636 Double_t qx3(0), qy3(0); // for psi3
637 Double_t qx2a(0), qy2a(0), qx2b(0), qy2b(0); // for psi2 a and b
638 Double_t qx3a(0), qy3a(0), qx3b(0), qy3b(0); // for psi3 a and b
640 Float_t excludeInEta[] = {-999, -999};
641 if(fExcludeLeadingJetsFromFit > 0 ) { // remove the leading jet from ep estimate
642 AliEmcalJet* leadingJet[] = {0x0, 0x0};
643 static Int_t lJets[9999] = {-1};
644 GetSortedArray(lJets, fJets);
645 for(Int_t i(0); i < fJets->GetEntriesFast(); i++) { // get the two leading jets
646 if (1 + i > fJets->GetEntriesFast()) break;
647 leadingJet[0] = static_cast<AliEmcalJet*>(fJets->At(lJets[i]));
648 leadingJet[1] = static_cast<AliEmcalJet*>(fJets->At(lJets[i+1]));
649 if(PassesCuts(leadingJet[0]) && PassesCuts(leadingJet[1])) break;
651 if(leadingJet[0] && leadingJet[1]) {
652 for(Int_t i(0); i < 2; i++) excludeInEta[i] = leadingJet[i]->Eta();
655 Int_t iTracks(fTracks->GetEntriesFast());
656 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
657 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
658 if(!PassesCuts(track) || track->Pt() < .15 || track->Pt() > 5.) continue;
659 if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta[0]) < fJetRadius*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - fJetRadius - fJetMaxEta ) > 0 )) continue;
661 qx2+= TMath::Cos(2.*track->Phi());
662 qy2+= TMath::Sin(2.*track->Phi());
663 qx3+= TMath::Cos(3.*track->Phi());
664 qy3+= TMath::Sin(3.*track->Phi());
665 if(track->Eta() < 0) { // A side, negative eta
666 qx2a+= TMath::Cos(2.*track->Phi());
667 qy2a+= TMath::Sin(2.*track->Phi());
668 qx3a+= TMath::Cos(3.*track->Phi());
669 qy3a+= TMath::Sin(3.*track->Phi());
670 } else { // B side, positive eta
671 qx2b+= TMath::Cos(2.*track->Phi());
672 qy2b+= TMath::Sin(2.*track->Phi());
673 qx3b+= TMath::Cos(3.*track->Phi());
674 qy3b+= TMath::Sin(3.*track->Phi());
678 tpc[0] = .5*TMath::ATan2(qy2, qx2);
679 tpc[1] = (1./3.)*TMath::ATan2(qy3, qx3);
680 tpc[2] = .5*TMath::ATan2(qy2a, qx2a);
681 tpc[3] = .5*TMath::ATan2(qy2b, qx2b);
682 tpc[4] = (1./3.)*TMath::ATan2(qy3a, qx3a);
683 tpc[5] = (1./3.)*TMath::ATan2(qy3b, qx3b);
685 //_____________________________________________________________________________
686 void AliAnalysisTaskRhoVnModulation::CalculateEventPlaneResolution(Double_t vzero[2][2], Double_t* tpc) const
688 // fill the profiles for the resolution parameters
689 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
690 fProfV2Resolution[fInCentralitySelection]->Fill(0., TMath::Cos(2.*(tpc[2] - tpc[3])));
691 fProfV2Resolution[fInCentralitySelection]->Fill(1., TMath::Cos(2.*(tpc[3] - tpc[2])));
692 fProfV2Resolution[fInCentralitySelection]->Fill(2., TMath::Cos(2.*(vzero[0][0] - vzero[1][0])));
693 fProfV2Resolution[fInCentralitySelection]->Fill(3., TMath::Cos(2.*(vzero[1][0] - vzero[0][0])));
694 fProfV2Resolution[fInCentralitySelection]->Fill(4., TMath::Cos(2.*(vzero[0][0] - tpc[0])));
695 fProfV2Resolution[fInCentralitySelection]->Fill(5., TMath::Cos(2.*(tpc[0] - vzero[0][0])));
696 fProfV2Resolution[fInCentralitySelection]->Fill(6., TMath::Cos(2.*(vzero[1][0] - tpc[0])));
697 fProfV2Resolution[fInCentralitySelection]->Fill(7., TMath::Cos(2.*(tpc[0] - vzero[1][0])));
698 fProfV3Resolution[fInCentralitySelection]->Fill(0., TMath::Cos(3.*(tpc[2] - tpc[3])));
699 fProfV3Resolution[fInCentralitySelection]->Fill(1., TMath::Cos(3.*(tpc[3] - tpc[2])));
700 fProfV3Resolution[fInCentralitySelection]->Fill(2., TMath::Cos(3.*(vzero[0][0] - vzero[1][0])));
701 fProfV3Resolution[fInCentralitySelection]->Fill(3., TMath::Cos(3.*(vzero[1][0] - vzero[0][0])));
702 fProfV3Resolution[fInCentralitySelection]->Fill(4., TMath::Cos(3.*(vzero[0][0] - tpc[0])));
703 fProfV3Resolution[fInCentralitySelection]->Fill(5., TMath::Cos(3.*(tpc[0] - vzero[0][0])));
704 fProfV3Resolution[fInCentralitySelection]->Fill(6., TMath::Cos(3.*(vzero[1][0] - tpc[0])));
705 fProfV3Resolution[fInCentralitySelection]->Fill(7., TMath::Cos(3.*(tpc[0] - vzero[1][0])));
707 //_____________________________________________________________________________
708 void AliAnalysisTaskRhoVnModulation::CalculateRandomCone(Float_t &pt, Float_t &eta, Float_t &phi,
709 AliEmcalJet* jet, Bool_t randomize) const
712 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
713 pt = 0; eta = 0; phi = 0;
714 Float_t etaJet(999), phiJet(999), dJet(999); // no jet: same as jet very far away
715 if(jet) { // if a leading jet is given, use its kinematic properties
719 // force the random cones to at least be within detector acceptance
720 Float_t minPhi(fJetMinPhi), maxPhi(fJetMaxPhi);
721 if(maxPhi > TMath::TwoPi()) maxPhi = TMath::TwoPi();
722 if(minPhi < 0 ) minPhi = 0;
723 Float_t diffRcRJR(TMath::Abs(fRandomConeRadius-fJetRadius));
724 // construct a random cone and see if it's far away enough from the leading jet
725 Int_t attempts(1000);
728 eta = gRandom->Uniform(fJetMinEta+diffRcRJR, fJetMaxEta-diffRcRJR);
729 phi = gRandom->Uniform(minPhi, maxPhi);
731 dJet = TMath::Sqrt((etaJet-eta)*(etaJet-eta)+(phiJet-phi)*(phiJet-phi));
732 if(dJet > fMinDisanceRCtoLJ) break;
733 else if (attempts == 0) {
734 printf(" > No random cone after 1000 tries, giving up ... !\n");
739 Int_t iTracks(fTracks->GetEntriesFast());
740 for(Int_t i(0); i < iTracks; i++) {
741 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
742 if(!PassesCuts(track)) continue;
743 Float_t etaTrack(track->Eta()), phiTrack(track->Phi()), ptTrack(track->Pt());
744 // if requested, randomize eta and phi to destroy any correlated fluctuations
746 etaTrack = gRandom->Uniform(fTrackMinEta, fTrackMaxEta);
747 phiTrack = gRandom->Uniform(minPhi, maxPhi);
749 // get distance from cone
750 if(TMath::Abs(phiTrack-phi) > TMath::Abs(phiTrack - phi + TMath::TwoPi())) phiTrack+=TMath::TwoPi();
751 if(TMath::Abs(phiTrack-phi) > TMath::Abs(phiTrack - phi - TMath::TwoPi())) phiTrack-=TMath::TwoPi();
752 if(TMath::Sqrt(TMath::Abs((etaTrack-eta)*(etaTrack-eta)+(phiTrack-phi)*(phiTrack-phi))) <= fRandomConeRadius) pt+=ptTrack;
756 //_____________________________________________________________________________
757 Bool_t AliAnalysisTaskRhoVnModulation::CorrectRho(Double_t* params, Double_t psi2, Double_t psi3, Double_t psi2b, Double_t psi3b)
759 // get rho' -> rho(phi)
760 // two routines are available
761 // [1] fitting a fourier expansion to the de/dphi distribution
762 // [2] getting vn from a fourier series around dn/dphi (see below for info)
763 // this function will return kTRUE if the fit passes a set of quality criteria
764 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
765 TString detector("");
766 switch (fDetectorType) {
767 case kTPC : detector+="TPC";
769 case kTPCSUB : detector+="kTPCSUB";
771 case kVZEROA : detector+="VZEROA";
773 case kVZEROC : detector+="VZEROC";
777 Int_t iTracks(fTracks->GetEntriesFast());
778 Double_t excludeInEta[] = {-999, -999};
779 Double_t excludeInPhi[] = {-999, -999};
780 Double_t excludeInPt[] = {-999, -999};
781 if(iTracks <= 0 || RhoVal() <= 0 ) return kFALSE; // no use fitting an empty event ...
782 if(fExcludeLeadingJetsFromFit > 0 ) {
783 AliEmcalJet* leadingJet[] = {0x0, 0x0};
784 static Int_t lJets[9999] = {-1};
785 GetSortedArray(lJets, fJets);
786 for(Int_t i(0); i < fJets->GetEntriesFast(); i++) { // get the two leading jets
787 if (1 + i > fJets->GetEntriesFast()) break;
788 leadingJet[0] = static_cast<AliEmcalJet*>(fJets->At(lJets[i]));
789 leadingJet[1] = static_cast<AliEmcalJet*>(fJets->At(lJets[i+1]));
790 if(PassesCuts(leadingJet[0]) && PassesCuts(leadingJet[1])) break;
792 if(leadingJet[0] && leadingJet[1]) {
793 for(Int_t i(0); i < 2; i++) {
794 excludeInEta[i] = leadingJet[i]->Eta();
795 excludeInPhi[i] = leadingJet[i]->Phi();
796 excludeInPt[i] = leadingJet[i]->Pt();
800 fHistSwap->Reset(); // clear the histogram
802 if(fRebinSwapHistoOnTheFly) {
803 if(fNAcceptedTracks < 49) fNAcceptedTracks = 49; // avoid aliasing effects
804 _tempSwap = TH1F("_tempSwap", "_tempSwap", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), 0, TMath::TwoPi());
806 else _tempSwap = *fHistSwap; // now _tempSwap holds the desired histo
807 for(Int_t i(0); i < iTracks; i++) {
808 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
809 if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta[0]) < fJetRadius*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - fJetRadius - fJetMaxEta ) > 0 )) continue;
810 if(!PassesCuts(track) || track->Pt() > 5 || track->Pt() < 0.15) continue;
811 if(fDetectorType == kTPCSUB && psi2 > -1000 && track->Eta() < 0 ) continue;
812 else if (fDetectorType == kTPCSUB && psi2 < -1000 && track->Eta() > 0 ) continue;
813 if(fUsePtWeight) _tempSwap.Fill(track->Phi(), track->Pt());
814 else _tempSwap.Fill(track->Phi());
816 for(Int_t i(0); i < _tempSwap.GetXaxis()->GetNbins(); i++) _tempSwap.SetBinError(1+i, TMath::Sqrt(_tempSwap.GetBinContent(1+i)));
817 fFitModulation->SetParameter(0, RhoVal());
818 switch (fFitModulationType) {
819 case kNoFit : { fFitModulation->FixParameter(0, RhoVal() );
822 fFitModulation->FixParameter(4, psi2);
825 fFitModulation->FixParameter(4, psi3);
828 fFitModulation->FixParameter(4, psi2);
829 fFitModulation->FixParameter(6, psi3);
831 case kFourierSeries : {
832 // in this approach, an explicit calculation will be made of vn = sqrt(xn^2+yn^2)
833 // where x[y] = Integrate[r(phi)cos[sin](n phi)dphi, 0, 2pi]
834 Double_t cos2(0), sin2(0), cos3(0), sin3(0), sumPt(0);
835 for(Int_t i(0); i < iTracks; i++) {
836 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
837 if(!PassesCuts(track) || track->Pt() > 5 || track->Pt() < 0.15) continue;
838 sumPt += track->Pt();
839 cos2 += track->Pt()*TMath::Cos(2*PhaseShift(track->Phi()-psi2));
840 sin2 += track->Pt()*TMath::Sin(2*PhaseShift(track->Phi()-psi2));
841 cos3 += track->Pt()*TMath::Cos(3*PhaseShift(track->Phi()-psi3));
842 sin3 += track->Pt()*TMath::Sin(3*PhaseShift(track->Phi()-psi3));
844 fFitModulation->SetParameter(3, TMath::Sqrt(cos2*cos2+sin2*sin2)/RhoVal());
845 fFitModulation->SetParameter(4, psi2);
846 fFitModulation->SetParameter(6, psi3);
847 fFitModulation->SetParameter(7, TMath::Sqrt(cos3*cos3+sin3*sin3)/RhoVal());
849 case kIntegratedFlow : {
850 // use v2 and v3 values from an earlier iteration over the data
851 fFitModulation->FixParameter(3, fUserSuppliedV2->GetBinContent(fUserSuppliedV2->GetXaxis()->FindBin(fCent)));
852 fFitModulation->FixParameter(4, psi2);
853 fFitModulation->FixParameter(6, psi3);
854 fFitModulation->FixParameter(7, fUserSuppliedV3->GetBinContent(fUserSuppliedV3->GetXaxis()->FindBin(fCent)));
855 return kTRUE; // no fit is performed
859 if(fDetectorType == kTPCSUB && psi2 > -1000 ) { // do the magic for the subevent case
860 Double_t v2(fFitModulation->GetParameter(3)), v3(fFitModulation->GetParameter(7));
861 CorrectRho(params, -9999, -9999, psi2b, psi3b);
862 v2 += fFitModulation->GetParameter(3);
863 v3 += fFitModulation->GetParameter(7);
864 fFitModulation->SetParameter(3, v2/2.);
865 fFitModulation->SetParameter(7, v3/3.);
867 _tempSwap.Fit(fFitModulation, fFitModulationOptions.Data(), "", 0, TMath::TwoPi());
868 // the quality of the fit is evaluated from 1 - the cdf of the chi square distribution
869 Double_t CDF(1.-ChiSquareCDF(fFitModulation->GetNDF(), fFitModulation->GetChisquare()));
870 // Double_t PDF(ChiSquarePDF(fFitModulation->GetNDF(), fFitModulation->GetChisquare()));
871 fHistPvalueCDF->Fill(CDF);
872 // fHistPvaluePDF->Fill(PDF);
873 if(CDF > fMinPvalue && CDF < fMaxPvalue && ( fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) > 0)) { // fit quality
874 // for LOCAL didactic purposes, save the best and the worst fits
875 // this routine can produce a lot of output histograms (it's not memory 'safe') and will not work on GRID
876 // since the output will become unmergeable (i.e. different nodes may produce conflicting output)
877 switch (fRunModeType) {
879 if(fRandom->Uniform(0, 100) > fPercentageOfFits) break;
880 static Int_t didacticCounterBest(0);
881 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
882 TF1* didactifFit = (TF1*)fFitModulation->Clone(Form("fit_%i_CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
883 didacticProfile->GetListOfFunctions()->Add(didactifFit);
884 fOutputListGood->Add(didacticProfile);
885 didacticCounterBest++;
886 TH2F* didacticSurface = BookTH2F(Form("surface_%s", didacticProfile->GetName()), "#phi", "#eta", 50, 0, TMath::TwoPi(), 50, -1, 1, -1, kFALSE);
887 for(Int_t i(0); i < iTracks; i++) {
888 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
889 if(PassesCuts(track)) {
890 if(fUsePtWeight) didacticSurface->Fill(track->Phi(), track->Eta(), track->Pt());
891 else didacticSurface->Fill(track->Phi(), track->Eta());
894 if(fExcludeLeadingJetsFromFit) { // visualize the excluded region
895 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);
896 f2->SetParameters(excludeInPt[0]/3.,excludeInPhi[0],.1,excludeInEta[0],.1);
897 didacticSurface->GetListOfFunctions()->Add(f2);
898 TF2 *f3 = new TF2(Form("%s_NLJ", didacticSurface->GetName()),"[0]*TMath::Gaus(x,[1],[2])*TMath::Gaus(y,[3],[4])", 0, TMath::TwoPi(), -1, 1);
899 f3->SetParameters(excludeInPt[1]/3.,excludeInPhi[1],.1,excludeInEta[1],.1);
900 f3->SetLineColor(kGreen);
901 didacticSurface->GetListOfFunctions()->Add(f3);
903 fOutputListGood->Add(didacticSurface);
907 } else { // if the fit is of poor quality revert to the original rho estimate
908 switch (fRunModeType) { // again see if we want to save the fit
910 static Int_t didacticCounterWorst(0);
911 if(fRandom->Uniform(0, 100) > fPercentageOfFits) break;
912 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data() ));
913 TF1* didactifFit = (TF1*)fFitModulation->Clone(Form("fit_%i_p_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data()));
914 didacticProfile->GetListOfFunctions()->Add(didactifFit);
915 fOutputListBad->Add(didacticProfile);
916 didacticCounterWorst++;
920 switch (fFitModulationType) {
921 case kNoFit : break; // nothing to do
922 case kUser : break; // FIXME not implemented yet
923 case kCombined : fFitModulation->SetParameter(7, 0); // no break
924 case kFourierSeries : fFitModulation->SetParameter(7, 0); // no break
925 default : { // needs to be done if there was a poor fit
926 fFitModulation->SetParameter(3, 0);
927 fFitModulation->SetParameter(0, RhoVal());
930 return kFALSE; // return false if the fit is rejected
932 for(Int_t i(0); i < fFitModulation->GetNpar(); i++) params[i] = fFitModulation->GetParameter(i);
935 //_____________________________________________________________________________
936 Bool_t AliAnalysisTaskRhoVnModulation::PassesCuts(AliVEvent* event)
939 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
940 if(!event) return kFALSE;
941 if(TMath::Abs(InputEvent()->GetPrimaryVertex()->GetZ()) > 10.) return kFALSE;
942 // aod and esd specific checks
945 AliESDEvent* esdEvent = static_cast<AliESDEvent*>(InputEvent());
946 if( (!esdEvent) || (TMath::Abs(esdEvent->GetPrimaryVertexSPD()->GetZ() - esdEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE;
949 AliAODEvent* aodEvent = static_cast<AliAODEvent*>(InputEvent());
950 if( (!aodEvent) || (TMath::Abs(aodEvent->GetPrimaryVertexSPD()->GetZ() - aodEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE;
954 fCent = InputEvent()->GetCentrality()->GetCentralityPercentile("V0M");
955 if(fCent <= fCentralityClasses->At(0) || fCent >= fCentralityClasses->At(fCentralityClasses->GetSize()-1) || TMath::Abs(fCent-InputEvent()->GetCentrality()->GetCentralityPercentile("TRK")) > 5.) return kFALSE;
956 // determine centrality class
957 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
958 if(fCent >= fCentralityClasses->At(i) && fCent <= fCentralityClasses->At(1+i)) {
959 fInCentralitySelection = i;
962 if(fExplicitOutlierCut == 2010 || fExplicitOutlierCut == 2011) {
963 if(!PassesCuts(fExplicitOutlierCut)) return kFALSE;
965 if(fFillQAHistograms) FillQAHistograms(event);
968 //_____________________________________________________________________________
969 Bool_t AliAnalysisTaskRhoVnModulation::PassesCuts(Int_t year)
971 // additional centrality cut based on relation between tpc and global multiplicity
972 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
973 AliAODEvent* event(dynamic_cast<AliAODEvent*>(InputEvent()));
974 if(!event) return kFALSE;
975 Int_t multTPC(0), multGlob(0), nTracks(InputEvent()->GetNumberOfTracks());
976 for(Int_t iTracks = 0; iTracks < nTracks; iTracks++) {
977 AliAODTrack* track = event->GetTrack(iTracks);
979 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
980 if (track->TestFilterBit(1) && track->Chi2perNDF() > 0.2) multTPC++;
981 if (!track->TestFilterBit(16) || track->Chi2perNDF() < 0.1) continue;
982 Double_t b[2] = {-99., -99.};
983 Double_t bCov[3] = {-99., -99., -99.};
984 if (track->PropagateToDCA(event->GetPrimaryVertex(), event->GetMagneticField(), 100., b, bCov) && TMath::Abs(b[0]) < 0.3 && TMath::Abs(b[1]) < 0.3) multGlob++;
986 if(year == 2010 && multTPC > (-40.3+1.22*multGlob) && multTPC < (32.1+1.59*multGlob)) return kTRUE;
987 if(year == 2011 && multTPC > (-36.73 + 1.48*multGlob) && multTPC < (62.87 + 1.78*multGlob)) return kTRUE;
990 //_____________________________________________________________________________
991 Bool_t AliAnalysisTaskRhoVnModulation::PassesCuts(const AliVCluster* cluster) const
994 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
995 if(!cluster) return kFALSE;
998 //_____________________________________________________________________________
999 void AliAnalysisTaskRhoVnModulation::FillHistogramsAfterSubtraction(Double_t vzero[2][2], Double_t* tpc) const
1002 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1003 FillTrackHistograms();
1004 /* FillClusterHistograms(); */
1005 FillJetHistograms(vzero, tpc);
1006 /* FillCorrectedClusterHistograms(); */
1007 FillEventPlaneHistograms(vzero, tpc);
1008 FillRhoHistograms();
1009 switch (fDetectorType) { // determine the detector type for the rho fit
1010 case kTPC : { FillDeltaPtHistograms(tpc[0], tpc[1]); } break;
1011 case kTPCSUB : { FillDeltaPtHistograms(tpc[2], tpc[4]);
1012 FillDeltaPtHistograms(tpc[3], tpc[5]); } break;
1013 case kVZEROA : { FillDeltaPtHistograms(vzero[0][0], vzero[0][1]); } break;
1014 case kVZEROC : { FillDeltaPtHistograms(vzero[1][0], vzero[1][1]); } break;
1017 FillDeltaPhiHistograms(vzero, tpc);
1019 //_____________________________________________________________________________
1020 void AliAnalysisTaskRhoVnModulation::FillTrackHistograms() const
1022 // fill track histograms
1023 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1024 Int_t iTracks(fTracks->GetEntriesFast());
1025 for(Int_t i(0); i < iTracks; i++) {
1026 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1027 if(!PassesCuts(track)) continue;
1028 fHistPicoTrackPt[fInCentralitySelection]->Fill(track->Pt());
1029 if(fFillQAHistograms) FillQAHistograms(track);
1032 //_____________________________________________________________________________
1033 void AliAnalysisTaskRhoVnModulation::FillClusterHistograms() const
1035 // fill cluster histograms
1036 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1037 /* Int_t iClusters(fCaloClusters->GetEntriesFast());
1038 for(Int_t i(0); i < iClusters; i++) {
1039 AliVCluster* cluster = static_cast<AliVCluster*>(fCaloClusters->At(iClusters));
1040 if (!PassesCuts(cluster)) continue;
1041 TLorentzVector clusterLorentzVector;
1042 cluster->GetMomentum(clusterLorentzVector, const_cast<Double_t*>(fVertex));
1043 fHistClusterPt[fInCentralitySelection]->Fill(clusterLorentzVector.Pt());
1044 fHistClusterEta[fInCentralitySelection]->Fill(clusterLorentzVector.Eta());
1045 fHistClusterPhi[fInCentralitySelection]->Fill(clusterLorentzVector.Phi());
1049 //_____________________________________________________________________________
1050 void AliAnalysisTaskRhoVnModulation::FillCorrectedClusterHistograms() const
1052 // fill clusters after hadronic correction FIXME implement
1053 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1055 //_____________________________________________________________________________
1056 void AliAnalysisTaskRhoVnModulation::FillEventPlaneHistograms(Double_t vzero[2][2], Double_t* tpc) const
1058 // fill event plane histograms
1059 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1060 fHistPsiControl->Fill(0.5, vzero[0][0]); // vzero a psi2
1061 fHistPsiControl->Fill(1.5, vzero[1][0]); // vzero c psi2
1062 fHistPsiControl->Fill(2.5, tpc[0]); // tpc psi 2
1063 fHistPsiControl->Fill(3.5, tpc[2]); // tpc sub a psi 2
1064 fHistPsiControl->Fill(4.5, tpc[3]); // tpc sub b psi 2
1065 fHistPsiControl->Fill(5.5, vzero[0][1]); // vzero a psi3
1066 fHistPsiControl->Fill(6.5, vzero[1][1]); // vzero b psi3
1067 fHistPsiControl->Fill(7.5, tpc[1]); // tpc psi 3
1068 fHistPsiControl->Fill(8.5, tpc[4]); // tpc sub a psi3
1069 fHistPsiControl->Fill(9.5, tpc[5]); // tpc sub b psi3
1070 fHistPsiVZEROA->Fill(vzero[0][0]);
1071 fHistPsiVZEROC->Fill(vzero[1][0]);
1072 fHistPsiTPC->Fill(tpc[0]);
1073 fHistPsiTPCSUBA->Fill(tpc[2]);
1074 fHistPsiTPCSUBB->Fill(tpc[3]);
1075 fHistPsiSpread->Fill(0.5, TMath::Abs(vzero[0][0]-vzero[1][0]));
1076 fHistPsiSpread->Fill(1.5, TMath::Abs(vzero[0][0]-tpc[0]));
1077 fHistPsiSpread->Fill(2.5, TMath::Abs(vzero[1][0]-tpc[0]));
1078 fHistPsiSpread->Fill(3.5, TMath::Abs(tpc[2]-tpc[3]));
1080 //_____________________________________________________________________________
1081 void AliAnalysisTaskRhoVnModulation::FillRhoHistograms() const
1083 // fill rho histograms
1084 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1085 fHistRhoPackage[fInCentralitySelection]->Fill(RhoVal()); // save the rho estimate from the emcal jet package
1086 // get multiplicity FIXME inefficient
1087 Int_t iTracks(fTracks->GetEntriesFast()), mult(0), iJets(fJets->GetEntriesFast());
1088 for(Int_t i(0); i < iTracks; i ++) { if(PassesCuts(static_cast<AliVTrack*>(fTracks->At(i)))) mult++; }
1089 Double_t rho(RhoVal(TMath::Pi(), TMath::Pi(), fRho->GetVal()));
1090 fHistRho[fInCentralitySelection]->Fill(rho);
1091 fHistRhoVsMult->Fill(mult, rho);
1092 fHistRhoVsCent->Fill(fCent, rho);
1093 for(Int_t i(0); i < iJets; i++) {
1094 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
1095 if(!PassesCuts(jet)) continue;
1096 fHistRhoAVsMult->Fill(mult, rho * jet->Area());
1097 fHistRhoAVsCent->Fill(fCent, rho * jet->Area());
1100 //_____________________________________________________________________________
1101 void AliAnalysisTaskRhoVnModulation::FillDeltaPtHistograms(Double_t psi2, Double_t psi3) const
1103 // fill delta pt histograms
1104 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1105 Int_t i(0), maxCones(20);
1106 AliEmcalJet* leadingJet(0x0);
1107 static Int_t sJets[9999] = {-1};
1108 GetSortedArray(sJets, fJets);
1109 do { // get the leading jet
1110 leadingJet = static_cast<AliEmcalJet*>(fJets->At(sJets[i]));
1113 while (!PassesCuts(leadingJet)&&i<fJets->GetEntriesFast());
1114 if(!leadingJet && fDebug > 0) printf(" > failed to retrieve leading jet ! < \n");
1115 const Float_t areaRC = fRandomConeRadius*fRandomConeRadius*TMath::Pi();
1116 // we're retrieved the leading jet, now get a random cone
1117 for(i = 0; i < maxCones; i++) {
1118 Float_t pt(0), eta(0), phi(0);
1119 // get a random cone without constraints on leading jet position
1120 CalculateRandomCone(pt, eta, phi, 0x0);
1122 fHistRCPhiEta[fInCentralitySelection]->Fill(phi, eta);
1123 fHistRhoVsRCPt[fInCentralitySelection]->Fill(pt, RhoVal(phi, fJetRadius, fRho->GetVal())*areaRC);
1124 fHistRCPt[fInCentralitySelection]->Fill(pt);
1125 fHistDeltaPtDeltaPhi2[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1126 fHistDeltaPtDeltaPhi3[fInCentralitySelection]->Fill(PhaseShift(phi-psi3, 3.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1128 // get a random cone excluding leading jet area
1129 CalculateRandomCone(pt, eta, phi, leadingJet);
1131 fHistRCPhiEtaExLJ[fInCentralitySelection]->Fill(phi, eta);
1132 fHistRhoVsRCPtExLJ[fInCentralitySelection]->Fill(pt, RhoVal(phi, fJetRadius, fRho->GetVal())*areaRC);
1133 fHistRCPtExLJ[fInCentralitySelection]->Fill(pt);
1134 fHistDeltaPtDeltaPhi2ExLJ[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1135 fHistDeltaPtDeltaPhi3ExLJ[fInCentralitySelection]->Fill(PhaseShift(phi-psi3, 3.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1137 // get a random cone in an event with randomized phi and eta
1138 /* CalculateRandomCone(pt, eta, phi, 0x0, kTRUE);
1140 fHistRCPhiEtaRand[fInCentralitySelection]->Fill(phi, eta);
1141 fHistRhoVsRCPtRand[fInCentralitySelection]->Fill(pt, RhoVal(phi, fJetRadius, fRho->GetVal())*areaRC);
1142 fHistRCPtRand[fInCentralitySelection]->Fill(pt);
1143 fHistDeltaPtDeltaPhi2Rand[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1144 fHistDeltaPtDeltaPhi3Rand[fInCentralitySelection]->Fill(PhaseShift(phi-psi3, 3.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1148 //_____________________________________________________________________________
1149 void AliAnalysisTaskRhoVnModulation::FillJetHistograms(Double_t vzero[2][2], Double_t* tpc) const
1151 // fill jet histograms
1152 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1153 Int_t iJets(fJets->GetEntriesFast());
1154 for(Int_t i(0); i < iJets; i++) {
1155 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
1156 if(PassesCuts(jet)) {
1157 Double_t pt(jet->Pt()), area(jet->Area()), eta(jet->Eta()), phi(jet->Phi());
1158 Double_t rho(RhoVal(phi, fJetRadius, fRho->GetVal()));
1159 fHistJetPtRaw[fInCentralitySelection]->Fill(pt);
1160 fHistJetPt[fInCentralitySelection]->Fill(pt-area*rho);
1161 fHistJetEtaPhi[fInCentralitySelection]->Fill(eta, phi);
1162 fHistJetPtArea[fInCentralitySelection]->Fill(pt-area*rho, area);
1163 fHistJetPsiTPCPt[fInCentralitySelection]->Fill(PhaseShift(phi-tpc[0], 2.), pt-area*rho);
1164 fHistJetPsiVZEROAPt[fInCentralitySelection]->Fill(PhaseShift(phi-vzero[0][0], 2.), pt-area*rho);
1165 fHistJetPsiVZEROCPt[fInCentralitySelection]->Fill(PhaseShift(phi-vzero[1][0], 2.), pt-area*rho);
1166 fHistJetPtConstituents[fInCentralitySelection]->Fill(pt-area*rho, jet->Nch());
1167 fHistJetEtaRho[fInCentralitySelection]->Fill(eta, pt/area);
1168 if(fSetPtSub) jet->SetPtSub(pt-area*rho);
1170 else { // if the jet is rejected, excluded it for the flow analysis
1171 if(fSetPtSub) jet->SetPtSub(-999.);
1175 //_____________________________________________________________________________
1176 void AliAnalysisTaskRhoVnModulation::FillDeltaPhiHistograms(Double_t vzero[2][2], Double_t* tpc) const
1178 // fill phi minus psi histograms
1179 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1181 Int_t iTracks(fTracks->GetEntriesFast());
1182 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
1183 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
1184 if(!PassesCuts(track)) continue;
1185 fHistDeltaPhi2VZEROA[fInCentralitySelection]->Fill(PhaseShift(track->Phi()-vzero[0][0], 2.));
1186 fHistDeltaPhi2VZEROC[fInCentralitySelection]->Fill(PhaseShift(track->Phi()-vzero[1][0], 2.));
1187 fHistDeltaPhi2TPC[fInCentralitySelection]->Fill(PhaseShift(track->Phi()-tpc[0], 2.));
1188 fHistDeltaPhi3VZEROA[fInCentralitySelection]->Fill(PhaseShift(track->Phi()-vzero[0][1], 3.));
1189 fHistDeltaPhi3VZEROC[fInCentralitySelection]->Fill(PhaseShift(track->Phi()-vzero[1][1], 3.));
1190 fHistDeltaPhi3TPC[fInCentralitySelection]->Fill(PhaseShift(track->Phi()-tpc[1], 3.));
1194 //_____________________________________________________________________________
1195 void AliAnalysisTaskRhoVnModulation::FillQAHistograms(AliVTrack* vtrack) const
1197 // fill qa histograms for pico tracks
1199 AliPicoTrack* track = static_cast<AliPicoTrack*>(vtrack);
1200 fHistRunnumbersPhi->Fill(fMappedRunNumber, track->Phi());
1201 fHistRunnumbersEta->Fill(fMappedRunNumber, track->Eta());
1202 Int_t type((int)(track->GetTrackType()));
1205 fHistPicoCat1[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1208 fHistPicoCat2[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1211 fHistPicoCat3[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1216 //_____________________________________________________________________________
1217 void AliAnalysisTaskRhoVnModulation::FillQAHistograms(AliVEvent* vevent)
1219 // fill qa histograms for events
1221 fHistVertexz->Fill(vevent->GetPrimaryVertex()->GetZ());
1222 fHistCentrality->Fill(fCent);
1223 Int_t runNumber(InputEvent()->GetRunNumber());
1224 Int_t runs[] = {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, 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};
1225 for(fMappedRunNumber = 0; fMappedRunNumber < 64; fMappedRunNumber++) {
1226 if(runs[fMappedRunNumber]==runNumber) break;
1229 //_____________________________________________________________________________
1230 void AliAnalysisTaskRhoVnModulation::Terminate(Option_t *)
1233 switch (fRunModeType) {
1235 printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1236 if(fFillQAHistograms) {
1237 Int_t runs[] = {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, 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};
1238 for(Int_t i(0); i < 64; i++) {
1239 fHistRunnumbersPhi->GetXaxis()->SetBinLabel(i+1, Form("%i", runs[i]));
1240 fHistRunnumbersEta->GetXaxis()->SetBinLabel(i+1, Form("%i", runs[i]));
1242 fHistRunnumbersPhi->GetXaxis()->SetBinLabel(65, "undetermined");
1243 fHistRunnumbersEta->GetXaxis()->SetBinLabel(65, "undetermined");
1245 AliAnalysisTaskRhoVnModulation::Dump();
1246 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));
1251 //_____________________________________________________________________________