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
43 #include <AliAnalysisTask.h>
44 #include <AliAnalysisManager.h>
45 #include <AliCentrality.h>
46 #include <AliVVertex.h>
47 #include <AliESDEvent.h>
48 #include <AliAODEvent.h>
49 #include <AliAODTrack.h>
50 // emcal jet framework includes
51 #include <AliPicoTrack.h>
52 #include <AliEmcalJet.h>
53 #include <AliRhoParameter.h>
54 #include <AliLocalRhoParameter.h>
55 #include <AliAnalysisTaskRhoVnModulation.h>
57 class AliAnalysisTaskRhoVnModulation;
60 ClassImp(AliAnalysisTaskRhoVnModulation)
62 AliAnalysisTaskRhoVnModulation::AliAnalysisTaskRhoVnModulation() : AliAnalysisTaskEmcalJet("AliAnalysisTaskRhoVnModulation", kTRUE),
63 fDebug(0), fRunToyMC(kFALSE), fLocalInit(0), fAttachToEvent(kTRUE), fSemiCentralInclusive(kFALSE), fFillHistograms(kTRUE), fFillQAHistograms(kTRUE), fReduceBinsXByFactor(-1.), fReduceBinsYByFactor(-1.), fNoEventWeightsForQC(kTRUE), fCentralityClasses(0), fPtBinsHybrids(0), fPtBinsJets(0), fUserSuppliedV2(0), fUserSuppliedV3(0), fUserSuppliedR2(0), fUserSuppliedR3(0), fTracksCont(0), fJetsCont(0), fLeadingJet(0), fUseScaledRho(0), fNAcceptedTracks(0), fNAcceptedTracksQCn(0), fFitModulationType(kNoFit), fFitGoodnessTest(kChi2Poisson), fQCRecovery(kTryFit), fUsePtWeight(kTRUE), fUsePtWeightErrorPropagation(kTRUE), fDetectorType(kTPC), 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), fNameJetClones(0), fNamePicoTrackClones(0), fNameRho(0), fLocalJetMinEta(-10), fLocalJetMaxEta(-10), fLocalJetMinPhi(-10), fLocalJetMaxPhi(-10), 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), fMinDisanceRCtoLJ(0), fRandomConeRadius(-1.), fMaxCones(-1), fAbsVnHarmonics(kTRUE), fExcludeLeadingJetsFromFit(1.), fRebinSwapHistoOnTheFly(kTRUE), fPercentageOfFits(10.), fUseV0EventPlaneFromHeader(kTRUE), fExplicitOutlierCut(-1), fMinLeadingHadronPt(0), fSubtractJetPt(kFALSE), 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) {
64 for(Int_t i(0); i < 10; i++) {
65 fProfV2Resolution[i] = 0;
66 fProfV3Resolution[i] = 0;
67 fHistPicoTrackPt[i] = 0;
68 fHistPicoTrackMult[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 fHistJetPsi2Pt[i] = 0;
102 fHistJetPsi3Pt[i] = 0;
104 // default constructor
106 //_____________________________________________________________________________
107 AliAnalysisTaskRhoVnModulation::AliAnalysisTaskRhoVnModulation(const char* name, runModeType type) : AliAnalysisTaskEmcalJet(name, kTRUE),
108 fDebug(0), fRunToyMC(kFALSE), fLocalInit(0), fAttachToEvent(kTRUE), fSemiCentralInclusive(kFALSE), fFillHistograms(kTRUE), fFillQAHistograms(kTRUE), fReduceBinsXByFactor(-1.), fReduceBinsYByFactor(-1.), fNoEventWeightsForQC(kTRUE), fCentralityClasses(0), fPtBinsHybrids(0), fPtBinsJets(0), fUserSuppliedV2(0), fUserSuppliedV3(0), fUserSuppliedR2(0), fUserSuppliedR3(0), fTracksCont(0), fJetsCont(0), fLeadingJet(0), fUseScaledRho(0), fNAcceptedTracks(0), fNAcceptedTracksQCn(0), fFitModulationType(kNoFit), fFitGoodnessTest(kChi2Poisson), fQCRecovery(kTryFit), fUsePtWeight(kTRUE), fUsePtWeightErrorPropagation(kTRUE), fDetectorType(kTPC), 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), fNameJetClones(0), fNamePicoTrackClones(0), fNameRho(0), fLocalJetMinEta(-10), fLocalJetMaxEta(-10), fLocalJetMinPhi(-10), fLocalJetMaxPhi(-10), 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), fMinDisanceRCtoLJ(0), fRandomConeRadius(-1.), fMaxCones(-1), fAbsVnHarmonics(kTRUE), fExcludeLeadingJetsFromFit(1.), fRebinSwapHistoOnTheFly(kTRUE), fPercentageOfFits(10.), fUseV0EventPlaneFromHeader(kTRUE), fExplicitOutlierCut(-1), fMinLeadingHadronPt(0), fSubtractJetPt(kFALSE), 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) {
109 for(Int_t i(0); i < 10; i++) {
110 fProfV2Resolution[i] = 0;
111 fProfV3Resolution[i] = 0;
112 fHistPicoTrackPt[i] = 0;
113 fHistPicoTrackMult[i] = 0;
114 fHistPicoCat1[i] = 0;
115 fHistPicoCat2[i] = 0;
116 fHistPicoCat3[i] = 0;
117 /* fHistClusterPt[i] = 0; */
118 /* fHistClusterPhi[i] = 0; */
119 /* fHistClusterEta[i] = 0; */
120 /* fHistClusterCorrPt[i] = 0; */
121 /* fHistClusterCorrPhi[i] = 0; */
122 /* fHistClusterCorrEta[i] = 0; */
123 fHistRhoPackage[i] = 0;
125 fHistRCPhiEta[i] = 0;
126 fHistRhoVsRCPt[i] = 0;
128 fHistDeltaPtDeltaPhi2[i] = 0;
129 fHistDeltaPtDeltaPhi3[i] = 0;
130 fHistRCPhiEtaExLJ[i] = 0;
131 fHistRhoVsRCPtExLJ[i] = 0;
132 fHistRCPtExLJ[i] = 0;
133 fHistDeltaPtDeltaPhi2ExLJ[i] = 0;
134 fHistDeltaPtDeltaPhi3ExLJ[i] = 0;
135 /* fHistRCPhiEtaRand[i] = 0; */
136 /* fHistRhoVsRCPtRand[i] = 0; */
137 /* fHistRCPtRand[i] = 0; */
138 /* fHistDeltaPtDeltaPhi2Rand[i] = 0; */
139 /* fHistDeltaPtDeltaPhi3Rand[i] = 0; */
140 fHistJetPtRaw[i] = 0;
142 fHistJetEtaPhi[i] = 0;
143 fHistJetPtArea[i] = 0;
144 fHistJetPtConstituents[i] = 0;
145 fHistJetEtaRho[i] = 0;
146 fHistJetPsi2Pt[i] = 0;
147 fHistJetPsi3Pt[i] = 0;
150 DefineInput(0, TChain::Class());
151 DefineOutput(1, TList::Class());
152 switch (fRunModeType) {
154 gStyle->SetOptFit(1);
155 DefineOutput(2, TList::Class());
156 DefineOutput(3, TList::Class());
158 default: fDebug = -1; // suppress debug info explicitely when not running locally
160 switch (fCollisionType) {
162 fFitModulationType = kNoFit;
166 if(fLocalRhoName=="") fLocalRhoName = Form("LocalRhoFrom_%s", GetName());
168 //_____________________________________________________________________________
169 AliAnalysisTaskRhoVnModulation::~AliAnalysisTaskRhoVnModulation()
172 if(fOutputList) delete fOutputList;
173 if(fOutputListGood) delete fOutputListGood;
174 if(fOutputListBad) delete fOutputListBad;
175 if(fFitModulation) delete fFitModulation;
176 if(fHistSwap) delete fHistSwap;
177 if(fCentralityClasses) delete fCentralityClasses;
178 if(fFitControl) delete fFitControl;
180 //_____________________________________________________________________________
181 void AliAnalysisTaskRhoVnModulation::ExecOnce()
184 fLocalRho = new AliLocalRhoParameter(fLocalRhoName.Data(), 0);
186 if(!(InputEvent()->FindListObject(fLocalRho->GetName()))) {
187 InputEvent()->AddObject(fLocalRho);
189 AliFatal(Form("%s: Container with name %s already present. Aborting", GetName(), fLocalRho->GetName()));
192 AliAnalysisTaskEmcalJet::ExecOnce(); // init the base class
193 AliAnalysisTaskEmcalJet::SetVzRange(-1.*fAbsVertexZ, fAbsVertexZ);
195 // unscaled rho has been retrieved by the parent class, now we retrieve rho scaled
196 fRho = dynamic_cast<AliRhoParameter*>(InputEvent()->FindListObject(Form("%s_Scaled", fRho->GetName())));
198 AliFatal(Form("%s: Couldn't find container for scaled rho. Aborting !", GetName()));
201 if(!GetJetContainer()) AliFatal(Form("%s: Couldn't find jet container. Aborting !", GetName()));
203 //_____________________________________________________________________________
204 Bool_t AliAnalysisTaskRhoVnModulation::InitializeAnalysis()
206 // initialize the anaysis
207 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
208 if(fRandomConeRadius <= 0) fRandomConeRadius = GetJetContainer()->GetJetRadius();
209 if(fMaxCones <= 0) fMaxCones = TMath::Nint(1.8*TMath::TwoPi()/(TMath::Pi()*fRandomConeRadius*fRandomConeRadius));
210 if(fLocalJetMinEta > -10 && fLocalJetMaxEta > -10) GetJetContainer()->SetJetEtaLimits(fLocalJetMinEta, fLocalJetMaxEta);
211 if(fLocalJetMinPhi > -10 && fLocalJetMaxPhi > -10) GetJetContainer()->SetJetPhiLimits(fLocalJetMinPhi, fLocalJetMaxPhi);
212 if(fMinDisanceRCtoLJ==0) fMinDisanceRCtoLJ = .5*GetJetRadius();
213 if(dynamic_cast<AliAODEvent*>(InputEvent())) fDataType = kAOD; // determine the datatype
214 else if(dynamic_cast<AliESDEvent*>(InputEvent())) fDataType = kESD;
215 fHistAnalysisSummary->SetBinContent(36, (int)fDataType);
216 if(!fRandom) fRandom = new TRandom3(0); // get a randomized if one hasn't been user-supplied
217 switch (fFitModulationType) {
218 case kNoFit : { SetModulationFit(new TF1("fix_kNoFit", "[0]", 0, TMath::TwoPi())); } break;
220 SetModulationFit(new TF1("fit_kV2", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
221 fFitModulation->SetParameter(0, 0.); // normalization
222 fFitModulation->SetParameter(3, 0.2); // v2
223 fFitModulation->FixParameter(1, 1.); // constant
224 fFitModulation->FixParameter(2, 2.); // constant
227 SetModulationFit(new TF1("fit_kV3", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
228 fFitModulation->SetParameter(0, 0.); // normalization
229 fFitModulation->SetParameter(3, 0.2); // v3
230 fFitModulation->FixParameter(1, 1.); // constant
231 fFitModulation->FixParameter(2, 3.); // constant
233 default : { // for the combined fit, the 'direct fourier series' or the user supplied vn values we use v2 and v3
234 SetModulationFit(new TF1("fit_kCombined", "[0]*([1]+[2]*([3]*TMath::Cos([2]*(x-[4]))+[7]*TMath::Cos([5]*(x-[6]))))", 0, TMath::TwoPi()));
235 fFitModulation->SetParameter(0, 0.); // normalization
236 fFitModulation->SetParameter(3, 0.2); // v2
237 fFitModulation->FixParameter(1, 1.); // constant
238 fFitModulation->FixParameter(2, 2.); // constant
239 fFitModulation->FixParameter(5, 3.); // constant
240 fFitModulation->SetParameter(7, 0.2); // v3
243 switch (fRunModeType) {
244 case kGrid : { fFitModulationOptions += "N0"; } break;
247 FillAnalysisSummaryHistogram();
250 //_____________________________________________________________________________
251 TH1F* AliAnalysisTaskRhoVnModulation::BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Int_t c, Bool_t append)
253 // book a TH1F 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 ) bins = TMath::Nint(bins/fReduceBinsXByFactor);
256 if(!fOutputList) return 0x0;
258 if(c!=-1) { // format centrality dependent histograms accordingly
259 name = Form("%s_%i", name, c);
260 title += Form("_%i-%i", (int)(fCentralityClasses->At(c)), (int)(fCentralityClasses->At((1+c))));
262 title += Form(";%s;[counts]", x);
263 TH1F* histogram = new TH1F(name, title.Data(), bins, min, max);
265 if(append) fOutputList->Add(histogram);
268 //_____________________________________________________________________________
269 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)
271 // book a TH2F and connect it to the output container
272 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
273 if(fReduceBinsXByFactor > 0 ) binsx = TMath::Nint(binsx/fReduceBinsXByFactor);
274 if(fReduceBinsYByFactor > 0 ) binsy = TMath::Nint(binsy/fReduceBinsYByFactor);
275 if(!fOutputList) return 0x0;
277 if(c!=-1) { // format centrality dependent histograms accordingly
278 name = Form("%s_%i", name, c);
279 title += Form("_%i-%i", (int)fCentralityClasses->At(c), (int)(fCentralityClasses->At((1+c))));
281 title += Form(";%s;%s", x, y);
282 TH2F* histogram = new TH2F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy);
284 if(append) fOutputList->Add(histogram);
287 //_____________________________________________________________________________
288 void AliAnalysisTaskRhoVnModulation::UserCreateOutputObjects()
290 // create output objects
291 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
292 fOutputList = new TList();
293 fOutputList->SetOwner(kTRUE);
294 if(!fCentralityClasses) { // classes must be defined at this point
295 Double_t c[] = {0., 20., 40., 60., 80., 100.};
296 fCentralityClasses = new TArrayD(sizeof(c)/sizeof(c[0]), c);
299 fHistCentrality = BookTH1F("fHistCentrality", "centrality", 102, -2, 100);
300 fHistVertexz = BookTH1F("fHistVertexz", "vertex z (cm)", 100, -12, 12);
302 // pico track kinematics
303 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
304 fHistPicoTrackPt[i] = BookTH1F("fHistPicoTrackPt", "p_{t} [GeV/c]", 100, 0, 100, i);
305 fHistPicoTrackMult[i] = BookTH1F("fHistPicoTrackMult", "multiplicity", 100, 0, 5000, i);
306 if(fFillQAHistograms) {
307 fHistPicoCat1[i] = BookTH2F("fHistPicoCat1", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
308 fHistPicoCat2[i] = BookTH2F("fHistPicoCat2", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
309 fHistPicoCat3[i] = BookTH2F("fHistPicoCat3", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
312 /* fHistClusterPt[i] = BookTH1F("fHistClusterPt", "p_{t} [GeV/c]", 100, 0, 100, i); */
313 /* fHistClusterPhi[i] = BookTH1F("fHistClusterPhi", "#phi", 100, 0, TMath::TwoPi(), i); */
314 /* fHistClusterEta[i] = BookTH1F("fHistClusterEta", "#eta", 100, -5, 5); */
316 // emcal kinematics after hadronic correction
317 /* fHistClusterCorrPt[i] = BookTH1F("fHistClusterCorrPt", "p_{t} [GeV/c]", 100, 0, 100, i); */
318 /* fHistClusterCorrPhi[i] = BookTH1F("fHistClusterCorrPhi", "#phi", 100, 0, TMath::TwoPi(), i); */
319 /* fHistClusterCorrEta[i] = BookTH1F("fHistClusterCorrEta", "#eta", 100, -5, 5, i); */
322 if(fFillQAHistograms) {
323 // event plane estimates and quality
324 fHistPsiControl = new TProfile("fHistPsiControl", "fHistPsiControl", 10, 0, 10);
325 fHistPsiControl->Sumw2();
326 fHistPsiSpread = new TProfile("fHistPsiSpread", "fHistPsiSpread", 4, 0, 4);
327 fHistPsiSpread->Sumw2();
328 fHistPsiControl->GetXaxis()->SetBinLabel(1, "<#Psi_{2, VZEROA}>");
329 fHistPsiControl->GetXaxis()->SetBinLabel(2, "<#Psi_{2, VZEROC}>");
330 fHistPsiControl->GetXaxis()->SetBinLabel(3, "<#Psi_{2, TPC}>");
331 fHistPsiControl->GetXaxis()->SetBinLabel(4, "<#Psi_{2, TPC, #eta < 0}>");
332 fHistPsiControl->GetXaxis()->SetBinLabel(5, "<#Psi_{2, TPC, #eta > 0}>");
333 fHistPsiControl->GetXaxis()->SetBinLabel(6, "<#Psi_{3, VZEROA}>");
334 fHistPsiControl->GetXaxis()->SetBinLabel(7, "<#Psi_{3, VZEROC}>");
335 fHistPsiControl->GetXaxis()->SetBinLabel(8, "<#Psi_{3, TPC}>");
336 fHistPsiControl->GetXaxis()->SetBinLabel(9, "<#Psi_{3, TPC, #eta < 0}>");
337 fHistPsiControl->GetXaxis()->SetBinLabel(10, "<#Psi_{3, TPC, #eta > 0}>");
338 fHistPsiSpread->GetXaxis()->SetBinLabel(1, "<#Psi_{2, VZEROA} - #Psi_{2, VZEROC}>");
339 fHistPsiSpread->GetXaxis()->SetBinLabel(2, "<#Psi_{2, VZEROC} - #Psi_{2, TPC}>");
340 fHistPsiSpread->GetXaxis()->SetBinLabel(3, "<#Psi_{2, VZEROC} - #Psi_{2, TPC}>");
341 fHistPsiSpread->GetXaxis()->SetBinLabel(4, "<#Psi_{2, TPC, #eta < 0} - #Psi_{2, TPC, #eta > 0}>");
342 fOutputList->Add(fHistPsiControl);
343 fOutputList->Add(fHistPsiSpread);
344 fHistPsiVZEROA = BookTH1F("fHistPsiVZEROA", "#Psi_{VZEROA}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
345 fHistPsiVZEROC = BookTH1F("fHistPsiVZEROC", "#Psi_{VZEROC}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
346 fHistPsiVZERO = BookTH1F("fHistPsiVZERO", "#Psi_{VZERO}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
347 fHistPsiTPC = BookTH1F("fHistPsiTPC", "#Psi_{TPC}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
348 fHistPsiVZEROAV0M = BookTH2F("fHistPsiVZEROAV0M", "V0M", "#Psi_{2, VZEROA}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
349 fHistPsiVZEROCV0M = BookTH2F("fHistPsiVZEROCV0M", "V0M", "#Psi_{2, VZEROC}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
350 fHistPsiVZEROVV0M = BookTH2F("fHistPsiVZEROV0M", "V0M", "#Psi_{2, VZERO}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
351 fHistPsiTPCiV0M = BookTH2F("fHistPsiTPCV0M", "V0M", "#Psi_{2, TRK}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
352 fHistPsiVZEROATRK = BookTH2F("fHistPsiVZEROATRK", "TRK", "#Psi_{2, VZEROA}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
353 fHistPsiVZEROCTRK = BookTH2F("fHistPsiVZEROCTRK", "TRK", "#Psi_{2, VZEROC}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
354 fHistPsiVZEROTRK = BookTH2F("fHistPsiVZEROTRK", "TRK", "#Psi_{2, VZERO}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
355 fHistPsiTPCTRK = BookTH2F("fHistPsiTPCTRK", "TRK", "#Psi_{2, TRK}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
358 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i ++) {
359 fHistRhoPackage[i] = BookTH1F("fHistRhoPackage", "#rho [GeV/c]", 100, 0, 150, i);
360 fHistRho[i] = BookTH1F("fHistRho", "#rho [GeV/c]", 100, 0, 150, i);
362 fHistRhoVsMult = BookTH2F("fHistRhoVsMult", "multiplicity", "#rho [GeV/c]", 100, 0, 4000, 100, 0, 250);
363 fHistRhoVsCent = BookTH2F("fHistRhoVsCent", "centrality", "#rho [GeV/c]", 100, 0, 100, 100, 0, 250);
364 fHistRhoAVsMult = BookTH2F("fHistRhoAVsMult", "multiplicity", "#rho * A (jet) [GeV/c]", 100, 0, 4000, 100, 0, 50);
365 fHistRhoAVsCent = BookTH2F("fHistRhoAVsCent", "centrality", "#rho * A (jet) [GeV/c]", 100, 0, 100, 100, 0, 50);
367 TString detector("");
368 switch (fDetectorType) {
369 case kTPC : detector+="TPC";
371 case kVZEROA : detector+="VZEROA";
373 case kVZEROC : detector+="VZEROC";
375 case kVZEROComb : detector+="VZEROComb";
379 // delta pt distributions
380 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i ++) {
381 if(fFillQAHistograms) fHistRCPhiEta[i] = BookTH2F("fHistRCPhiEta", "#phi (RC)", "#eta (RC)", 40, 0, TMath::TwoPi(), 40, -1, 1, i);
382 fHistRhoVsRCPt[i] = BookTH2F("fHistRhoVsRCPt", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i);
383 fHistRCPt[i] = BookTH1F("fHistRCPt", "p_{t} (RC) [GeV/c]", 130, -20, 150, i);
384 if(fFillQAHistograms) fHistRCPhiEtaExLJ[i] = BookTH2F("fHistRCPhiEtaExLJ", "#phi (RC)", "#eta (RC)", 40, 0, TMath::TwoPi(), 40, -1, 1, i);
385 fHistDeltaPtDeltaPhi2[i] = BookTH2F("fHistDeltaPtDeltaPhi2", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
386 fHistDeltaPtDeltaPhi3[i] = BookTH2F("fHistDeltaPtDeltaPhi3", Form("#phi - #Psi_{3, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::TwoPi()/3., 400, -70, 130, i);
387 fHistRhoVsRCPtExLJ[i] = BookTH2F("fHistRhoVsRCPtExLJ", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i);
388 fHistRCPtExLJ[i] = BookTH1F("fHistRCPtExLJ", "p_{t} (RC) [GeV/c]", 130, -20, 150, i);
389 /* fHistRCPhiEtaRand[i] = BookTH2F("fHistRCPhiEtaRand", "#phi (RC)", "#eta (RC)", 100, 0, TMath::TwoPi(), 100, -1, 1, i); */
390 fHistDeltaPtDeltaPhi2ExLJ[i] = BookTH2F("fHistDeltaPtDeltaPhi2ExLJ", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
391 fHistDeltaPtDeltaPhi3ExLJ[i] = BookTH2F("fHistDeltaPtDeltaPhi3ExLJ", Form("#phi - #Psi_{3, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::TwoPi()/3., 400, -70, 130, i);
392 /* fHistRhoVsRCPtRand[i] = BookTH2F("fHistRhoVsRCPtRand", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i); */
393 /* fHistRCPtRand[i] = BookTH1F("fHistRCPtRand", "p_{t} (RC) [GeV/c]", 130, -20, 150, i); */
394 /* fHistDeltaPtDeltaPhi2Rand[i] = BookTH2F("fHistDeltaPtDeltaPhi2Rand", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::Pi(), 100, -50, 100, i); */
395 /* fHistDeltaPtDeltaPhi3Rand[i] = BookTH2F("fHistDeltaPtDeltaPhi3Rand", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::TwoPi()/3., 100, -50, 100, i); */
396 // jet histograms (after kinematic cuts)
397 fHistJetPtRaw[i] = BookTH1F("fHistJetPtRaw", "p_{t, jet} RAW [GeV/c]", 200, -50, 150, i);
398 fHistJetPt[i] = BookTH1F("fHistJetPt", "p_{t, jet} [GeV/c]", 350, -100, 250, i);
399 if(fFillQAHistograms) fHistJetEtaPhi[i] = BookTH2F("fHistJetEtaPhi", "#eta", "#phi", 100, -1, 1, 100, 0, TMath::TwoPi(), i);
400 fHistJetPtArea[i] = BookTH2F("fHistJetPtArea", "p_{t, jet} [GeV/c]", "Area", 175, -100, 250, 30, 0, 0.9, i);
401 fHistJetPtConstituents[i] = BookTH2F("fHistJetPtConstituents", "p_{t, jet} [GeV/c]", "Area", 350, -100, 250, 60, 0, 150, i);
402 fHistJetEtaRho[i] = BookTH2F("fHistJetEtaRho", "#eta", "#rho", 100, -1, 1, 100, 0, 300, i);
403 // in plane and out of plane spectra
404 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);
405 fHistJetPsi3Pt[i] = BookTH2F("fHistJetPsi3Pt", Form("#phi_{jet} - #Psi_{3, %s}", detector.Data()), "p_{t, jet} [GeV/c]", 40, 0., TMath::TwoPi()/3., 350, -100, 250, i);
406 // profiles for all correlator permutations which are necessary to calculate each second and third order event plane resolution
407 fProfV2Resolution[i] = new TProfile(Form("fProfV2Resolution_%i", i), Form("fProfV2Resolution_%i", i), 11, -0.5, 10.5);
408 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(3, "<cos(2(#Psi_{VZEROA} - #Psi_{VZEROC}))>");
409 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(4, "<cos(2(#Psi_{VZEROC} - #Psi_{VZEROA}))>");
410 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(5, "<cos(2(#Psi_{VZEROA} - #Psi_{TPC}))>");
411 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(6, "<cos(2(#Psi_{TPC} - #Psi_{VZEROA}))>");
412 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(7, "<cos(2(#Psi_{VZEROC} - #Psi_{TPC}))>");
413 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(8, "<cos(2(#Psi_{TPC} - #Psi_{VZEROC}))>");
414 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(9, "<cos(2(#Psi_{VZERO} - #Psi_{TPC_A}))>");
415 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(10, "<cos(2(#Psi_{VZERO} - #Psi_{TPC_B}))>");
416 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(11, "<cos(2(#Psi_{TPC_A} - #Psi_{TPC_B}))>");
417 fOutputList->Add(fProfV2Resolution[i]);
418 fProfV3Resolution[i] = new TProfile(Form("fProfV3Resolution_%i", i), Form("fProfV3Resolution_%i", i), 11, -0.5, 10.5);
419 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(3, "<cos(3(#Psi_{VZEROA} - #Psi_{VZEROC}))>");
420 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(4, "<cos(3(#Psi_{VZEROC} - #Psi_{VZEROA}))>");
421 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(5, "<cos(3(#Psi_{VZEROA} - #Psi_{TPC}))>");
422 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(6, "<cos(3(#Psi_{TPC} - #Psi_{VZEROA}))>");
423 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(7, "<cos(3(#Psi_{VZEROC} - #Psi_{TPC}))>");
424 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(8, "<cos(3(#Psi_{TPC} - #Psi_{VZEROC}))>");
425 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(9, "<cos(3(#Psi_{VZERO} - #Psi_{TPC_A}))>");
426 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(10, "<cos(3(#Psi_{VZERO} - #Psi_{TPC_B}))>");
427 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(11, "<cos(3(#Psi_{TPC_A} - #Psi_{TPC_B}))>");
428 fOutputList->Add(fProfV3Resolution[i]);
431 Float_t temp[fCentralityClasses->GetSize()];
432 for(Int_t i(0); i < fCentralityClasses->GetSize(); i++) temp[i] = fCentralityClasses->At(i);
433 fProfV2 = new TProfile("fProfV2", "fProfV2", fCentralityClasses->GetSize()-1, temp);
434 fProfV3 = new TProfile("fProfV3", "fProfV3", fCentralityClasses->GetSize()-1, temp);
435 fOutputList->Add(fProfV2);
436 fOutputList->Add(fProfV3);
437 switch (fFitModulationType) {
439 fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
440 fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
441 fOutputList->Add(fProfV2Cumulant);
442 fOutputList->Add(fProfV3Cumulant);
445 fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
446 fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
447 fOutputList->Add(fProfV2Cumulant);
448 fOutputList->Add(fProfV3Cumulant);
452 // for the histograms initialized below, binning is fixed to runnumbers or flags
453 fReduceBinsXByFactor = 1;
454 fReduceBinsYByFactor = 1;
455 if(fFillQAHistograms) {
456 fHistRunnumbersEta = new TH2F("fHistRunnumbersEta", "fHistRunnumbersEta", 100, -.5, 99.5, 100, -1.1, 1.1);
457 fHistRunnumbersEta->Sumw2();
458 fOutputList->Add(fHistRunnumbersEta);
459 fHistRunnumbersPhi = new TH2F("fHistRunnumbersPhi", "fHistRunnumbersPhi", 100, -.5, 99.5, 100, -0.2, TMath::TwoPi()+0.2);
460 fHistRunnumbersPhi->Sumw2();
461 fOutputList->Add(fHistRunnumbersPhi);
463 fHistAnalysisSummary = BookTH1F("fHistAnalysisSummary", "flag", 52, -0.5, 52.5);
464 fHistSwap = new TH1F("fHistSwap", "fHistSwap", 20, 0, TMath::TwoPi());
465 if(fUsePtWeight) fHistSwap->Sumw2();
467 if(fUserSuppliedV2) fOutputList->Add(fUserSuppliedV2);
468 if(fUserSuppliedV3) fOutputList->Add(fUserSuppliedV3);
469 if(fUserSuppliedR2) fOutputList->Add(fUserSuppliedR2);
470 if(fUserSuppliedR3) fOutputList->Add(fUserSuppliedR3);
471 // increase readability of output list
473 // cdf and pdf of chisquare distribution
474 fHistPvalueCDF = BookTH1F("fHistPvalueCDF", "CDF #chi^{2}", 50, 0, 1);
475 fHistPvalueCDFCent = BookTH2F("fHistPvalueCDFCent", "centrality", "p-value", 40, 0, 100, 40, 0, 1);
476 fHistChi2Cent = BookTH2F("fHistChi2Cent", "centrality", "#tilde{#chi^{2}}", 100, 0, 100, 100, 0, 5);
477 fHistPChi2 = BookTH2F("fHistPChi2", "p-value", "#tilde{#chi^{2}}", 1000, 0, 1, 100, 0, 5);
478 fHistKolmogorovTest = BookTH1F("fHistKolmogorovTest", "KolmogorovTest", 50, 0, 1);
479 fHistKolmogorovTestCent = BookTH2F("fHistKolmogorovTestCent", "centrality", "Kolmogorov p", 40, 0, 100, 45, 0, 1);
480 fHistPvalueCDFROOT = BookTH1F("fHistPvalueCDFROOT", "CDF #chi^{2} ROOT", 50, 0, 1);
481 fHistPvalueCDFROOTCent = BookTH2F("fHistPvalueCDFROOTCent", "centrality", "p-value ROOT", 40, 0, 100, 45, 0, 1);
482 fHistChi2ROOTCent = BookTH2F("fHistChi2ROOTCent", "centrality", "#tilde{#chi^{2}}", 40, 0, 100, 45, 0, 5);
483 fHistPChi2Root = BookTH2F("fHistPChi2Root", "p-value", "#tilde{#chi^{2}} ROOT", 1000, 0, 1, 100, 0, 5);
484 fHistPKolmogorov = BookTH2F("fHistPKolmogorov", "p-value", "kolmogorov p",40, 0, 1, 40, 0, 1);
485 fHistRhoStatusCent = BookTH2F("fHistRhoStatusCent", "centrality", "status [-1=lin was better, 0=ok, 1 = failed]", 101, -1, 100, 3, -1.5, 1.5);
487 PostData(1, fOutputList);
489 switch (fRunModeType) {
491 fOutputListGood = new TList();
492 fOutputListGood->SetOwner(kTRUE);
493 fOutputListBad = new TList();
494 fOutputListBad->SetOwner(kTRUE);
495 PostData(2, fOutputListGood);
496 PostData(3, fOutputListBad);
501 // get the containers
502 fTracksCont = GetParticleContainer("Tracks");
503 fJetsCont = GetJetContainer("Jets");
505 //_____________________________________________________________________________
506 Bool_t AliAnalysisTaskRhoVnModulation::Run()
508 // user exec: execute once for each event
509 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
510 if(!fTracks||!fJets||!fRho) return kFALSE;
511 if(!fLocalInit) fLocalInit = InitializeAnalysis();
512 // reject the event if expected data is missing
513 if(!PassesCuts(InputEvent())) return kFALSE;
514 fLeadingJet = GetLeadingJet(); // store the leading jet
516 fLocalRho->SetVal(fRho->GetVal());
517 // [0][0] psi2a [1,0] psi2c
518 // [0][1] psi3a [1,1] psi3c
519 Double_t vzero[2][2];
520 CalculateEventPlaneVZERO(vzero);
521 /* for the combined vzero event plane
523 * not fully implmemented yet, use with caution ! */
524 Double_t vzeroComb[2];
525 CalculateEventPlaneCombinedVZERO(vzeroComb);
528 CalculateEventPlaneTPC(tpc);
529 Double_t psi2(-1), psi3(-1);
530 // arrays which will hold the fit parameters
531 switch (fDetectorType) { // determine the detector type for the rho fit
532 case kTPC : { psi2 = tpc[0]; psi3 = tpc[1]; } break;
533 case kVZEROA : { psi2 = vzero[0][0]; psi3 = vzero[0][1]; } break;
534 case kVZEROC : { psi2 = vzero[1][0]; psi3 = vzero[1][1]; } break;
535 case kVZEROComb : { psi2 = vzeroComb[0]; psi3 = vzeroComb[1];} break;
538 switch (fFitModulationType) { // do the fits
540 switch (fCollisionType) {
541 case kPythia : { // background is zero for pp jets
542 fFitModulation->FixParameter(0, 0);
543 fLocalRho->SetVal(0);
546 fFitModulation->FixParameter(0, fLocalRho->GetVal());
550 case kV2 : { // only v2
551 if(CorrectRho(psi2, psi3)) {
552 fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
553 if(fUserSuppliedR2) {
554 Double_t r(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
555 if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
557 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
560 case kV3 : { // only v3
561 if(CorrectRho(psi2, psi3)) {
562 if(fUserSuppliedR3) {
563 Double_t r(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
564 if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
566 fProfV3->Fill(fCent, fFitModulation->GetParameter(3));
567 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
570 case kQC2 : { // qc2 analysis
571 if(CorrectRho(psi2, psi3)) {
572 if(fUserSuppliedR2 && fUserSuppliedR3) {
573 // note for the qc method, resolution is REVERSED to go back to v2obs
574 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
575 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
576 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
577 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)*r3);
579 if (fUsePtWeight) { // use weighted weights
580 Double_t dQCnM11 = (fNoEventWeightsForQC) ? 1. : QCnM11();
581 fProfV2->Fill(fCent, fFitModulation->GetParameter(3), dQCnM11);
582 fProfV3->Fill(fCent, fFitModulation->GetParameter(7), dQCnM11);
584 Double_t dQCnM = (fNoEventWeightsForQC) ? 2. : QCnM();
585 fProfV2->Fill(fCent, fFitModulation->GetParameter(3), dQCnM*(dQCnM-1));
586 fProfV3->Fill(fCent, fFitModulation->GetParameter(7), dQCnM*(dQCnM-1));
588 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
592 if(CorrectRho(psi2, psi3)) {
593 if(fUserSuppliedR2 && fUserSuppliedR3) {
594 // note for the qc method, resolution is REVERSED to go back to v2obs
595 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
596 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
597 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
598 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)*r3);
600 if (fUsePtWeight) { // use weighted weights
601 fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM1111()*/);
602 fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM1111()*/);
604 fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
605 fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
608 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
611 if(CorrectRho(psi2, psi3)) {
612 if(fUserSuppliedR2 && fUserSuppliedR3) {
613 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
614 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
615 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r2);
616 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)/r3);
618 fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
619 fProfV3->Fill(fCent, fFitModulation->GetParameter(7));
620 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
624 // if all went well, update the local rho parameter
625 fLocalRho->SetLocalRho(fFitModulation);
626 // fill a number of histograms
627 if(fFillHistograms) FillHistogramsAfterSubtraction(psi2, psi3, vzero, vzeroComb, tpc);
628 if(fFillQAHistograms) FillQAHistograms(InputEvent());
629 // send the output to the connected output container
630 PostData(1, fOutputList);
631 switch (fRunModeType) {
633 PostData(2, fOutputListGood);
634 PostData(3, fOutputListBad);
641 //_____________________________________________________________________________
642 void AliAnalysisTaskRhoVnModulation::CalculateEventPlaneVZERO(Double_t vzero[2][2]) const
644 // get the vzero event plane
645 if(fUseV0EventPlaneFromHeader) { // use the vzero from the header
646 Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0), h(0);
647 vzero[0][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 2, a, b);
648 vzero[1][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 2, c, d);
649 vzero[0][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 3, e, f);
650 vzero[1][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 3, g, h);
653 // grab the vzero event plane without recentering
654 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
655 Double_t qxa2(0), qya2(0), qxc2(0), qyc2(0); // for psi2
656 Double_t qxa3(0), qya3(0), qxc3(0), qyc3(0); // for psi3
657 for(Int_t iVZERO(0); iVZERO < 64; iVZERO++) {
658 Double_t phi(TMath::PiOver4()*(.5+iVZERO%8)), /* eta(0), */ weight(InputEvent()->GetVZEROEqMultiplicity(iVZERO));
659 // (iVZERO<32) ? eta = -3.45+.5*(iVZERO/8) : eta = 4.8-.6*((iVZERO/8)-4);
661 qxa2 += weight*TMath::Cos(2.*phi);
662 qya2 += weight*TMath::Sin(2.*phi);
663 qxa3 += weight*TMath::Cos(3.*phi);
664 qya3 += weight*TMath::Sin(3.*phi);
667 qxc2 += weight*TMath::Cos(2.*phi);
668 qyc2 += weight*TMath::Sin(2.*phi);
669 qxc3 += weight*TMath::Cos(3.*phi);
670 qyc3 += weight*TMath::Sin(3.*phi);
673 vzero[0][0] = .5*TMath::ATan2(qya2, qxa2);
674 vzero[1][0] = .5*TMath::ATan2(qyc2, qxc2);
675 vzero[0][1] = (1./3.)*TMath::ATan2(qya3, qxa3);
676 vzero[1][1] = (1./3.)*TMath::ATan2(qyc3, qxc3);
678 //_____________________________________________________________________________
679 void AliAnalysisTaskRhoVnModulation::CalculateEventPlaneTPC(Double_t* tpc)
681 // grab the TPC event plane
682 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
683 fNAcceptedTracks = 0; // reset the track counter
684 Double_t qx2(0), qy2(0); // for psi2
685 Double_t qx3(0), qy3(0); // for psi3
687 Float_t excludeInEta = -999;
688 if(fExcludeLeadingJetsFromFit > 0 ) { // remove the leading jet from ep estimate
689 if(fLeadingJet) excludeInEta = fLeadingJet->Eta();
691 Int_t iTracks(fTracks->GetEntriesFast());
692 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
693 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
694 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
695 if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta) < GetJetContainer()->GetJetRadius()*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - GetJetContainer()->GetJetRadius() - GetJetContainer()->GetJetEtaMax() ) > 0 )) continue;
697 qx2+= TMath::Cos(2.*track->Phi());
698 qy2+= TMath::Sin(2.*track->Phi());
699 qx3+= TMath::Cos(3.*track->Phi());
700 qy3+= TMath::Sin(3.*track->Phi());
703 tpc[0] = .5*TMath::ATan2(qy2, qx2);
704 tpc[1] = (1./3.)*TMath::ATan2(qy3, qx3);
706 //_____________________________________________________________________________
707 void AliAnalysisTaskRhoVnModulation::CalculateEventPlaneCombinedVZERO(Double_t* comb) const
709 // grab the combined vzero event plane
710 // if(fUseV0EventPlaneFromHeader) { // use the vzero from the header
711 Double_t a(0), b(0), c(0), d(0);
712 comb[0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 10, 2, a, b);
713 comb[1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 10, 3, c, d);
715 // Double_t qx2a(0), qy2a(0), qx2c(0), qy2c(0), qx3a(0), qy3a(0), qx3c(0), qy3c(0);
716 // InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 2, qx2a, qy2a);
717 // InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 2, qx2c, qy2c);
718 // InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 3, qx3a, qy3a);
719 // InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 3, qx3c, qy3c);
720 // FIXME the rest of this function isn't impelmented yet (as of 01-07-2013)
721 // Double_t chi2A(-1), chi2C(-1), chi3A(-1), chi3C(-1); // get chi from the resolution
722 // Double_t qx2(chi2A*chi2A*qx2a+chi2C*chi2C*qx2c);
723 // Double_t qy2(chi2A*chi2A*qy2a+chi2C*chi2C*qy2c);
724 // Double_t qx3(chi3A*chi3A*qx3a+chi3C*chi3C*qx3c);
725 // Double_t qy3(chi3A*chi3A*qy3a+chi3C*chi3C*qy3c);
726 // comb[0] = .5*TMath::ATan2(qy2, qx2);
727 // comb[1] = (1./3.)*TMath::ATan2(qy3, qx3);
730 //_____________________________________________________________________________
731 void AliAnalysisTaskRhoVnModulation::CalculateEventPlaneResolution(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc)
733 // fill the profiles for the resolution parameters
734 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
735 fProfV2Resolution[fInCentralitySelection]->Fill(2., TMath::Cos(2.*(vzero[0][0] - vzero[1][0])));
736 fProfV2Resolution[fInCentralitySelection]->Fill(3., TMath::Cos(2.*(vzero[1][0] - vzero[0][0])));
737 fProfV2Resolution[fInCentralitySelection]->Fill(4., TMath::Cos(2.*(vzero[0][0] - tpc[0])));
738 fProfV2Resolution[fInCentralitySelection]->Fill(5., TMath::Cos(2.*(tpc[0] - vzero[0][0])));
739 fProfV2Resolution[fInCentralitySelection]->Fill(6., TMath::Cos(2.*(vzero[1][0] - tpc[0])));
740 fProfV2Resolution[fInCentralitySelection]->Fill(7., TMath::Cos(2.*(tpc[0] - vzero[1][0])));
741 fProfV3Resolution[fInCentralitySelection]->Fill(2., TMath::Cos(3.*(vzero[0][0] - vzero[1][0])));
742 fProfV3Resolution[fInCentralitySelection]->Fill(3., TMath::Cos(3.*(vzero[1][0] - vzero[0][0])));
743 fProfV3Resolution[fInCentralitySelection]->Fill(4., TMath::Cos(3.*(vzero[0][0] - tpc[0])));
744 fProfV3Resolution[fInCentralitySelection]->Fill(5., TMath::Cos(3.*(tpc[0] - vzero[0][0])));
745 fProfV3Resolution[fInCentralitySelection]->Fill(6., TMath::Cos(3.*(vzero[1][0] - tpc[0])));
746 fProfV3Resolution[fInCentralitySelection]->Fill(7., TMath::Cos(3.*(tpc[0] - vzero[1][0])));
747 // for the resolution of the combined vzero event plane, use two tpc halves as uncorrelated subdetectors
748 Double_t qx2a(0), qy2a(0); // for psi2a, negative eta
749 Double_t qx3a(0), qy3a(0); // for psi3a, negative eta
750 Double_t qx2b(0), qy2b(0); // for psi2a, positive eta
751 Double_t qx3b(0), qy3b(0); // for psi3a, positive eta
753 Int_t iTracks(fTracks->GetEntriesFast());
754 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
755 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
756 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
757 if(track->Eta() < 0 ) {
758 qx2a+= TMath::Cos(2.*track->Phi());
759 qy2a+= TMath::Sin(2.*track->Phi());
760 qx3a+= TMath::Cos(3.*track->Phi());
761 qy3a+= TMath::Sin(3.*track->Phi());
762 } else if (track->Eta() > 0) {
763 qx2b+= TMath::Cos(2.*track->Phi());
764 qy2b+= TMath::Sin(2.*track->Phi());
765 qx3b+= TMath::Cos(3.*track->Phi());
766 qy3b+= TMath::Sin(3.*track->Phi());
770 Double_t tpca2(.5*TMath::ATan2(qy2a, qx2a));
771 Double_t tpca3((1./3.)*TMath::ATan2(qy3a, qx3a));
772 Double_t tpcb2(.5*TMath::ATan2(qy2b, qx2b));
773 Double_t tpcb3((1./3.)*TMath::ATan2(qy3b, qx3b));
774 fProfV2Resolution[fInCentralitySelection]->Fill(8., TMath::Cos(2.*(vzeroComb[0] - tpca2)));
775 fProfV2Resolution[fInCentralitySelection]->Fill(9., TMath::Cos(2.*(vzeroComb[0] - tpcb2)));
776 fProfV2Resolution[fInCentralitySelection]->Fill(10., TMath::Cos(2.*(tpca2 - tpcb2)));
777 fProfV3Resolution[fInCentralitySelection]->Fill(8., TMath::Cos(3.*(vzeroComb[1] - tpca3)));
778 fProfV3Resolution[fInCentralitySelection]->Fill(9., TMath::Cos(3.*(vzeroComb[1] - tpcb3)));
779 fProfV3Resolution[fInCentralitySelection]->Fill(10., TMath::Cos(3.*(tpca3 - tpcb3)));
781 //_____________________________________________________________________________
782 Double_t AliAnalysisTaskRhoVnModulation::CalculateEventPlaneChi(Double_t resEP) const
784 // Get Chi from EP resolution (PRC 58 1671)
785 Double_t chi(2.), delta (1.);
786 for (Int_t i(0); i < 15; i++) {
787 chi = ((TMath::Sqrt(TMath::Pi()/2.)/2.)*chi*exp(-chi*chi/4.)*(TMath::BesselI0(chi*chi/4.)+TMath::BesselI1(chi* chi/4.)) < resEP) ? chi+delta : chi-delta;
792 //_____________________________________________________________________________
793 void AliAnalysisTaskRhoVnModulation::CalculateRandomCone(Float_t &pt, Float_t &eta, Float_t &phi,
794 AliEmcalJet* jet) const
797 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
798 pt = 0; eta = 0; phi = 0;
799 Float_t etaJet(999), phiJet(999), dJet(999); // no jet: same as jet very far away
800 if(jet) { // if a leading jet is given, use its kinematic properties
804 // the random cone acceptance has to equal the jet acceptance
805 // this also insures safety when runnnig on the semi-good tpc runs for 11h data,
806 // where jet acceptance is adjusted to reduced acceptance - hence random cone acceptance as well
807 Float_t minPhi(GetJetContainer()->GetJetPhiMin()), maxPhi(GetJetContainer()->GetJetPhiMax());
808 if(maxPhi > TMath::TwoPi()) maxPhi = TMath::TwoPi();
809 if(minPhi < 0 ) minPhi = 0;
810 Float_t diffRcRJR(TMath::Abs(fRandomConeRadius-GetJetContainer()->GetJetRadius()));
811 // construct a random cone and see if it's far away enough from the leading jet
812 Int_t attempts(1000);
815 eta = gRandom->Uniform(GetJetContainer()->GetJetEtaMin()+diffRcRJR, GetJetContainer()->GetJetEtaMax()-diffRcRJR);
816 phi = gRandom->Uniform(minPhi, maxPhi);
818 dJet = TMath::Sqrt((etaJet-eta)*(etaJet-eta)+(phiJet-phi)*(phiJet-phi));
819 if(dJet > fMinDisanceRCtoLJ) break;
820 else if (attempts == 0) {
821 printf(" > No random cone after 1000 tries, giving up ... !\n");
826 AliVParticle* track = fTracksCont->GetNextAcceptParticle(0);
828 Float_t etaTrack(track->Eta()), phiTrack(track->Phi()), ptTrack(track->Pt());
829 // get distance from cone
830 if(TMath::Abs(phiTrack-phi) > TMath::Abs(phiTrack - phi + TMath::TwoPi())) phiTrack+=TMath::TwoPi();
831 if(TMath::Abs(phiTrack-phi) > TMath::Abs(phiTrack - phi - TMath::TwoPi())) phiTrack-=TMath::TwoPi();
832 if(TMath::Sqrt(TMath::Abs((etaTrack-eta)*(etaTrack-eta)+(phiTrack-phi)*(phiTrack-phi))) <= fRandomConeRadius) pt+=ptTrack;
833 track = fTracksCont->GetNextAcceptParticle();
837 //_____________________________________________________________________________
838 Double_t AliAnalysisTaskRhoVnModulation::CalculateQC2(Int_t harm) {
839 // get the second order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho
840 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
841 Double_t reQ(0), imQ(0), modQ(0), M11(0), M(0);
842 if(fUsePtWeight) { // for the weighted 2-nd order q-cumulant
843 QCnQnk(harm, 1, reQ, imQ); // get the weighted 2-nd order q-vectors
844 modQ = reQ*reQ+imQ*imQ; // get abs Q-squared
845 M11 = QCnM11(); // equals S2,1 - S1,2
846 return (M11 > 0) ? ((modQ - QCnS(1,2))/M11) : -999;
847 } // else return the non-weighted 2-nd order q-cumulant
848 QCnQnk(harm, 0, reQ, imQ); // get the non-weighted 2-nd order q-vectors
849 modQ = reQ*reQ+imQ*imQ; // get abs Q-squared
851 return (M > 1) ? (modQ - M)/(M*(M-1)) : -999;
853 //_____________________________________________________________________________
854 Double_t AliAnalysisTaskRhoVnModulation::CalculateQC4(Int_t harm) {
855 // get the fourth order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho
856 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
857 Double_t reQn1(0), imQn1(0), reQ2n2(0), imQ2n2(0), reQn3(0), imQn3(0), M1111(0), M(0);
858 Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0); // terms of the calculation
859 if(fUsePtWeight) { // for the weighted 4-th order q-cumulant
860 QCnQnk(harm, 1, reQn1, imQn1);
861 QCnQnk(harm*2, 2, reQ2n2, imQ2n2);
862 QCnQnk(harm, 3, reQn3, imQn3);
863 // fill in the terms ...
864 a = (reQn1*reQn1+imQn1*imQn1)*(reQn1*reQn1+imQn1*imQn1);
865 b = reQ2n2*reQ2n2 + imQ2n2*imQ2n2;
866 c = -2.*(reQ2n2*reQn1*reQn1-reQ2n2*imQn1*imQn1+2.*imQ2n2*reQn1*imQn1);
867 d = 8.*(reQn3*reQn1+imQn3*imQn1);
868 e = -4.*QCnS(1,2)*(reQn1*reQn1+imQn1*imQn1);
872 return (M1111 > 0) ? (a+b+c+d+e+f+g)/M1111 : -999;
873 } // else return the unweighted case
874 Double_t reQn(0), imQn(0), reQ2n(0), imQ2n(0);
875 QCnQnk(harm, 0, reQn, imQn);
876 QCnQnk(harm*2, 0, reQ2n, imQ2n);
877 // fill in the terms ...
879 if(M < 4) return -999;
880 a = (reQn*reQn+imQn*imQn)*(reQn*reQn+imQn*imQn);
881 b = reQ2n*reQ2n + imQ2n*imQ2n;
882 c = -2.*(reQ2n*reQn*reQn-reQ2n*imQn*imQn+2.*imQ2n*reQn*imQn);
883 e = -4.*(M-2)*(reQn*reQn+imQn*imQn);
885 return (a+b+c+e+f)/(M*(M-1)*(M-2)*(M-3));
887 //_____________________________________________________________________________
888 void AliAnalysisTaskRhoVnModulation::QCnQnk(Int_t n, Int_t k, Double_t &reQ, Double_t &imQ) {
889 // get the weighted n-th order q-vector, pass real and imaginary part as reference
890 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
892 fNAcceptedTracksQCn = 0;
893 Int_t iTracks(fTracks->GetEntriesFast());
894 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
895 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
896 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
897 fNAcceptedTracksQCn++;
898 // for the unweighted case, k equals zero and the weight doesn't contribute to the equation below
899 reQ += TMath::Power(track->Pt(), k) * TMath::Cos(((double)n)*track->Phi());
900 imQ += TMath::Power(track->Pt(), k) * TMath::Sin(((double)n)*track->Phi());
903 //_____________________________________________________________________________
904 void AliAnalysisTaskRhoVnModulation::QCnDiffentialFlowVectors(
905 TClonesArray* pois, TArrayD* ptBins, Bool_t vpart, Double_t* repn, Double_t* impn,
906 Double_t *mp, Double_t *reqn, Double_t *imqn, Double_t* mq, Int_t n)
908 // get unweighted differential flow vectors
909 Int_t iPois(pois->GetEntriesFast());
911 for(Int_t i(0); i < iPois; i++) {
912 for(Int_t ptBin(0); ptBin < ptBins->GetSize()-1; ptBin++) {
913 AliVTrack* poi = static_cast<AliVTrack*>(pois->At(i));
914 if(PassesCuts(poi)) {
915 if(poi->Pt() >= ptBins->At(ptBin) && poi->Pt() < ptBins->At(ptBin+1)) {
916 // fill the flow vectors assuming that all poi's are in the rp selection (true by design)
917 repn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
918 impn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
920 reqn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
921 imqn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
928 for(Int_t i(0); i < iPois; i++) {
929 for(Int_t ptBin(0); ptBin < ptBins->GetSize()-1; ptBin++) {
930 AliEmcalJet* poi = static_cast<AliEmcalJet*>(pois->At(i));
931 if(PassesCuts(poi)) {
932 Double_t pt(poi->Pt()-poi->Area()*fLocalRho->GetLocalVal(poi->Phi(), GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
933 if(pt >= ptBins->At(ptBin) && pt < ptBins->At(ptBin+1)) {
934 repn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
935 impn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
936 mp[ptBin]++; // qn isn't filled, no overlap between poi's and rp's
943 //_____________________________________________________________________________
944 Double_t AliAnalysisTaskRhoVnModulation::QCnS(Int_t i, Int_t j) {
945 // get the weighted ij-th order autocorrelation correction
946 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
947 if(!fTracks || i <= 0 || j <= 0) return -999;
948 Int_t iTracks(fTracks->GetEntriesFast());
950 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
951 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
952 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
953 Sij+=TMath::Power(track->Pt(), j);
955 return TMath::Power(Sij, i);
957 //_____________________________________________________________________________
958 Double_t AliAnalysisTaskRhoVnModulation::QCnM() {
959 // get multiplicity for unweighted q-cumulants. function QCnQnk should be called first
960 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
961 return (Double_t) fNAcceptedTracksQCn;
963 //_____________________________________________________________________________
964 Double_t AliAnalysisTaskRhoVnModulation::QCnM11() {
965 // get multiplicity weights for the weighted two particle cumulant
966 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
967 return (QCnS(2,1) - QCnS(1,2));
969 //_____________________________________________________________________________
970 Double_t AliAnalysisTaskRhoVnModulation::QCnM1111() {
971 // get multiplicity weights for the weighted four particle cumulant
972 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
973 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));
975 //_____________________________________________________________________________
976 Bool_t AliAnalysisTaskRhoVnModulation::QCnRecovery(Double_t psi2, Double_t psi3) {
977 // decides how to deal with the situation where c2 or c3 is negative
978 // returns kTRUE depending on whether or not a modulated rho is used for the jet background
979 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
980 if(TMath::AreEqualAbs(fFitModulation->GetParameter(3), .0, 1e-10) && TMath::AreEqualAbs(fFitModulation->GetParameter(7), .0,1e-10)) {
981 fFitModulation->SetParameter(7, 0);
982 fFitModulation->SetParameter(3, 0);
983 fFitModulation->SetParameter(0, fLocalRho->GetVal());
984 return kTRUE; // v2 and v3 have physical null values
986 switch (fQCRecovery) {
987 case kFixedRho : { // roll back to the original rho
988 fFitModulation->SetParameter(7, 0);
989 fFitModulation->SetParameter(3, 0);
990 fFitModulation->SetParameter(0, fLocalRho->GetVal());
991 return kFALSE; // rho is forced to be fixed
994 Double_t c2(fFitModulation->GetParameter(3));
995 Double_t c3(fFitModulation->GetParameter(7));
996 if( c2 < 0 ) c2 = -1.*TMath::Sqrt(-1.*c2);
997 if( c3 < 0 ) c3 = -1.*TMath::Sqrt(-1.*c3);
998 fFitModulation->SetParameter(3, c2);
999 fFitModulation->SetParameter(7, c3);
1000 return kTRUE; // is this a physical quantity ?
1003 fitModulationType tempType(fFitModulationType); // store temporarily
1004 fFitModulationType = kCombined;
1005 fFitModulation->SetParameter(7, 0);
1006 fFitModulation->SetParameter(3, 0);
1007 Bool_t pass(CorrectRho(psi2, psi3)); // do the fit and all quality checks
1008 fFitModulationType = tempType; // roll back for next event
1011 default : return kFALSE;
1015 //_____________________________________________________________________________
1016 Bool_t AliAnalysisTaskRhoVnModulation::CorrectRho(Double_t psi2, Double_t psi3)
1018 // get rho' -> rho(phi)
1019 // two routines are available, both can be used with or without pt weights
1020 // [1] get vn from q-cumulants or as an integrated value from a user supplied histogram
1021 // in case of cumulants, both cumulants and vn values are stored. in both cases, v2 and v3
1022 // are expected. a check is performed to see if rho has no negative local minimum
1023 // for full description, see Phys. Rev. C 83, 044913
1024 // since the cn distribution has negative values, vn = sqrt(cn) can be imaginary sometimes
1025 // in this case one can either roll back to the 'original' rixed rho, do a fit for vn or take use
1026 // vn = - sqrt(|cn|)
1027 // [2] fitting a fourier expansion to the de/dphi distribution
1028 // the fit can be done with either v2, v3 or a combination.
1029 // in all cases, a cut can be made on the p-value of the chi-squared value of the fit
1030 // and a check can be performed to see if rho has no negative local minimum
1031 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1032 Int_t freeParams(2); // free parameters of the fit (for NDF)
1033 switch (fFitModulationType) { // for approaches where no fitting is required
1035 fFitModulation->FixParameter(4, psi2);
1036 fFitModulation->FixParameter(6, psi3);
1037 fFitModulation->FixParameter(3, CalculateQC2(2)); // set here with cn, vn = sqrt(cn)
1038 fFitModulation->FixParameter(7, CalculateQC2(3));
1039 // first fill the histos of the raw cumulant distribution
1040 if (fUsePtWeight) { // use weighted weights
1041 Double_t dQCnM11 = (fNoEventWeightsForQC) ? 1. : QCnM11();
1042 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), dQCnM11);
1043 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), dQCnM11);
1045 Double_t dQCnM = (fNoEventWeightsForQC) ? 2. : QCnM();
1046 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), dQCnM*(dQCnM-1));
1047 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), dQCnM*(dQCnM-1));
1049 // then see if one of the cn value is larger than zero and vn is readily available
1050 if(fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) {
1051 fFitModulation->FixParameter(3, TMath::Sqrt(fFitModulation->GetParameter(3)));
1052 fFitModulation->FixParameter(7, TMath::Sqrt(fFitModulation->GetParameter(7)));
1053 } else if (!QCnRecovery(psi2, psi3)) return kFALSE; // try to recover the cumulant, this will set v2 and v3
1054 if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { // general check
1055 fFitModulation->SetParameter(7, 0);
1056 fFitModulation->SetParameter(3, 0);
1057 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1063 fFitModulation->FixParameter(4, psi2);
1064 fFitModulation->FixParameter(6, psi3);
1065 fFitModulation->FixParameter(3, CalculateQC4(2)); // set here with cn, vn = sqrt(cn)
1066 fFitModulation->FixParameter(7, CalculateQC4(3));
1067 // first fill the histos of the raw cumulant distribution
1068 if (fUsePtWeight) { // use weighted weights
1069 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
1070 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
1072 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
1073 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
1075 // then see if one of the cn value is larger than zero and vn is readily available
1076 if(fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) {
1077 fFitModulation->FixParameter(3, TMath::Sqrt(fFitModulation->GetParameter(3)));
1078 fFitModulation->FixParameter(7, TMath::Sqrt(fFitModulation->GetParameter(7)));
1079 } else if (!QCnRecovery(psi2, psi3)) return kFALSE; // try to recover the cumulant, this will set v2 and v3
1080 if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { // general check
1081 fFitModulation->SetParameter(7, 0);
1082 fFitModulation->SetParameter(3, 0);
1083 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1087 case kIntegratedFlow : {
1088 // use v2 and v3 values from an earlier iteration over the data
1089 fFitModulation->FixParameter(3, fUserSuppliedV2->GetBinContent(fUserSuppliedV2->GetXaxis()->FindBin(fCent)));
1090 fFitModulation->FixParameter(4, psi2);
1091 fFitModulation->FixParameter(6, psi3);
1092 fFitModulation->FixParameter(7, fUserSuppliedV3->GetBinContent(fUserSuppliedV3->GetXaxis()->FindBin(fCent)));
1093 if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) {
1094 fFitModulation->SetParameter(7, 0);
1095 fFitModulation->SetParameter(3, 0);
1096 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1103 TString detector("");
1104 switch (fDetectorType) {
1105 case kTPC : detector+="TPC";
1107 case kVZEROA : detector+="VZEROA";
1109 case kVZEROC : detector+="VZEROC";
1111 case kVZEROComb : detector+="VZEROComb";
1115 Int_t iTracks(fTracks->GetEntriesFast());
1116 Double_t excludeInEta = -999;
1117 Double_t excludeInPhi = -999;
1118 Double_t excludeInPt = -999;
1119 if(iTracks <= 0 || fLocalRho->GetVal() <= 0 ) return kFALSE; // no use fitting an empty event ...
1120 if(fExcludeLeadingJetsFromFit > 0 ) {
1122 excludeInEta = fLeadingJet->Eta();
1123 excludeInPhi = fLeadingJet->Phi();
1124 excludeInPt = fLeadingJet->Pt();
1127 // check the acceptance of the track selection that will be used
1128 // if one uses e.g. semi-good tpc tracks, accepance in phi is reduced to 0 < phi < 4
1129 // the defaults (-10 < phi < 10) which accept all, are then overwritten
1130 Double_t lowBound(0.), upBound(TMath::TwoPi()); // bounds for fit
1131 if(GetParticleContainer()->GetParticlePhiMin() > lowBound) lowBound = GetParticleContainer()->GetParticlePhiMin();
1132 if(GetParticleContainer()->GetParticlePhiMax() < upBound) upBound = GetParticleContainer()->GetParticlePhiMax();
1134 fHistSwap->Reset(); // clear the histogram
1135 TH1F _tempSwap; // on stack for quick access
1136 TH1F _tempSwapN; // on stack for quick access, bookkeeping histogram
1137 if(fRebinSwapHistoOnTheFly) {
1138 if(fNAcceptedTracks < 49) fNAcceptedTracks = 49; // avoid aliasing effects
1139 _tempSwap = TH1F("_tempSwap", "_tempSwap", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), lowBound, upBound);
1140 if(fUsePtWeightErrorPropagation) _tempSwapN = TH1F("_tempSwapN", "_tempSwapN", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), lowBound, upBound);
1141 if(fUsePtWeight) _tempSwap.Sumw2();
1143 else _tempSwap = *fHistSwap; // now _tempSwap holds the desired histo
1144 // non poissonian error when using pt weights
1145 Double_t totalpts(0.), totalptsquares(0.), totalns(0.);
1146 for(Int_t i(0); i < iTracks; i++) {
1147 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1148 if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta) < GetJetContainer()->GetJetRadius()*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - GetJetContainer()->GetJetRadius() - GetJetContainer()->GetJetEtaMax() ) > 0 )) continue;
1149 if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
1151 _tempSwap.Fill(track->Phi(), track->Pt());
1152 if(fUsePtWeightErrorPropagation) {
1153 totalpts += track->Pt();
1154 totalptsquares += track->Pt()*track->Pt();
1156 _tempSwapN.Fill(track->Phi());
1159 else _tempSwap.Fill(track->Phi());
1161 if(fUsePtWeight && fUsePtWeightErrorPropagation) {
1162 // in the case of pt weights overwrite the poissonian error estimate which is assigned by root by a more sophisticated appraoch
1163 // 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
1164 // 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
1165 // event, for the latter use a poissonian estimate. the two contrubitions are assumed to be uncorrelated
1166 if(totalns < 1) return kFALSE; // not one track passes the cuts
1167 for(Int_t l = 0; l < _tempSwap.GetNbinsX(); l++) {
1168 if(_tempSwapN.GetBinContent(l+1) == 0) {
1169 _tempSwap.SetBinContent(l+1,0);
1170 _tempSwap.SetBinError(l+1,0);
1173 Double_t vartimesnsq = totalptsquares*totalns - totalpts*totalpts;
1174 Double_t variance = vartimesnsq/(totalns*(totalns-1.));
1175 Double_t SDOMSq = variance / _tempSwapN.GetBinContent(l+1);
1176 Double_t SDOMSqOverMeanSq = SDOMSq * _tempSwapN.GetBinContent(l+1) * _tempSwapN.GetBinContent(l+1) / (_tempSwapN.GetBinContent(l+1) * _tempSwapN.GetBinContent(l+1));
1177 Double_t poissonfrac = 1./_tempSwapN.GetBinContent(l+1);
1178 Double_t vartotalfrac = SDOMSqOverMeanSq + poissonfrac;
1179 Double_t vartotal = vartotalfrac * _tempSwap.GetBinContent(l+1) * _tempSwap.GetBinContent(l+1);
1180 if(vartotal > 0.0001) _tempSwap.SetBinError(l+1,TMath::Sqrt(vartotal));
1182 _tempSwap.SetBinContent(l+1,0);
1183 _tempSwap.SetBinError(l+1,0);
1189 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1190 switch (fFitModulationType) {
1192 fFitModulation->FixParameter(0, fLocalRho->GetVal() );
1196 fFitModulation->FixParameter(4, psi2);
1200 fFitModulation->FixParameter(4, psi3);
1204 fFitModulation->FixParameter(4, psi2);
1205 fFitModulation->FixParameter(6, psi3);
1208 case kFourierSeries : {
1209 // in this approach, an explicit calculation will be made of vn = sqrt(xn^2+yn^2)
1210 // where x[y] = Integrate[r(phi)cos[sin](n phi)dphi, 0, 2pi]
1211 Double_t cos2(0), sin2(0), cos3(0), sin3(0), sumPt(0);
1212 for(Int_t i(0); i < iTracks; i++) {
1213 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1214 if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
1215 sumPt += track->Pt();
1216 cos2 += track->Pt()*TMath::Cos(2*PhaseShift(track->Phi()-psi2));
1217 sin2 += track->Pt()*TMath::Sin(2*PhaseShift(track->Phi()-psi2));
1218 cos3 += track->Pt()*TMath::Cos(3*PhaseShift(track->Phi()-psi3));
1219 sin3 += track->Pt()*TMath::Sin(3*PhaseShift(track->Phi()-psi3));
1221 fFitModulation->SetParameter(3, TMath::Sqrt(cos2*cos2+sin2*sin2)/fLocalRho->GetVal());
1222 fFitModulation->SetParameter(4, psi2);
1223 fFitModulation->SetParameter(6, psi3);
1224 fFitModulation->SetParameter(7, TMath::Sqrt(cos3*cos3+sin3*sin3)/fLocalRho->GetVal());
1229 // toy mc, just here to check procedure, azimuthal profile is filled from hypothesis so p-value distribution should be flat
1230 Int_t _bins = _tempSwap.GetXaxis()->GetNbins();
1231 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());
1232 _tempFit->SetParameter(0, fFitModulation->GetParameter(0)); // normalization
1233 _tempFit->SetParameter(3, 0.1); // v2
1234 _tempFit->FixParameter(1, 1.); // constant
1235 _tempFit->FixParameter(2, 2.); // constant
1236 _tempFit->FixParameter(5, 3.); // constant
1237 _tempFit->FixParameter(4, fFitModulation->GetParameter(4));
1238 _tempFit->FixParameter(6, fFitModulation->GetParameter(6));
1239 _tempFit->SetParameter(7, 0.1); // v3
1240 _tempSwap.Reset(); // rese bin content
1241 for(int _binsI = 0; _binsI < _bins*_bins; _binsI++) _tempSwap.Fill(_tempFit->GetRandom());
1243 _tempSwap.Fit(fFitModulation, fFitModulationOptions.Data(), "", lowBound, upBound);
1244 // the quality of the fit is evaluated from 1 - the cdf of the chi square distribution
1245 // three methods are available, all with their drawbacks. all are stored, one is selected to do the cut
1246 Int_t NDF(_tempSwap.GetXaxis()->GetNbins()-freeParams);
1247 if(NDF == 0) return kFALSE;
1248 Double_t CDF(1.-ChiSquareCDF(NDF, ChiSquare(_tempSwap, fFitModulation)));
1249 Double_t CDFROOT(1.-ChiSquareCDF(NDF, fFitModulation->GetChisquare()));
1250 Double_t CDFKolmogorov(KolmogorovTest(_tempSwap, fFitModulation));
1251 // fill the values and centrality correlation (redundant but easy on the eyes)
1252 fHistPvalueCDF->Fill(CDF);
1253 fHistPvalueCDFCent->Fill(fCent, CDF);
1254 fHistPvalueCDFROOT->Fill(CDFROOT);
1255 fHistPvalueCDFROOTCent->Fill(fCent, CDFROOT);
1256 fHistKolmogorovTest->Fill(CDFKolmogorov);
1257 fHistChi2ROOTCent->Fill(fCent, fFitModulation->GetChisquare()/((float)NDF));
1258 fHistChi2Cent->Fill(fCent, ChiSquare(_tempSwap, fFitModulation)/((float)NDF));
1259 fHistKolmogorovTestCent->Fill(fCent, CDFKolmogorov);
1260 fHistPChi2Root->Fill(CDFROOT, fFitModulation->GetChisquare()/((float)NDF));
1261 fHistPChi2->Fill(CDF, ChiSquare(_tempSwap, fFitModulation)/((float)NDF));
1262 fHistPKolmogorov->Fill(CDF, CDFKolmogorov);
1264 // variable CDF is used for making cuts, so we fill it with the selected p-value
1265 switch (fFitGoodnessTest) {
1269 case kChi2Poisson : break; // CDF is already CDF
1270 case kKolmogorov : {
1271 CDF = CDFKolmogorov;
1277 // as an additional quality check, see if fitting a control fit has a higher significance
1278 _tempSwap.Fit(fFitControl, fFitModulationOptions.Data(), "", lowBound, upBound);
1279 Double_t CDFControl(-1.);
1280 switch (fFitGoodnessTest) {
1282 CDFControl = 1.-ChiSquareCDF(fFitControl->GetNDF(), fFitModulation->GetChisquare());
1284 case kChi2Poisson : {
1285 CDFControl = 1.-ChiSquareCDF(fFitControl->GetNDF(), ChiSquare(_tempSwap, fFitModulation));
1287 case kKolmogorov : {
1288 CDFControl = KolmogorovTest(_tempSwap, fFitControl);
1292 if(CDFControl > CDF) {
1293 CDF = -1.; // control fit is more significant, so throw out the 'old' fit
1294 fHistRhoStatusCent->Fill(fCent, -1);
1297 if(CDF >= fMinPvalue && CDF <= fMaxPvalue && ( fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) > 0)) { // fit quality. not that although with limited acceptance the fit is performed on just
1298 // part of phase space, the requirement that energy desntiy is larger than zero is applied
1299 // to the FULL spectrum
1300 fHistRhoStatusCent->Fill(fCent, 0.);
1301 // for LOCAL didactic purposes, save the best and the worst fits
1302 // this routine can produce a lot of output histograms (it's not memory 'safe') and will not work on GRID
1303 // since the output will become unmergeable (i.e. different nodes may produce conflicting output)
1304 switch (fRunModeType) {
1306 if(fRandom->Uniform(0, 100) > fPercentageOfFits) break;
1307 static Int_t didacticCounterBest(0);
1308 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
1309 TF1* didacticFit = (TF1*)fFitModulation->Clone(Form("fit_%i_CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
1310 switch(fFitModulationType) {
1312 // to make a nice picture also plot the separate components (v2 and v3) of the fit
1313 // only done for cobined fit where there are actually components to split ...
1314 TF1* v2(new TF1("dfit_kV2", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
1315 v2->SetParameter(0, didacticFit->GetParameter(0)); // normalization
1316 v2->SetParameter(3, didacticFit->GetParameter(3)); // v2
1317 v2->FixParameter(1, 1.); // constant
1318 v2->FixParameter(2, 2.); // constant
1319 v2->FixParameter(4, didacticFit->GetParameter(4)); // psi2
1320 v2->SetLineColor(kGreen);
1321 didacticProfile->GetListOfFunctions()->Add(v2);
1322 TF1* v3(new TF1("dfit_kV3", "[0]*([1]+[2]*[3]*TMath::Cos([5]*(x-[4])))", 0, TMath::TwoPi()));
1323 v3->SetParameter(0, didacticFit->GetParameter(0)); // normalization
1324 v3->SetParameter(3, didacticFit->GetParameter(7)); // v3
1325 v3->FixParameter(1, 1.); // constant
1326 v3->FixParameter(2, 2.); // constant
1327 v3->FixParameter(4, didacticFit->GetParameter(6)); // psi3
1328 v3->FixParameter(5, 3.); // constant
1329 v3->SetLineColor(kYellow);
1330 didacticProfile->GetListOfFunctions()->Add(v3);
1334 didacticProfile->GetListOfFunctions()->Add(didacticFit);
1335 fOutputListGood->Add(didacticProfile);
1336 didacticCounterBest++;
1337 TH2F* didacticSurface = BookTH2F(Form("surface_%s", didacticProfile->GetName()), "#phi", "#eta", 50, 0, TMath::TwoPi(), 50, -1, 1, -1, kFALSE);
1338 for(Int_t i(0); i < iTracks; i++) {
1339 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1340 if(PassesCuts(track)) {
1341 if(fUsePtWeight) didacticSurface->Fill(track->Phi(), track->Eta(), track->Pt());
1342 else didacticSurface->Fill(track->Phi(), track->Eta());
1345 if(fExcludeLeadingJetsFromFit) { // visualize the excluded region
1346 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);
1347 f2->SetParameters(excludeInPt/3.,excludeInPhi,.1,excludeInEta,.1);
1348 didacticSurface->GetListOfFunctions()->Add(f2);
1350 fOutputListGood->Add(didacticSurface);
1354 } else { // if the fit is of poor quality revert to the original rho estimate
1355 switch (fRunModeType) { // again see if we want to save the fit
1357 static Int_t didacticCounterWorst(0);
1358 if(fRandom->Uniform(0, 100) > fPercentageOfFits) break;
1359 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data() ));
1360 TF1* didacticFit = (TF1*)fFitModulation->Clone(Form("fit_%i_p_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data()));
1361 didacticProfile->GetListOfFunctions()->Add(didacticFit);
1362 fOutputListBad->Add(didacticProfile);
1363 didacticCounterWorst++;
1367 switch (fFitModulationType) {
1368 case kNoFit : break; // nothing to do
1369 case kCombined : fFitModulation->SetParameter(7, 0); // no break
1370 case kFourierSeries : fFitModulation->SetParameter(7, 0); // no break
1371 default : { // needs to be done if there was a poor fit
1372 fFitModulation->SetParameter(3, 0);
1373 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1376 if(CDF > -.5) fHistRhoStatusCent->Fill(fCent, 1.);
1377 return kFALSE; // return false if the fit is rejected
1381 //_____________________________________________________________________________
1382 Bool_t AliAnalysisTaskRhoVnModulation::PassesCuts(AliVEvent* event)
1385 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1386 if(!event || !AliAnalysisTaskEmcal::IsEventSelected()) return kFALSE;
1387 if(fSemiCentralInclusive && ! (event->GetTriggerMask() & (ULong64_t(1)<<7))) return kFALSE;
1388 if(TMath::Abs(InputEvent()->GetPrimaryVertex()->GetZ()) > 10.) return kFALSE;
1389 // aod and esd specific checks
1390 switch (fDataType) {
1392 AliESDEvent* esdEvent = static_cast<AliESDEvent*>(InputEvent());
1393 if( (!esdEvent) || (TMath::Abs(esdEvent->GetPrimaryVertexSPD()->GetZ() - esdEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE;
1396 AliAODEvent* aodEvent = static_cast<AliAODEvent*>(InputEvent());
1397 if( (!aodEvent) || (TMath::Abs(aodEvent->GetPrimaryVertexSPD()->GetZ() - aodEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE;
1401 fCent = InputEvent()->GetCentrality()->GetCentralityPercentile("V0M");
1402 if(fCent <= fCentralityClasses->At(0) || fCent >= fCentralityClasses->At(fCentralityClasses->GetSize()-1) || TMath::Abs(fCent-InputEvent()->GetCentrality()->GetCentralityPercentile("TRK")) > 5.) return kFALSE;
1403 // determine centrality class
1404 fInCentralitySelection = -1;
1405 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
1406 if(fCent >= fCentralityClasses->At(i) && fCent <= fCentralityClasses->At(1+i)) {
1407 fInCentralitySelection = i;
1410 if(fInCentralitySelection<0) return kFALSE; // should be null op
1411 if(fExplicitOutlierCut == 2010 || fExplicitOutlierCut == 2011) {
1412 if(!PassesCuts(fExplicitOutlierCut)) return kFALSE;
1414 if(fRho->GetVal() <= 0 ) return kFALSE;
1415 if(fTracks->GetEntries() < 1) return kFALSE;
1416 if(fRunNumber != InputEvent()->GetRunNumber()) {
1417 fRunNumber = InputEvent()->GetRunNumber(); // set the current run number
1418 if(fDebug > 0 ) printf("> NEW RUNNNUMBER DETECTED: %i <\n", fRunNumber);
1419 // determine the run number to see if the track and jet cuts should be refreshed for semi-good TPC runs
1420 Int_t semiGoodTPCruns[] = {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};
1421 for(Int_t i(0); i < (int)(sizeof(semiGoodTPCruns)/sizeof(semiGoodTPCruns[0])); i++) {
1422 if(semiGoodTPCruns[i] == fRunNumber) { // run is semi-good
1423 if(fDebug > 0) printf(" > SEMI-GOOD TPC DETECTED, adjusting acceptance accordingly ... <\n");
1424 AliAnalysisTaskEmcalJet::SetJetPhiLimits(fSemiGoodJetMinPhi, fSemiGoodJetMaxPhi);
1425 AliAnalysisTaskEmcal::SetTrackPhiLimits(fSemiGoodTrackMinPhi, fSemiGoodTrackMaxPhi);
1431 //_____________________________________________________________________________
1432 Bool_t AliAnalysisTaskRhoVnModulation::PassesCuts(Int_t year)
1434 // additional centrality cut based on relation between tpc and global multiplicity
1435 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1436 AliAODEvent* event(dynamic_cast<AliAODEvent*>(InputEvent()));
1437 if(!event) return kFALSE;
1438 Int_t multTPC(0), multGlob(0), nTracks(InputEvent()->GetNumberOfTracks());
1439 for(Int_t iTracks = 0; iTracks < nTracks; iTracks++) {
1440 AliAODTrack* track = event->GetTrack(iTracks);
1441 if(!track) continue;
1442 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
1443 if (track->TestFilterBit(1) && track->Chi2perNDF() > 0.2) multTPC++;
1444 if (!track->TestFilterBit(16) || track->Chi2perNDF() < 0.1) continue;
1445 Double_t b[2] = {-99., -99.};
1446 Double_t bCov[3] = {-99., -99., -99.};
1447 if (track->PropagateToDCA(event->GetPrimaryVertex(), event->GetMagneticField(), 100., b, bCov) && TMath::Abs(b[0]) < 0.3 && TMath::Abs(b[1]) < 0.3) multGlob++;
1449 if(year == 2010 && multTPC > (-40.3+1.22*multGlob) && multTPC < (32.1+1.59*multGlob)) return kTRUE;
1450 if(year == 2011 && multTPC > (-36.73 + 1.48*multGlob) && multTPC < (62.87 + 1.78*multGlob)) return kTRUE;
1453 //_____________________________________________________________________________
1454 Bool_t AliAnalysisTaskRhoVnModulation::PassesCuts(const AliVCluster* cluster) const
1457 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1458 if(!cluster) return kFALSE;
1461 //_____________________________________________________________________________
1462 void AliAnalysisTaskRhoVnModulation::FillHistogramsAfterSubtraction(Double_t psi2, Double_t psi3, Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc)
1465 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1466 FillTrackHistograms();
1467 /* FillClusterHistograms(); */
1468 FillJetHistograms(psi2, psi3);
1469 /* FillCorrectedClusterHistograms(); */
1470 if(fFillQAHistograms) FillEventPlaneHistograms(vzero, vzeroComb, tpc);
1471 FillRhoHistograms();
1472 FillDeltaPtHistograms(psi2, psi3);
1474 //_____________________________________________________________________________
1475 void AliAnalysisTaskRhoVnModulation::FillTrackHistograms() const
1477 // fill track histograms
1478 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1479 Int_t iTracks(fTracks->GetEntriesFast()), iAcceptedTracks(0);
1480 for(Int_t i(0); i < iTracks; i++) {
1481 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1482 if(!PassesCuts(track)) continue;
1484 fHistPicoTrackPt[fInCentralitySelection]->Fill(track->Pt());
1485 if(fFillQAHistograms) FillQAHistograms(track);
1487 fHistPicoTrackMult[fInCentralitySelection]->Fill(iAcceptedTracks);
1489 //_____________________________________________________________________________
1490 void AliAnalysisTaskRhoVnModulation::FillClusterHistograms() const
1492 // fill cluster histograms
1493 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1494 Int_t iClusters(fCaloClusters->GetEntriesFast());
1495 for(Int_t i(0); i < iClusters; i++) {
1496 AliVCluster* cluster = static_cast<AliVCluster*>(fCaloClusters->At(iClusters));
1497 if (!PassesCuts(cluster)) continue;
1498 TLorentzVector clusterLorentzVector;
1499 cluster->GetMomentum(clusterLorentzVector, const_cast<Double_t*>(fVertex));
1500 //fHistClusterPt[fInCentralitySelection]->Fill(clusterLorentzVector.Pt());
1501 //fHistClusterEta[fInCentralitySelection]->Fill(clusterLorentzVector.Eta());
1502 //fHistClusterPhi[fInCentralitySelection]->Fill(clusterLorentzVector.Phi());
1506 //_____________________________________________________________________________
1507 void AliAnalysisTaskRhoVnModulation::FillCorrectedClusterHistograms() const
1509 // fill clusters after hadronic correction FIXME implement
1510 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1512 //_____________________________________________________________________________
1513 void AliAnalysisTaskRhoVnModulation::FillEventPlaneHistograms(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc) const
1515 // fill event plane histograms
1516 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1517 fHistPsiControl->Fill(0.5, vzero[0][0]); // vzero a psi2
1518 fHistPsiControl->Fill(1.5, vzero[1][0]); // vzero c psi2
1519 fHistPsiControl->Fill(2.5, tpc[0]); // tpc psi 2
1520 fHistPsiControl->Fill(5.5, vzero[0][1]); // vzero a psi3
1521 fHistPsiControl->Fill(6.5, vzero[1][1]); // vzero b psi3
1522 fHistPsiControl->Fill(7.5, tpc[1]); // tpc psi 3
1523 fHistPsiVZEROA->Fill(vzero[0][0]);
1524 fHistPsiVZEROC->Fill(vzero[1][0]);
1525 fHistPsiVZERO->Fill(vzeroComb[0]);
1526 fHistPsiTPC->Fill(tpc[0]);
1527 fHistPsiSpread->Fill(0.5, TMath::Abs(vzero[0][0]-vzero[1][0]));
1528 fHistPsiSpread->Fill(1.5, TMath::Abs(vzero[0][0]-tpc[0]));
1529 fHistPsiSpread->Fill(2.5, TMath::Abs(vzero[1][0]-tpc[0]));
1530 // event plane vs centrality QA histo's to check recentering
1531 Double_t TRK(InputEvent()->GetCentrality()->GetCentralityPercentile("TRK"));
1532 Double_t V0M(InputEvent()->GetCentrality()->GetCentralityPercentile("V0M"));
1533 fHistPsiVZEROAV0M->Fill(V0M, vzero[0][0]);
1534 fHistPsiVZEROCV0M->Fill(V0M, vzero[1][0]);
1535 fHistPsiVZEROVV0M->Fill(V0M, vzeroComb[0]);
1536 fHistPsiTPCiV0M->Fill(V0M, tpc[0]);
1537 fHistPsiVZEROATRK->Fill(TRK, vzero[0][0]);
1538 fHistPsiVZEROCTRK->Fill(TRK, vzero[1][0]);
1539 fHistPsiVZEROTRK->Fill(TRK, vzeroComb[0]);
1540 fHistPsiTPCTRK->Fill(TRK, tpc[0]);
1542 //_____________________________________________________________________________
1543 void AliAnalysisTaskRhoVnModulation::FillRhoHistograms()
1545 // fill rho histograms
1546 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1547 fHistRhoPackage[fInCentralitySelection]->Fill(fLocalRho->GetVal()); // save the rho estimate from the emcal jet package
1548 // get multiplicity FIXME inefficient
1549 Int_t iJets(fJets->GetEntriesFast());
1550 Double_t rho(fLocalRho->GetLocalVal(TMath::Pi(), TMath::Pi(), fLocalRho->GetVal()));
1551 fHistRho[fInCentralitySelection]->Fill(rho);
1552 fHistRhoVsMult->Fill(fTracks->GetEntries(), rho);
1553 fHistRhoVsCent->Fill(fCent, rho);
1554 for(Int_t i(0); i < iJets; i++) {
1555 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
1556 if(!PassesCuts(jet)) continue;
1557 fHistRhoAVsMult->Fill(fTracks->GetEntries(), rho * jet->Area());
1558 fHistRhoAVsCent->Fill(fCent, rho * jet->Area());
1561 //_____________________________________________________________________________
1562 void AliAnalysisTaskRhoVnModulation::FillDeltaPtHistograms(Double_t psi2, Double_t psi3) const
1564 // fill delta pt histograms
1565 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1567 const Float_t areaRC = fRandomConeRadius*fRandomConeRadius*TMath::Pi();
1568 // we're retrieved the leading jet, now get a random cone
1569 for(i = 0; i < fMaxCones; i++) {
1570 Float_t pt(0), eta(0), phi(0);
1571 // get a random cone without constraints on leading jet position
1572 CalculateRandomCone(pt, eta, phi, 0x0);
1574 if(fFillQAHistograms) fHistRCPhiEta[fInCentralitySelection]->Fill(phi, eta);
1575 fHistRhoVsRCPt[fInCentralitySelection]->Fill(pt, fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal())*areaRC);
1576 fHistRCPt[fInCentralitySelection]->Fill(pt);
1577 fHistDeltaPtDeltaPhi2[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1578 fHistDeltaPtDeltaPhi3[fInCentralitySelection]->Fill(PhaseShift(phi-psi3, 3.), pt - areaRC*fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1580 // get a random cone excluding leading jet area
1581 CalculateRandomCone(pt, eta, phi, fLeadingJet);
1583 if(fFillQAHistograms) fHistRCPhiEtaExLJ[fInCentralitySelection]->Fill(phi, eta);
1584 fHistRhoVsRCPtExLJ[fInCentralitySelection]->Fill(pt, fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal())*areaRC);
1585 fHistRCPtExLJ[fInCentralitySelection]->Fill(pt);
1586 fHistDeltaPtDeltaPhi2ExLJ[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1587 fHistDeltaPtDeltaPhi3ExLJ[fInCentralitySelection]->Fill(PhaseShift(phi-psi3, 3.), pt - areaRC*fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1591 //_____________________________________________________________________________
1592 void AliAnalysisTaskRhoVnModulation::FillJetHistograms(Double_t psi2, Double_t psi3)
1594 // fill jet histograms
1595 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1596 Int_t iJets(fJets->GetEntriesFast());
1597 for(Int_t i(0); i < iJets; i++) {
1598 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
1599 if(PassesCuts(jet)) {
1600 Double_t pt(jet->Pt()), area(jet->Area()), eta(jet->Eta()), phi(jet->Phi());
1601 Double_t rho(fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1602 fHistJetPtRaw[fInCentralitySelection]->Fill(pt);
1603 fHistJetPt[fInCentralitySelection]->Fill(pt-area*rho);
1604 if(fFillQAHistograms) fHistJetEtaPhi[fInCentralitySelection]->Fill(eta, phi);
1605 fHistJetPtArea[fInCentralitySelection]->Fill(pt-area*rho, area);
1606 fHistJetPsi2Pt[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt-area*rho);
1607 fHistJetPsi3Pt[fInCentralitySelection]->Fill(PhaseShift(phi-psi3, 3.), pt-area*rho);
1608 fHistJetPtConstituents[fInCentralitySelection]->Fill(pt-area*rho, jet->Nch());
1609 fHistJetEtaRho[fInCentralitySelection]->Fill(eta, pt/area);
1610 if(fSubtractJetPt) jet->SetPtSub(pt-area*rho); // if requested, save the subtracted jet pt
1611 } else if(fSubtractJetPt) jet->SetPtSub(-999.);
1614 //_____________________________________________________________________________
1615 void AliAnalysisTaskRhoVnModulation::FillQAHistograms(AliVTrack* vtrack) const
1617 // fill qa histograms for pico tracks
1619 AliPicoTrack* track = static_cast<AliPicoTrack*>(vtrack);
1620 fHistRunnumbersPhi->Fill(fMappedRunNumber, track->Phi());
1621 fHistRunnumbersEta->Fill(fMappedRunNumber, track->Eta());
1622 Int_t type((int)(track->GetTrackType()));
1625 fHistPicoCat1[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1628 fHistPicoCat2[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1631 fHistPicoCat3[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1636 //_____________________________________________________________________________
1637 void AliAnalysisTaskRhoVnModulation::FillQAHistograms(AliVEvent* vevent)
1639 // fill qa histograms for events
1641 fHistVertexz->Fill(vevent->GetPrimaryVertex()->GetZ());
1642 fHistCentrality->Fill(fCent);
1643 Int_t runNumber(InputEvent()->GetRunNumber());
1644 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};
1645 for(fMappedRunNumber = 0; fMappedRunNumber < 64; fMappedRunNumber++) {
1646 if(runs[fMappedRunNumber]==runNumber) break;
1649 //_____________________________________________________________________________
1650 void AliAnalysisTaskRhoVnModulation::FillAnalysisSummaryHistogram() const
1652 // fill the analysis summary histrogram, saves all relevant analysis settigns
1653 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1654 fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fJetRadius");
1655 fHistAnalysisSummary->SetBinContent(2, GetJetContainer()->GetJetRadius());
1656 fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fJetEtaMin");
1657 fHistAnalysisSummary->SetBinContent(3, GetJetContainer()->GetJetEtaMin());
1658 fHistAnalysisSummary->GetXaxis()->SetBinLabel(4, "fJetEtaMax");
1659 fHistAnalysisSummary->SetBinContent(4, GetJetContainer()->GetJetEtaMax());
1660 fHistAnalysisSummary->GetXaxis()->SetBinLabel(5, "fJetPhiMin");
1661 fHistAnalysisSummary->SetBinContent(5, GetJetContainer()->GetJetPhiMin());
1662 fHistAnalysisSummary->GetXaxis()->SetBinLabel(6, "fJetPhiMax");
1663 fHistAnalysisSummary->SetBinContent(6, GetJetContainer()->GetJetPhiMin());
1664 fHistAnalysisSummary->GetXaxis()->SetBinLabel(16, "fForceBeamType");
1665 fHistAnalysisSummary->SetBinContent(16, fForceBeamType);
1666 fHistAnalysisSummary->GetXaxis()->SetBinLabel(17, "fMinCent");
1667 fHistAnalysisSummary->SetBinContent(17, fMinCent);
1668 fHistAnalysisSummary->GetXaxis()->SetBinLabel(18, "fMaxCent");
1669 fHistAnalysisSummary->SetBinContent(18, fMaxCent);
1670 fHistAnalysisSummary->GetXaxis()->SetBinLabel(19, "fMinVz");
1671 fHistAnalysisSummary->SetBinContent(19, fMinVz);
1672 fHistAnalysisSummary->GetXaxis()->SetBinLabel(20, "fMaxVz");
1673 fHistAnalysisSummary->SetBinContent(20, fMaxVz);
1674 fHistAnalysisSummary->GetXaxis()->SetBinLabel(21, "fOffTrigger");
1675 fHistAnalysisSummary->SetBinContent(21, fOffTrigger);
1676 fHistAnalysisSummary->GetXaxis()->SetBinLabel(33, "fRandomConeRadius");
1677 fHistAnalysisSummary->SetBinContent(33, fRandomConeRadius);
1678 fHistAnalysisSummary->GetXaxis()->SetBinLabel(34, "fitModulationType");
1679 fHistAnalysisSummary->SetBinContent(34, (int)fFitModulationType);
1680 fHistAnalysisSummary->GetXaxis()->SetBinLabel(35, "runModeType");
1681 fHistAnalysisSummary->SetBinContent(35, (int)fRunModeType);
1682 fHistAnalysisSummary->GetXaxis()->SetBinLabel(36, "data type");
1683 fHistAnalysisSummary->SetBinContent(36, (int)fDataType);
1684 fHistAnalysisSummary->GetXaxis()->SetBinLabel(37, "iterator");
1685 fHistAnalysisSummary->SetBinContent(37, 1.);
1686 fHistAnalysisSummary->GetXaxis()->SetBinLabel(38, "fMinPvalue");
1687 fHistAnalysisSummary->SetBinContent(38, fMinPvalue);
1688 fHistAnalysisSummary->GetXaxis()->SetBinLabel(39, "fMaxPvalue");
1689 fHistAnalysisSummary->SetBinContent(39, fMaxPvalue);
1690 fHistAnalysisSummary->GetXaxis()->SetBinLabel(40, "fExcludeLeadingJetsFromFit");
1691 fHistAnalysisSummary->SetBinContent(40, fExcludeLeadingJetsFromFit);
1692 fHistAnalysisSummary->GetXaxis()->SetBinLabel(41, "fRebinSwapHistoOnTheFly");
1693 fHistAnalysisSummary->SetBinContent(41, (int)fRebinSwapHistoOnTheFly);
1694 fHistAnalysisSummary->GetXaxis()->SetBinLabel(42, "fUsePtWeight");
1695 fHistAnalysisSummary->SetBinContent(42, (int)fUsePtWeight);
1696 fHistAnalysisSummary->GetXaxis()->SetBinLabel(43, "fMinLeadingHadronPt");
1697 fHistAnalysisSummary->SetBinContent(43, fMinLeadingHadronPt);
1698 fHistAnalysisSummary->GetXaxis()->SetBinLabel(44, "fExplicitOutlierCut");
1699 fHistAnalysisSummary->SetBinContent(44, fExplicitOutlierCut);
1700 fHistAnalysisSummary->GetXaxis()->SetBinLabel(45, "fLocalJetMinEta");
1701 fHistAnalysisSummary->SetBinContent(45,fLocalJetMinEta );
1702 fHistAnalysisSummary->GetXaxis()->SetBinLabel(46, "fLocalJetMaxEta");
1703 fHistAnalysisSummary->SetBinContent(46, fLocalJetMaxEta);
1704 fHistAnalysisSummary->GetXaxis()->SetBinLabel(47, "fLocalJetMinPhi");
1705 fHistAnalysisSummary->SetBinContent(47, fLocalJetMinPhi);
1706 fHistAnalysisSummary->GetXaxis()->SetBinLabel(48, "fLocalJetMaxPhi");
1707 fHistAnalysisSummary->SetBinContent(48, fLocalJetMaxPhi);
1708 fHistAnalysisSummary->GetXaxis()->SetBinLabel(49, "fSoftTrackMinPt");
1709 fHistAnalysisSummary->SetBinContent(49, fSoftTrackMinPt);
1710 fHistAnalysisSummary->GetXaxis()->SetBinLabel(50, "fSoftTrackMaxPt");
1711 fHistAnalysisSummary->SetBinContent(50, fSoftTrackMaxPt);
1712 fHistAnalysisSummary->GetXaxis()->SetBinLabel(51, "fMaxCones");
1713 fHistAnalysisSummary->SetBinContent(51, fMaxCones);
1714 fHistAnalysisSummary->GetXaxis()->SetBinLabel(52, "fUseScaledRho");
1715 fHistAnalysisSummary->SetBinContent(52, fUseScaledRho);
1717 //_____________________________________________________________________________
1718 void AliAnalysisTaskRhoVnModulation::Terminate(Option_t *)
1721 switch (fRunModeType) {
1723 printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1724 if(fFillQAHistograms) {
1725 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};
1726 for(Int_t i(0); i < 64; i++) {
1727 fHistRunnumbersPhi->GetXaxis()->SetBinLabel(i+1, Form("%i", runs[i]));
1728 fHistRunnumbersEta->GetXaxis()->SetBinLabel(i+1, Form("%i", runs[i]));
1730 fHistRunnumbersPhi->GetXaxis()->SetBinLabel(65, "undetermined");
1731 fHistRunnumbersEta->GetXaxis()->SetBinLabel(65, "undetermined");
1733 AliAnalysisTaskRhoVnModulation::Dump();
1734 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));
1739 //_____________________________________________________________________________
1740 void AliAnalysisTaskRhoVnModulation::SetModulationFit(TF1* fit)
1742 // set modulation fit
1743 if (fFitModulation) delete fFitModulation;
1744 fFitModulation = fit;
1746 //_____________________________________________________________________________
1747 void AliAnalysisTaskRhoVnModulation::SetUseControlFit(Bool_t c)
1750 if (fFitControl) delete fFitControl;
1752 fFitControl = new TF1("controlFit", "pol0", 0, TMath::TwoPi());
1753 } else fFitControl = 0x0;
1755 //_____________________________________________________________________________
1756 TH1F* AliAnalysisTaskRhoVnModulation::GetResolutionFromOuptutFile(detectorType det, Int_t h, TArrayD* cen)
1758 // INTERFACE METHOD FOR OUTPUTFILE
1759 // get the detector resolution, user has ownership of the returned histogram
1761 printf(" > Please add fOutputList first < \n");
1765 (cen) ? r = new TH1F("R", "R", cen->GetSize()-1, cen->GetArray()) : r = new TH1F("R", "R", 10, 0, 10);
1766 if(!cen) r->GetXaxis()->SetTitle("number of centrality bin");
1767 r->GetYaxis()->SetTitle(Form("Resolution #Psi_{%i}", h));
1768 for(Int_t i(0); i < 10; i++) {
1769 TProfile* temp((TProfile*)fOutputList->FindObject(Form("fProfV%iResolution_%i", h, i)));
1771 Double_t a(temp->GetBinContent(3)), b(temp->GetBinContent(5)), c(temp->GetBinContent(7));
1772 Double_t d(temp->GetBinContent(9)), e(temp->GetBinContent(10)), f(temp->GetBinContent(11));
1773 Double_t _a(temp->GetBinError(3)), _b(temp->GetBinError(5)), _c(temp->GetBinError(7));
1774 Double_t _d(temp->GetBinError(9)), _e(temp->GetBinError(10)), _f(temp->GetBinError(11));
1775 if(a <= 0 || b <= 0 || c <= 0 || d <= 0 || e <= 0 || f <= 0) continue;
1778 r->SetBinContent(1+i, TMath::Sqrt((a*b)/c));
1779 if(i==0) r->SetNameTitle("VZEROA resolution", "VZEROA resolution");
1780 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
1783 r->SetBinContent(1+i, TMath::Sqrt((a*c)/b));
1784 if(i==0) r->SetNameTitle("VZEROC resolution", "VZEROC resolution");
1785 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
1788 r->SetBinContent(1+i, TMath::Sqrt((b*c)/a));
1789 if(i==0) r->SetNameTitle("TPC resolution", "TPC resolution");
1790 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
1793 r->SetBinContent(1+i, TMath::Sqrt((d*e)/f));
1794 if(i==0) r->SetNameTitle("VZEROComb resolution", "VZEROComb resolution");
1795 r->SetBinError(1+i, TMath::Sqrt(_d*_d+_e*_e+_f*_f));
1802 //_____________________________________________________________________________
1803 TH1F* AliAnalysisTaskRhoVnModulation::CorrectForResolutionDiff(TH1F* v, detectorType det, TArrayD* cen, Int_t c, Int_t h)
1805 // INTERFACE METHOD FOR OUTPUT FILE
1806 // correct the supplied differential vn histogram v for detector resolution
1807 TH1F* r(GetResolutionFromOuptutFile(det, h, cen));
1809 printf(" > Couldn't find resolution < \n");
1812 Double_t res(1./r->GetBinContent(1+r->FindBin(c)));
1813 TF1* line = new TF1("line", "pol0", 0, 200);
1814 line->SetParameter(0, res);
1818 //_____________________________________________________________________________
1819 TH1F* AliAnalysisTaskRhoVnModulation::CorrectForResolutionInt(TH1F* v, detectorType det, TArrayD* cen, Int_t h)
1821 // INTERFACE METHOD FOR OUTPUT FILE
1822 // correct the supplied intetrated vn histogram v for detector resolution
1823 // integrated vn must have the same centrality binning as the resolotion correction
1824 TH1F* r(GetResolutionFromOuptutFile(det, h, cen));
1828 //_____________________________________________________________________________
1829 TH1F* AliAnalysisTaskRhoVnModulation::GetDifferentialQC(TProfile* refCumulants, TProfile* diffCumlants, TArrayD* ptBins, Int_t h)
1831 // get differential QC
1832 Double_t r(refCumulants->GetBinContent(h-1)); // v2 reference flow
1833 if(r > 0) r = TMath::Sqrt(r);
1834 TH1F* qc = new TH1F(Form("QC2v%i", h), Form("QC2v%i", h), ptBins->GetSize()-1, ptBins->GetArray());
1835 Double_t a(0), b(0), c(0); // dummy variables
1836 for(Int_t i(0); i < ptBins->GetSize(); i++) {
1838 a = diffCumlants->GetBinContent(1+i);
1839 b = diffCumlants->GetBinError(1+i);
1841 qc->SetBinContent(1+i, c);
1842 (a <= 0 || b <= 0) ? qc->SetBinError(1+i, b) : qc->SetBinError(1+i, TMath::Sqrt(c*c*b*b/(a*a)));
1848 //_____________________________________________________________________________