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), fExpectedRuns(0), fExpectedSemiGoodRuns(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), fNameSmallRho(""), fCachedRho(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), fHistUndeterminedRunQA(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), fExpectedRuns(0), fExpectedSemiGoodRuns(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), fNameSmallRho(""), fCachedRho(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), fHistUndeterminedRunQA(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(fExpectedRuns) delete fExpectedRuns;
179 if(fExpectedSemiGoodRuns) delete fExpectedSemiGoodRuns;
180 if(fFitControl) delete fFitControl;
182 //_____________________________________________________________________________
183 void AliAnalysisTaskRhoVnModulation::ExecOnce()
186 fLocalRho = new AliLocalRhoParameter(fLocalRhoName.Data(), 0);
188 if(!(InputEvent()->FindListObject(fLocalRho->GetName()))) {
189 InputEvent()->AddObject(fLocalRho);
191 AliFatal(Form("%s: Container with name %s already present. Aborting", GetName(), fLocalRho->GetName()));
194 AliAnalysisTaskEmcalJet::ExecOnce(); // init the base class
195 AliAnalysisTaskEmcalJet::SetVzRange(-1.*fAbsVertexZ, fAbsVertexZ);
197 // unscaled rho has been retrieved by the parent class, now we retrieve rho scaled
198 fRho = dynamic_cast<AliRhoParameter*>(InputEvent()->FindListObject(Form("%s_Scaled", fRho->GetName())));
200 AliFatal(Form("%s: Couldn't find container for scaled rho. Aborting !", GetName()));
203 if(!GetJetContainer()) AliFatal(Form("%s: Couldn't find jet container. Aborting !", GetName()));
205 //_____________________________________________________________________________
206 Bool_t AliAnalysisTaskRhoVnModulation::InitializeAnalysis()
208 // initialize the anaysis
209 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
210 if(fRandomConeRadius <= 0) fRandomConeRadius = GetJetContainer()->GetJetRadius();
211 if(fMaxCones <= 0) fMaxCones = TMath::Nint(1.8*TMath::TwoPi()/(TMath::Pi()*fRandomConeRadius*fRandomConeRadius));
212 if(fLocalJetMinEta > -10 && fLocalJetMaxEta > -10) GetJetContainer()->SetJetEtaLimits(fLocalJetMinEta, fLocalJetMaxEta);
213 if(fLocalJetMinPhi > -10 && fLocalJetMaxPhi > -10) GetJetContainer()->SetJetPhiLimits(fLocalJetMinPhi, fLocalJetMaxPhi);
214 if(fMinDisanceRCtoLJ==0) fMinDisanceRCtoLJ = .5*GetJetRadius();
215 if(dynamic_cast<AliAODEvent*>(InputEvent())) fDataType = kAOD; // determine the datatype
216 else if(dynamic_cast<AliESDEvent*>(InputEvent())) fDataType = kESD;
217 fHistAnalysisSummary->SetBinContent(36, (int)fDataType);
218 if(!fRandom) fRandom = new TRandom3(0); // get a randomized if one hasn't been user-supplied
219 switch (fFitModulationType) {
220 case kNoFit : { SetModulationFit(new TF1("fix_kNoFit", "[0]", 0, TMath::TwoPi())); } break;
222 SetModulationFit(new TF1("fit_kV2", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
223 fFitModulation->SetParameter(0, 0.); // normalization
224 fFitModulation->SetParameter(3, 0.2); // v2
225 fFitModulation->FixParameter(1, 1.); // constant
226 fFitModulation->FixParameter(2, 2.); // constant
229 SetModulationFit(new TF1("fit_kV3", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
230 fFitModulation->SetParameter(0, 0.); // normalization
231 fFitModulation->SetParameter(3, 0.2); // v3
232 fFitModulation->FixParameter(1, 1.); // constant
233 fFitModulation->FixParameter(2, 3.); // constant
235 default : { // for the combined fit, the 'direct fourier series' or the user supplied vn values we use v2 and v3
236 SetModulationFit(new TF1("fit_kCombined", "[0]*([1]+[2]*([3]*TMath::Cos([2]*(x-[4]))+[7]*TMath::Cos([5]*(x-[6]))))", 0, TMath::TwoPi()));
237 fFitModulation->SetParameter(0, 0.); // normalization
238 fFitModulation->SetParameter(3, 0.2); // v2
239 fFitModulation->FixParameter(1, 1.); // constant
240 fFitModulation->FixParameter(2, 2.); // constant
241 fFitModulation->FixParameter(5, 3.); // constant
242 fFitModulation->SetParameter(7, 0.2); // v3
245 switch (fRunModeType) {
246 case kGrid : { fFitModulationOptions += "N0"; } break;
249 FillAnalysisSummaryHistogram();
252 //_____________________________________________________________________________
253 TH1F* AliAnalysisTaskRhoVnModulation::BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Int_t c, Bool_t append)
255 // book a TH1F and connect it to the output container
256 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
257 if(fReduceBinsXByFactor > 0 ) bins = TMath::Nint(bins/fReduceBinsXByFactor);
258 if(!fOutputList) return 0x0;
260 if(c!=-1) { // format centrality dependent histograms accordingly
261 name = Form("%s_%i", name, c);
262 title += Form("_%i-%i", (int)(fCentralityClasses->At(c)), (int)(fCentralityClasses->At((1+c))));
264 title += Form(";%s;[counts]", x);
265 TH1F* histogram = new TH1F(name, title.Data(), bins, min, max);
267 if(append) fOutputList->Add(histogram);
270 //_____________________________________________________________________________
271 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)
273 // book a TH2F and connect it to the output container
274 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
275 if(fReduceBinsXByFactor > 0 ) binsx = TMath::Nint(binsx/fReduceBinsXByFactor);
276 if(fReduceBinsYByFactor > 0 ) binsy = TMath::Nint(binsy/fReduceBinsYByFactor);
277 if(!fOutputList) return 0x0;
279 if(c!=-1) { // format centrality dependent histograms accordingly
280 name = Form("%s_%i", name, c);
281 title += Form("_%i-%i", (int)fCentralityClasses->At(c), (int)(fCentralityClasses->At((1+c))));
283 title += Form(";%s;%s", x, y);
284 TH2F* histogram = new TH2F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy);
286 if(append) fOutputList->Add(histogram);
289 //_____________________________________________________________________________
290 void AliAnalysisTaskRhoVnModulation::UserCreateOutputObjects()
292 // create output objects
293 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
294 fOutputList = new TList();
295 fOutputList->SetOwner(kTRUE);
296 if(!fCentralityClasses) { // classes must be defined at this point
297 Double_t c[] = {0., 20., 40., 60., 80., 100.};
298 fCentralityClasses = new TArrayD(sizeof(c)/sizeof(c[0]), c);
300 if(!fExpectedRuns) { // expected runs must be defined at this point
301 Int_t r[] = {167813, 167988, 168066, 168068, 168069, 168076, 168104, 168212, 168311, 168322, 168325, 168341, 168361, 168362, 168458, 168460, 168461, 168992, 169091, 169094, 169138, 169143, 169167, 169417, 169835, 169837, 169838, 169846, 169855, 169858, 169859, 169923, 169956, 170027, 170036, 170081, /* up till here original good TPC list */169975, 169981, 170038, 170040, 170083, 170084, 170085, 170088, 170089, 170091, 170152, 170155, 170159, 170163, 170193, 170195, 170203, 170204, 170205, 170228, 170230, 170264, 170268, 170269, 170270, 170306, 170308, 170309, /* original semi-good tpc list */169415, 169411, 169035, 168988, 168984, 168826, 168777, 168512, 168511, 168467, 168464, 168342, 168310, 168115, 168108, 168107, 167987, 167915, 167903, /*new runs, good according to RCT */ 169238, 169160, 169156, 169148, 169145, 169144 /* run swith missing OROC 8 but seem ok in QA */};
302 fExpectedRuns = new TArrayI(sizeof(r)/sizeof(r[0]), r);
304 if(!fExpectedSemiGoodRuns) {
305 Int_t r[] = {169975, 169981, 170038, 170040, 170083, 170084, 170085, 170088, 170089, 170091, 170152, 170155, 170159, 170163, 170193, 170195, 170203, 170204, 170205, 170228, 170230, 170264, 170268, 170269, 170270, 170306, 170308, 170309};
306 fExpectedSemiGoodRuns = new TArrayI(sizeof(r)/sizeof(r[0]), r);
309 fHistCentrality = BookTH1F("fHistCentrality", "centrality", 102, -2, 100);
310 fHistVertexz = BookTH1F("fHistVertexz", "vertex z (cm)", 100, -12, 12);
312 // pico track kinematics
313 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
314 fHistPicoTrackPt[i] = BookTH1F("fHistPicoTrackPt", "p_{t} [GeV/c]", 100, 0, 100, i);
315 fHistPicoTrackMult[i] = BookTH1F("fHistPicoTrackMult", "multiplicity", 100, 0, 5000, i);
316 if(fFillQAHistograms) {
317 fHistPicoCat1[i] = BookTH2F("fHistPicoCat1", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
318 fHistPicoCat2[i] = BookTH2F("fHistPicoCat2", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
319 fHistPicoCat3[i] = BookTH2F("fHistPicoCat3", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
322 /* fHistClusterPt[i] = BookTH1F("fHistClusterPt", "p_{t} [GeV/c]", 100, 0, 100, i); */
323 /* fHistClusterPhi[i] = BookTH1F("fHistClusterPhi", "#phi", 100, 0, TMath::TwoPi(), i); */
324 /* fHistClusterEta[i] = BookTH1F("fHistClusterEta", "#eta", 100, -5, 5); */
326 // emcal kinematics after hadronic correction
327 /* fHistClusterCorrPt[i] = BookTH1F("fHistClusterCorrPt", "p_{t} [GeV/c]", 100, 0, 100, i); */
328 /* fHistClusterCorrPhi[i] = BookTH1F("fHistClusterCorrPhi", "#phi", 100, 0, TMath::TwoPi(), i); */
329 /* fHistClusterCorrEta[i] = BookTH1F("fHistClusterCorrEta", "#eta", 100, -5, 5, i); */
332 if(fFillQAHistograms) {
333 // event plane estimates and quality
334 fHistPsiControl = new TProfile("fHistPsiControl", "fHistPsiControl", 10, 0, 10);
335 fHistPsiControl->Sumw2();
336 fHistPsiSpread = new TProfile("fHistPsiSpread", "fHistPsiSpread", 4, 0, 4);
337 fHistPsiSpread->Sumw2();
338 fHistPsiControl->GetXaxis()->SetBinLabel(1, "<#Psi_{2, VZEROA}>");
339 fHistPsiControl->GetXaxis()->SetBinLabel(2, "<#Psi_{2, VZEROC}>");
340 fHistPsiControl->GetXaxis()->SetBinLabel(3, "<#Psi_{2, TPC}>");
341 fHistPsiControl->GetXaxis()->SetBinLabel(4, "<#Psi_{2, TPC, #eta < 0}>");
342 fHistPsiControl->GetXaxis()->SetBinLabel(5, "<#Psi_{2, TPC, #eta > 0}>");
343 fHistPsiControl->GetXaxis()->SetBinLabel(6, "<#Psi_{3, VZEROA}>");
344 fHistPsiControl->GetXaxis()->SetBinLabel(7, "<#Psi_{3, VZEROC}>");
345 fHistPsiControl->GetXaxis()->SetBinLabel(8, "<#Psi_{3, TPC}>");
346 fHistPsiControl->GetXaxis()->SetBinLabel(9, "<#Psi_{3, TPC, #eta < 0}>");
347 fHistPsiControl->GetXaxis()->SetBinLabel(10, "<#Psi_{3, TPC, #eta > 0}>");
348 fHistPsiSpread->GetXaxis()->SetBinLabel(1, "<#Psi_{2, VZEROA} - #Psi_{2, VZEROC}>");
349 fHistPsiSpread->GetXaxis()->SetBinLabel(2, "<#Psi_{2, VZEROC} - #Psi_{2, TPC}>");
350 fHistPsiSpread->GetXaxis()->SetBinLabel(3, "<#Psi_{2, VZEROC} - #Psi_{2, TPC}>");
351 fHistPsiSpread->GetXaxis()->SetBinLabel(4, "<#Psi_{2, TPC, #eta < 0} - #Psi_{2, TPC, #eta > 0}>");
352 fOutputList->Add(fHistPsiControl);
353 fOutputList->Add(fHistPsiSpread);
354 fHistPsiVZEROA = BookTH1F("fHistPsiVZEROA", "#Psi_{VZEROA}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
355 fHistPsiVZEROC = BookTH1F("fHistPsiVZEROC", "#Psi_{VZEROC}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
356 fHistPsiVZERO = BookTH1F("fHistPsiVZERO", "#Psi_{VZERO}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
357 fHistPsiTPC = BookTH1F("fHistPsiTPC", "#Psi_{TPC}", 40, -.5*TMath::Pi(), .5*TMath::Pi());
358 fHistPsiVZEROAV0M = BookTH2F("fHistPsiVZEROAV0M", "V0M", "#Psi_{2, VZEROA}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
359 fHistPsiVZEROCV0M = BookTH2F("fHistPsiVZEROCV0M", "V0M", "#Psi_{2, VZEROC}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
360 fHistPsiVZEROVV0M = BookTH2F("fHistPsiVZEROV0M", "V0M", "#Psi_{2, VZERO}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
361 fHistPsiTPCiV0M = BookTH2F("fHistPsiTPCV0M", "V0M", "#Psi_{2, TRK}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
362 fHistPsiVZEROATRK = BookTH2F("fHistPsiVZEROATRK", "TRK", "#Psi_{2, VZEROA}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
363 fHistPsiVZEROCTRK = BookTH2F("fHistPsiVZEROCTRK", "TRK", "#Psi_{2, VZEROC}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
364 fHistPsiVZEROTRK = BookTH2F("fHistPsiVZEROTRK", "TRK", "#Psi_{2, VZERO}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
365 fHistPsiTPCTRK = BookTH2F("fHistPsiTPCTRK", "TRK", "#Psi_{2, TRK}", 60, 0, 60, 40, -.5*TMath::Pi(), .5*TMath::Pi());
368 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i ++) {
369 fHistRhoPackage[i] = BookTH1F("fHistRhoPackage", "#rho [GeV/c]", 100, 0, 150, i);
370 fHistRho[i] = BookTH1F("fHistRho", "#rho [GeV/c]", 100, 0, 150, i);
372 fHistRhoVsMult = BookTH2F("fHistRhoVsMult", "multiplicity", "#rho [GeV/c]", 100, 0, 4000, 100, 0, 250);
373 fHistRhoVsCent = BookTH2F("fHistRhoVsCent", "centrality", "#rho [GeV/c]", 100, 0, 100, 100, 0, 250);
374 fHistRhoAVsMult = BookTH2F("fHistRhoAVsMult", "multiplicity", "#rho * A (jet) [GeV/c]", 100, 0, 4000, 100, 0, 50);
375 fHistRhoAVsCent = BookTH2F("fHistRhoAVsCent", "centrality", "#rho * A (jet) [GeV/c]", 100, 0, 100, 100, 0, 50);
377 TString detector("");
378 switch (fDetectorType) {
379 case kTPC : detector+="TPC";
381 case kVZEROA : detector+="VZEROA";
383 case kVZEROC : detector+="VZEROC";
385 case kVZEROComb : detector+="VZEROComb";
389 // delta pt distributions
390 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i ++) {
391 if(fFillQAHistograms) fHistRCPhiEta[i] = BookTH2F("fHistRCPhiEta", "#phi (RC)", "#eta (RC)", 40, 0, TMath::TwoPi(), 40, -1, 1, i);
392 fHistRhoVsRCPt[i] = BookTH2F("fHistRhoVsRCPt", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i);
393 fHistRCPt[i] = BookTH1F("fHistRCPt", "p_{t} (RC) [GeV/c]", 130, -20, 150, i);
394 if(fFillQAHistograms) fHistRCPhiEtaExLJ[i] = BookTH2F("fHistRCPhiEtaExLJ", "#phi (RC)", "#eta (RC)", 40, 0, TMath::TwoPi(), 40, -1, 1, i);
395 fHistDeltaPtDeltaPhi2[i] = BookTH2F("fHistDeltaPtDeltaPhi2", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
396 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);
397 fHistRhoVsRCPtExLJ[i] = BookTH2F("fHistRhoVsRCPtExLJ", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i);
398 fHistRCPtExLJ[i] = BookTH1F("fHistRCPtExLJ", "p_{t} (RC) [GeV/c]", 130, -20, 150, i);
399 /* fHistRCPhiEtaRand[i] = BookTH2F("fHistRCPhiEtaRand", "#phi (RC)", "#eta (RC)", 100, 0, TMath::TwoPi(), 100, -1, 1, i); */
400 fHistDeltaPtDeltaPhi2ExLJ[i] = BookTH2F("fHistDeltaPtDeltaPhi2ExLJ", Form("#phi - #Psi_{2, %s}", detector.Data()), "#delta p_{t} [GeV/c]", 40, 0, TMath::Pi(), 400, -70, 130, i);
401 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);
402 /* fHistRhoVsRCPtRand[i] = BookTH2F("fHistRhoVsRCPtRand", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i); */
403 /* fHistRCPtRand[i] = BookTH1F("fHistRCPtRand", "p_{t} (RC) [GeV/c]", 130, -20, 150, i); */
404 /* fHistDeltaPtDeltaPhi2Rand[i] = BookTH2F("fHistDeltaPtDeltaPhi2Rand", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::Pi(), 100, -50, 100, i); */
405 /* fHistDeltaPtDeltaPhi3Rand[i] = BookTH2F("fHistDeltaPtDeltaPhi3Rand", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::TwoPi()/3., 100, -50, 100, i); */
406 // jet histograms (after kinematic cuts)
407 fHistJetPtRaw[i] = BookTH1F("fHistJetPtRaw", "p_{t, jet} RAW [GeV/c]", 200, -50, 150, i);
408 fHistJetPt[i] = BookTH1F("fHistJetPt", "p_{t, jet} [GeV/c]", 350, -100, 250, i);
409 if(fFillQAHistograms) fHistJetEtaPhi[i] = BookTH2F("fHistJetEtaPhi", "#eta", "#phi", 100, -1, 1, 100, 0, TMath::TwoPi(), i);
410 fHistJetPtArea[i] = BookTH2F("fHistJetPtArea", "p_{t, jet} [GeV/c]", "Area", 175, -100, 250, 30, 0, 0.9, i);
411 fHistJetPtConstituents[i] = BookTH2F("fHistJetPtConstituents", "p_{t, jet} [GeV/c]", "Area", 350, -100, 250, 60, 0, 150, i);
412 fHistJetEtaRho[i] = BookTH2F("fHistJetEtaRho", "#eta", "#rho", 100, -1, 1, 100, 0, 300, i);
413 // in plane and out of plane spectra
414 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);
415 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);
416 // profiles for all correlator permutations which are necessary to calculate each second and third order event plane resolution
417 fProfV2Resolution[i] = new TProfile(Form("fProfV2Resolution_%i", i), Form("fProfV2Resolution_%i", i), 11, -0.5, 10.5);
418 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(3, "<cos(2(#Psi_{VZEROA} - #Psi_{VZEROC}))>");
419 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(4, "<cos(2(#Psi_{VZEROC} - #Psi_{VZEROA}))>");
420 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(5, "<cos(2(#Psi_{VZEROA} - #Psi_{TPC}))>");
421 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(6, "<cos(2(#Psi_{TPC} - #Psi_{VZEROA}))>");
422 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(7, "<cos(2(#Psi_{VZEROC} - #Psi_{TPC}))>");
423 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(8, "<cos(2(#Psi_{TPC} - #Psi_{VZEROC}))>");
424 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(9, "<cos(2(#Psi_{VZERO} - #Psi_{TPC_A}))>");
425 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(10, "<cos(2(#Psi_{VZERO} - #Psi_{TPC_B}))>");
426 fProfV2Resolution[i]->GetXaxis()->SetBinLabel(11, "<cos(2(#Psi_{TPC_A} - #Psi_{TPC_B}))>");
427 fOutputList->Add(fProfV2Resolution[i]);
428 fProfV3Resolution[i] = new TProfile(Form("fProfV3Resolution_%i", i), Form("fProfV3Resolution_%i", i), 11, -0.5, 10.5);
429 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(3, "<cos(3(#Psi_{VZEROA} - #Psi_{VZEROC}))>");
430 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(4, "<cos(3(#Psi_{VZEROC} - #Psi_{VZEROA}))>");
431 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(5, "<cos(3(#Psi_{VZEROA} - #Psi_{TPC}))>");
432 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(6, "<cos(3(#Psi_{TPC} - #Psi_{VZEROA}))>");
433 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(7, "<cos(3(#Psi_{VZEROC} - #Psi_{TPC}))>");
434 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(8, "<cos(3(#Psi_{TPC} - #Psi_{VZEROC}))>");
435 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(9, "<cos(3(#Psi_{VZERO} - #Psi_{TPC_A}))>");
436 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(10, "<cos(3(#Psi_{VZERO} - #Psi_{TPC_B}))>");
437 fProfV3Resolution[i]->GetXaxis()->SetBinLabel(11, "<cos(3(#Psi_{TPC_A} - #Psi_{TPC_B}))>");
438 fOutputList->Add(fProfV3Resolution[i]);
441 Float_t temp[fCentralityClasses->GetSize()];
442 for(Int_t i(0); i < fCentralityClasses->GetSize(); i++) temp[i] = fCentralityClasses->At(i);
443 fProfV2 = new TProfile("fProfV2", "fProfV2", fCentralityClasses->GetSize()-1, temp);
444 fProfV3 = new TProfile("fProfV3", "fProfV3", fCentralityClasses->GetSize()-1, temp);
445 fOutputList->Add(fProfV2);
446 fOutputList->Add(fProfV3);
447 switch (fFitModulationType) {
449 fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
450 fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
451 fOutputList->Add(fProfV2Cumulant);
452 fOutputList->Add(fProfV3Cumulant);
455 fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
456 fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
457 fOutputList->Add(fProfV2Cumulant);
458 fOutputList->Add(fProfV3Cumulant);
462 // for the histograms initialized below, binning is fixed to runnumbers or flags
463 fReduceBinsXByFactor = 1;
464 fReduceBinsYByFactor = 1;
465 if(fFillQAHistograms) {
466 fHistRunnumbersEta = new TH2F("fHistRunnumbersEta", "fHistRunnumbersEta", fExpectedRuns->GetSize()+1, -.5, fExpectedRuns->GetSize()+.5, 100, -1.1, 1.1);
467 fHistRunnumbersEta->Sumw2();
468 fOutputList->Add(fHistRunnumbersEta);
469 fHistRunnumbersPhi = new TH2F("fHistRunnumbersPhi", "fHistRunnumbersPhi", fExpectedRuns->GetSize()+1, -.5, fExpectedRuns->GetSize()+.5, 100, -0.2, TMath::TwoPi()+0.2);
470 fHistRunnumbersPhi->Sumw2();
471 fOutputList->Add(fHistRunnumbersPhi);
472 for(Int_t i(0); i < fExpectedRuns->GetSize(); i++) {
473 fHistRunnumbersPhi->GetXaxis()->SetBinLabel(i+1, Form("%i", fExpectedRuns->At(i)));
474 fHistRunnumbersEta->GetXaxis()->SetBinLabel(i+1, Form("%i", fExpectedRuns->At(i)));
476 fHistRunnumbersPhi->GetXaxis()->SetBinLabel(fExpectedRuns->GetSize()+1, "undetermined");
477 fHistRunnumbersEta->GetXaxis()->SetBinLabel(fExpectedRuns->GetSize()+1, "undetermined");
479 fHistAnalysisSummary = BookTH1F("fHistAnalysisSummary", "flag", 54, -0.5, 54.5);
480 fHistSwap = new TH1F("fHistSwap", "fHistSwap", 20, 0, TMath::TwoPi());
481 if(fUsePtWeight) fHistSwap->Sumw2();
483 if(fUserSuppliedV2) fOutputList->Add(fUserSuppliedV2);
484 if(fUserSuppliedV3) fOutputList->Add(fUserSuppliedV3);
485 if(fUserSuppliedR2) fOutputList->Add(fUserSuppliedR2);
486 if(fUserSuppliedR3) fOutputList->Add(fUserSuppliedR3);
487 // increase readability of output list
489 // cdf and pdf of chisquare distribution
490 fHistPvalueCDF = BookTH1F("fHistPvalueCDF", "CDF #chi^{2}", 50, 0, 1);
491 fHistPvalueCDFCent = BookTH2F("fHistPvalueCDFCent", "centrality", "p-value", 40, 0, 100, 40, 0, 1);
492 fHistChi2Cent = BookTH2F("fHistChi2Cent", "centrality", "#tilde{#chi^{2}}", 100, 0, 100, 100, 0, 5);
493 fHistPChi2 = BookTH2F("fHistPChi2", "p-value", "#tilde{#chi^{2}}", 1000, 0, 1, 100, 0, 5);
494 fHistKolmogorovTest = BookTH1F("fHistKolmogorovTest", "KolmogorovTest", 50, 0, 1);
495 fHistKolmogorovTestCent = BookTH2F("fHistKolmogorovTestCent", "centrality", "Kolmogorov p", 40, 0, 100, 45, 0, 1);
496 fHistPvalueCDFROOT = BookTH1F("fHistPvalueCDFROOT", "CDF #chi^{2} ROOT", 50, 0, 1);
497 fHistPvalueCDFROOTCent = BookTH2F("fHistPvalueCDFROOTCent", "centrality", "p-value ROOT", 40, 0, 100, 45, 0, 1);
498 fHistChi2ROOTCent = BookTH2F("fHistChi2ROOTCent", "centrality", "#tilde{#chi^{2}}", 40, 0, 100, 45, 0, 5);
499 fHistPChi2Root = BookTH2F("fHistPChi2Root", "p-value", "#tilde{#chi^{2}} ROOT", 1000, 0, 1, 100, 0, 5);
500 fHistPKolmogorov = BookTH2F("fHistPKolmogorov", "p-value", "kolmogorov p",40, 0, 1, 40, 0, 1);
501 fHistRhoStatusCent = BookTH2F("fHistRhoStatusCent", "centrality", "status [-1=lin was better, 0=ok, 1 = failed]", 101, -1, 100, 3, -1.5, 1.5);
502 fHistUndeterminedRunQA = BookTH1F("fHistUndeterminedRunQA", "runnumber", 10, 0, 10);
504 PostData(1, fOutputList);
506 switch (fRunModeType) {
508 fOutputListGood = new TList();
509 fOutputListGood->SetOwner(kTRUE);
510 fOutputListBad = new TList();
511 fOutputListBad->SetOwner(kTRUE);
512 PostData(2, fOutputListGood);
513 PostData(3, fOutputListBad);
518 // get the containers
519 fTracksCont = GetParticleContainer("Tracks");
520 fJetsCont = GetJetContainer("Jets");
522 //_____________________________________________________________________________
523 Bool_t AliAnalysisTaskRhoVnModulation::Run()
525 // user exec: execute once for each event
526 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
527 if(!fTracks||!fJets||!fRho) return kFALSE;
528 if(!fLocalInit) fLocalInit = InitializeAnalysis();
529 // reject the event if expected data is missing
530 if(!PassesCuts(InputEvent())) return kFALSE;
531 fLeadingJet = GetLeadingJet(); // store the leading jet
533 fLocalRho->SetVal(fRho->GetVal());
534 // [0][0] psi2a [1,0] psi2c
535 // [0][1] psi3a [1,1] psi3c
536 Double_t vzero[2][2];
537 CalculateEventPlaneVZERO(vzero);
538 /* for the combined vzero event plane
540 * not fully implmemented yet, use with caution ! */
541 Double_t vzeroComb[2];
542 CalculateEventPlaneCombinedVZERO(vzeroComb);
545 CalculateEventPlaneTPC(tpc);
546 Double_t psi2(-1), psi3(-1);
547 // arrays which will hold the fit parameters
548 switch (fDetectorType) { // determine the detector type for the rho fit
549 case kTPC : { psi2 = tpc[0]; psi3 = tpc[1]; } break;
550 case kVZEROA : { psi2 = vzero[0][0]; psi3 = vzero[0][1]; } break;
551 case kVZEROC : { psi2 = vzero[1][0]; psi3 = vzero[1][1]; } break;
552 case kVZEROComb : { psi2 = vzeroComb[0]; psi3 = vzeroComb[1];} break;
555 switch (fFitModulationType) { // do the fits
557 switch (fCollisionType) {
558 case kPythia : { // background is zero for pp jets
559 fFitModulation->FixParameter(0, 0);
560 fLocalRho->SetVal(0);
563 fFitModulation->FixParameter(0, fLocalRho->GetVal());
567 case kV2 : { // only v2
568 if(CorrectRho(psi2, psi3)) {
569 fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
570 if(fUserSuppliedR2) {
571 Double_t r(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
572 if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
574 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
577 case kV3 : { // only v3
578 if(CorrectRho(psi2, psi3)) {
579 if(fUserSuppliedR3) {
580 Double_t r(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
581 if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
583 fProfV3->Fill(fCent, fFitModulation->GetParameter(3));
584 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
587 case kQC2 : { // qc2 analysis
588 if(CorrectRho(psi2, psi3)) {
589 if(fUserSuppliedR2 && fUserSuppliedR3) {
590 // note for the qc method, resolution is REVERSED to go back to v2obs
591 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
592 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
593 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
594 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)*r3);
596 if (fUsePtWeight) { // use weighted weights
597 Double_t dQCnM11 = (fNoEventWeightsForQC) ? 1. : QCnM11();
598 fProfV2->Fill(fCent, fFitModulation->GetParameter(3), dQCnM11);
599 fProfV3->Fill(fCent, fFitModulation->GetParameter(7), dQCnM11);
601 Double_t dQCnM = (fNoEventWeightsForQC) ? 2. : QCnM();
602 fProfV2->Fill(fCent, fFitModulation->GetParameter(3), dQCnM*(dQCnM-1));
603 fProfV3->Fill(fCent, fFitModulation->GetParameter(7), dQCnM*(dQCnM-1));
605 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
609 if(CorrectRho(psi2, psi3)) {
610 if(fUserSuppliedR2 && fUserSuppliedR3) {
611 // note for the qc method, resolution is REVERSED to go back to v2obs
612 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
613 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
614 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
615 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)*r3);
617 if (fUsePtWeight) { // use weighted weights
618 fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM1111()*/);
619 fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM1111()*/);
621 fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
622 fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
625 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
628 if(CorrectRho(psi2, psi3)) {
629 if(fUserSuppliedR2 && fUserSuppliedR3) {
630 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
631 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
632 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r2);
633 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)/r3);
635 fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
636 fProfV3->Fill(fCent, fFitModulation->GetParameter(7));
637 CalculateEventPlaneResolution(vzero, vzeroComb, tpc);
641 // if all went well, update the local rho parameter
642 fLocalRho->SetLocalRho(fFitModulation);
643 // fill a number of histograms. event qa needs to be filled first as it also determines the runnumber for the track qa
644 if(fFillQAHistograms) FillQAHistograms(InputEvent());
645 if(fFillHistograms) FillHistogramsAfterSubtraction(psi2, psi3, vzero, vzeroComb, tpc);
646 // send the output to the connected output container
647 PostData(1, fOutputList);
648 switch (fRunModeType) {
650 PostData(2, fOutputListGood);
651 PostData(3, fOutputListBad);
658 //_____________________________________________________________________________
659 void AliAnalysisTaskRhoVnModulation::CalculateEventPlaneVZERO(Double_t vzero[2][2]) const
661 // get the vzero event plane
662 if(fUseV0EventPlaneFromHeader) { // use the vzero from the header
663 Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0), h(0);
664 vzero[0][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 2, a, b);
665 vzero[1][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 2, c, d);
666 vzero[0][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 3, e, f);
667 vzero[1][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 3, g, h);
670 // grab the vzero event plane without recentering
671 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
672 Double_t qxa2(0), qya2(0), qxc2(0), qyc2(0); // for psi2
673 Double_t qxa3(0), qya3(0), qxc3(0), qyc3(0); // for psi3
674 for(Int_t iVZERO(0); iVZERO < 64; iVZERO++) {
675 Double_t phi(TMath::PiOver4()*(.5+iVZERO%8)), /* eta(0), */ weight(InputEvent()->GetVZEROEqMultiplicity(iVZERO));
676 // (iVZERO<32) ? eta = -3.45+.5*(iVZERO/8) : eta = 4.8-.6*((iVZERO/8)-4);
678 qxa2 += weight*TMath::Cos(2.*phi);
679 qya2 += weight*TMath::Sin(2.*phi);
680 qxa3 += weight*TMath::Cos(3.*phi);
681 qya3 += weight*TMath::Sin(3.*phi);
684 qxc2 += weight*TMath::Cos(2.*phi);
685 qyc2 += weight*TMath::Sin(2.*phi);
686 qxc3 += weight*TMath::Cos(3.*phi);
687 qyc3 += weight*TMath::Sin(3.*phi);
690 vzero[0][0] = .5*TMath::ATan2(qya2, qxa2);
691 vzero[1][0] = .5*TMath::ATan2(qyc2, qxc2);
692 vzero[0][1] = (1./3.)*TMath::ATan2(qya3, qxa3);
693 vzero[1][1] = (1./3.)*TMath::ATan2(qyc3, qxc3);
695 //_____________________________________________________________________________
696 void AliAnalysisTaskRhoVnModulation::CalculateEventPlaneTPC(Double_t* tpc)
698 // grab the TPC event plane
699 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
700 fNAcceptedTracks = 0; // reset the track counter
701 Double_t qx2(0), qy2(0); // for psi2
702 Double_t qx3(0), qy3(0); // for psi3
704 Float_t excludeInEta = -999;
705 if(fExcludeLeadingJetsFromFit > 0 ) { // remove the leading jet from ep estimate
706 if(fLeadingJet) excludeInEta = fLeadingJet->Eta();
708 Int_t iTracks(fTracks->GetEntriesFast());
709 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
710 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
711 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
712 if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta) < GetJetContainer()->GetJetRadius()*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - GetJetContainer()->GetJetRadius() - GetJetContainer()->GetJetEtaMax() ) > 0 )) continue;
714 qx2+= TMath::Cos(2.*track->Phi());
715 qy2+= TMath::Sin(2.*track->Phi());
716 qx3+= TMath::Cos(3.*track->Phi());
717 qy3+= TMath::Sin(3.*track->Phi());
720 tpc[0] = .5*TMath::ATan2(qy2, qx2);
721 tpc[1] = (1./3.)*TMath::ATan2(qy3, qx3);
723 //_____________________________________________________________________________
724 void AliAnalysisTaskRhoVnModulation::CalculateEventPlaneCombinedVZERO(Double_t* comb) const
726 // grab the combined vzero event plane
727 // if(fUseV0EventPlaneFromHeader) { // use the vzero from the header
728 Double_t a(0), b(0), c(0), d(0);
729 comb[0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 10, 2, a, b);
730 comb[1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 10, 3, c, d);
732 // Double_t qx2a(0), qy2a(0), qx2c(0), qy2c(0), qx3a(0), qy3a(0), qx3c(0), qy3c(0);
733 // InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 2, qx2a, qy2a);
734 // InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 2, qx2c, qy2c);
735 // InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 3, qx3a, qy3a);
736 // InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 3, qx3c, qy3c);
737 // FIXME the rest of this function isn't impelmented yet (as of 01-07-2013)
738 // Double_t chi2A(-1), chi2C(-1), chi3A(-1), chi3C(-1); // get chi from the resolution
739 // Double_t qx2(chi2A*chi2A*qx2a+chi2C*chi2C*qx2c);
740 // Double_t qy2(chi2A*chi2A*qy2a+chi2C*chi2C*qy2c);
741 // Double_t qx3(chi3A*chi3A*qx3a+chi3C*chi3C*qx3c);
742 // Double_t qy3(chi3A*chi3A*qy3a+chi3C*chi3C*qy3c);
743 // comb[0] = .5*TMath::ATan2(qy2, qx2);
744 // comb[1] = (1./3.)*TMath::ATan2(qy3, qx3);
747 //_____________________________________________________________________________
748 void AliAnalysisTaskRhoVnModulation::CalculateEventPlaneResolution(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc)
750 // fill the profiles for the resolution parameters
751 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
752 fProfV2Resolution[fInCentralitySelection]->Fill(2., TMath::Cos(2.*(vzero[0][0] - vzero[1][0])));
753 fProfV2Resolution[fInCentralitySelection]->Fill(3., TMath::Cos(2.*(vzero[1][0] - vzero[0][0])));
754 fProfV2Resolution[fInCentralitySelection]->Fill(4., TMath::Cos(2.*(vzero[0][0] - tpc[0])));
755 fProfV2Resolution[fInCentralitySelection]->Fill(5., TMath::Cos(2.*(tpc[0] - vzero[0][0])));
756 fProfV2Resolution[fInCentralitySelection]->Fill(6., TMath::Cos(2.*(vzero[1][0] - tpc[0])));
757 fProfV2Resolution[fInCentralitySelection]->Fill(7., TMath::Cos(2.*(tpc[0] - vzero[1][0])));
758 fProfV3Resolution[fInCentralitySelection]->Fill(2., TMath::Cos(3.*(vzero[0][0] - vzero[1][0])));
759 fProfV3Resolution[fInCentralitySelection]->Fill(3., TMath::Cos(3.*(vzero[1][0] - vzero[0][0])));
760 fProfV3Resolution[fInCentralitySelection]->Fill(4., TMath::Cos(3.*(vzero[0][0] - tpc[0])));
761 fProfV3Resolution[fInCentralitySelection]->Fill(5., TMath::Cos(3.*(tpc[0] - vzero[0][0])));
762 fProfV3Resolution[fInCentralitySelection]->Fill(6., TMath::Cos(3.*(vzero[1][0] - tpc[0])));
763 fProfV3Resolution[fInCentralitySelection]->Fill(7., TMath::Cos(3.*(tpc[0] - vzero[1][0])));
764 // for the resolution of the combined vzero event plane, use two tpc halves as uncorrelated subdetectors
765 Double_t qx2a(0), qy2a(0); // for psi2a, negative eta
766 Double_t qx3a(0), qy3a(0); // for psi3a, negative eta
767 Double_t qx2b(0), qy2b(0); // for psi2a, positive eta
768 Double_t qx3b(0), qy3b(0); // for psi3a, positive eta
770 Int_t iTracks(fTracks->GetEntriesFast());
771 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
772 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
773 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
774 if(track->Eta() < 0 ) {
775 qx2a+= TMath::Cos(2.*track->Phi());
776 qy2a+= TMath::Sin(2.*track->Phi());
777 qx3a+= TMath::Cos(3.*track->Phi());
778 qy3a+= TMath::Sin(3.*track->Phi());
779 } else if (track->Eta() > 0) {
780 qx2b+= TMath::Cos(2.*track->Phi());
781 qy2b+= TMath::Sin(2.*track->Phi());
782 qx3b+= TMath::Cos(3.*track->Phi());
783 qy3b+= TMath::Sin(3.*track->Phi());
787 Double_t tpca2(.5*TMath::ATan2(qy2a, qx2a));
788 Double_t tpca3((1./3.)*TMath::ATan2(qy3a, qx3a));
789 Double_t tpcb2(.5*TMath::ATan2(qy2b, qx2b));
790 Double_t tpcb3((1./3.)*TMath::ATan2(qy3b, qx3b));
791 fProfV2Resolution[fInCentralitySelection]->Fill(8., TMath::Cos(2.*(vzeroComb[0] - tpca2)));
792 fProfV2Resolution[fInCentralitySelection]->Fill(9., TMath::Cos(2.*(vzeroComb[0] - tpcb2)));
793 fProfV2Resolution[fInCentralitySelection]->Fill(10., TMath::Cos(2.*(tpca2 - tpcb2)));
794 fProfV3Resolution[fInCentralitySelection]->Fill(8., TMath::Cos(3.*(vzeroComb[1] - tpca3)));
795 fProfV3Resolution[fInCentralitySelection]->Fill(9., TMath::Cos(3.*(vzeroComb[1] - tpcb3)));
796 fProfV3Resolution[fInCentralitySelection]->Fill(10., TMath::Cos(3.*(tpca3 - tpcb3)));
798 //_____________________________________________________________________________
799 Double_t AliAnalysisTaskRhoVnModulation::CalculateEventPlaneChi(Double_t resEP) const
801 // Get Chi from EP resolution (PRC 58 1671)
802 Double_t chi(2.), delta (1.);
803 for (Int_t i(0); i < 15; i++) {
804 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;
809 //_____________________________________________________________________________
810 void AliAnalysisTaskRhoVnModulation::CalculateRandomCone(Float_t &pt, Float_t &eta, Float_t &phi,
811 AliEmcalJet* jet) const
814 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
815 pt = 0; eta = 0; phi = 0;
816 Float_t etaJet(999), phiJet(999), dJet(999); // no jet: same as jet very far away
817 if(jet) { // if a leading jet is given, use its kinematic properties
821 // the random cone acceptance has to equal the jet acceptance
822 // this also insures safety when runnnig on the semi-good tpc runs for 11h data,
823 // where jet acceptance is adjusted to reduced acceptance - hence random cone acceptance as well
824 Float_t minPhi(GetJetContainer()->GetJetPhiMin()), maxPhi(GetJetContainer()->GetJetPhiMax());
825 if(maxPhi > TMath::TwoPi()) maxPhi = TMath::TwoPi();
826 if(minPhi < 0 ) minPhi = 0;
827 Float_t diffRcRJR(TMath::Abs(fRandomConeRadius-GetJetContainer()->GetJetRadius()));
828 // construct a random cone and see if it's far away enough from the leading jet
829 Int_t attempts(1000);
832 eta = gRandom->Uniform(GetJetContainer()->GetJetEtaMin()+diffRcRJR, GetJetContainer()->GetJetEtaMax()-diffRcRJR);
833 phi = gRandom->Uniform(minPhi, maxPhi);
835 dJet = TMath::Sqrt((etaJet-eta)*(etaJet-eta)+(phiJet-phi)*(phiJet-phi));
836 if(dJet > fMinDisanceRCtoLJ) break;
837 else if (attempts == 0) {
838 printf(" > No random cone after 1000 tries, giving up ... !\n");
843 AliVParticle* track = fTracksCont->GetNextAcceptParticle(0);
845 Float_t etaTrack(track->Eta()), phiTrack(track->Phi()), ptTrack(track->Pt());
846 // get distance from cone
847 if(TMath::Abs(phiTrack-phi) > TMath::Abs(phiTrack - phi + TMath::TwoPi())) phiTrack+=TMath::TwoPi();
848 if(TMath::Abs(phiTrack-phi) > TMath::Abs(phiTrack - phi - TMath::TwoPi())) phiTrack-=TMath::TwoPi();
849 if(TMath::Sqrt(TMath::Abs((etaTrack-eta)*(etaTrack-eta)+(phiTrack-phi)*(phiTrack-phi))) <= fRandomConeRadius) pt+=ptTrack;
850 track = fTracksCont->GetNextAcceptParticle();
854 //_____________________________________________________________________________
855 Double_t AliAnalysisTaskRhoVnModulation::CalculateQC2(Int_t harm) {
856 // get the second order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho
857 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
858 Double_t reQ(0), imQ(0), modQ(0), M11(0), M(0);
859 if(fUsePtWeight) { // for the weighted 2-nd order q-cumulant
860 QCnQnk(harm, 1, reQ, imQ); // get the weighted 2-nd order q-vectors
861 modQ = reQ*reQ+imQ*imQ; // get abs Q-squared
862 M11 = QCnM11(); // equals S2,1 - S1,2
863 return (M11 > 0) ? ((modQ - QCnS(1,2))/M11) : -999;
864 } // else return the non-weighted 2-nd order q-cumulant
865 QCnQnk(harm, 0, reQ, imQ); // get the non-weighted 2-nd order q-vectors
866 modQ = reQ*reQ+imQ*imQ; // get abs Q-squared
868 return (M > 1) ? (modQ - M)/(M*(M-1)) : -999;
870 //_____________________________________________________________________________
871 Double_t AliAnalysisTaskRhoVnModulation::CalculateQC4(Int_t harm) {
872 // get the fourth order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho
873 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
874 Double_t reQn1(0), imQn1(0), reQ2n2(0), imQ2n2(0), reQn3(0), imQn3(0), M1111(0), M(0);
875 Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0); // terms of the calculation
876 if(fUsePtWeight) { // for the weighted 4-th order q-cumulant
877 QCnQnk(harm, 1, reQn1, imQn1);
878 QCnQnk(harm*2, 2, reQ2n2, imQ2n2);
879 QCnQnk(harm, 3, reQn3, imQn3);
880 // fill in the terms ...
881 a = (reQn1*reQn1+imQn1*imQn1)*(reQn1*reQn1+imQn1*imQn1);
882 b = reQ2n2*reQ2n2 + imQ2n2*imQ2n2;
883 c = -2.*(reQ2n2*reQn1*reQn1-reQ2n2*imQn1*imQn1+2.*imQ2n2*reQn1*imQn1);
884 d = 8.*(reQn3*reQn1+imQn3*imQn1);
885 e = -4.*QCnS(1,2)*(reQn1*reQn1+imQn1*imQn1);
889 return (M1111 > 0) ? (a+b+c+d+e+f+g)/M1111 : -999;
890 } // else return the unweighted case
891 Double_t reQn(0), imQn(0), reQ2n(0), imQ2n(0);
892 QCnQnk(harm, 0, reQn, imQn);
893 QCnQnk(harm*2, 0, reQ2n, imQ2n);
894 // fill in the terms ...
896 if(M < 4) return -999;
897 a = (reQn*reQn+imQn*imQn)*(reQn*reQn+imQn*imQn);
898 b = reQ2n*reQ2n + imQ2n*imQ2n;
899 c = -2.*(reQ2n*reQn*reQn-reQ2n*imQn*imQn+2.*imQ2n*reQn*imQn);
900 e = -4.*(M-2)*(reQn*reQn+imQn*imQn);
902 return (a+b+c+e+f)/(M*(M-1)*(M-2)*(M-3));
904 //_____________________________________________________________________________
905 void AliAnalysisTaskRhoVnModulation::QCnQnk(Int_t n, Int_t k, Double_t &reQ, Double_t &imQ) {
906 // get the weighted n-th order q-vector, pass real and imaginary part as reference
907 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
909 fNAcceptedTracksQCn = 0;
910 Int_t iTracks(fTracks->GetEntriesFast());
911 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
912 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
913 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
914 fNAcceptedTracksQCn++;
915 // for the unweighted case, k equals zero and the weight doesn't contribute to the equation below
916 reQ += TMath::Power(track->Pt(), k) * TMath::Cos(((double)n)*track->Phi());
917 imQ += TMath::Power(track->Pt(), k) * TMath::Sin(((double)n)*track->Phi());
920 //_____________________________________________________________________________
921 void AliAnalysisTaskRhoVnModulation::QCnDiffentialFlowVectors(
922 TClonesArray* pois, TArrayD* ptBins, Bool_t vpart, Double_t* repn, Double_t* impn,
923 Double_t *mp, Double_t *reqn, Double_t *imqn, Double_t* mq, Int_t n)
925 // get unweighted differential flow vectors
926 Int_t iPois(pois->GetEntriesFast());
928 for(Int_t i(0); i < iPois; i++) {
929 for(Int_t ptBin(0); ptBin < ptBins->GetSize()-1; ptBin++) {
930 AliVTrack* poi = static_cast<AliVTrack*>(pois->At(i));
931 if(PassesCuts(poi)) {
932 if(poi->Pt() >= ptBins->At(ptBin) && poi->Pt() < ptBins->At(ptBin+1)) {
933 // fill the flow vectors assuming that all poi's are in the rp selection (true by design)
934 repn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
935 impn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
937 reqn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
938 imqn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
945 for(Int_t i(0); i < iPois; i++) {
946 for(Int_t ptBin(0); ptBin < ptBins->GetSize()-1; ptBin++) {
947 AliEmcalJet* poi = static_cast<AliEmcalJet*>(pois->At(i));
948 if(PassesCuts(poi)) {
949 Double_t pt(poi->Pt()-poi->Area()*fLocalRho->GetLocalVal(poi->Phi(), GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
950 if(pt >= ptBins->At(ptBin) && pt < ptBins->At(ptBin+1)) {
951 repn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
952 impn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
953 mp[ptBin]++; // qn isn't filled, no overlap between poi's and rp's
960 //_____________________________________________________________________________
961 Double_t AliAnalysisTaskRhoVnModulation::QCnS(Int_t i, Int_t j) {
962 // get the weighted ij-th order autocorrelation correction
963 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
964 if(!fTracks || i <= 0 || j <= 0) return -999;
965 Int_t iTracks(fTracks->GetEntriesFast());
967 for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
968 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
969 if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
970 Sij+=TMath::Power(track->Pt(), j);
972 return TMath::Power(Sij, i);
974 //_____________________________________________________________________________
975 Double_t AliAnalysisTaskRhoVnModulation::QCnM() {
976 // get multiplicity for unweighted q-cumulants. function QCnQnk should be called first
977 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
978 return (Double_t) fNAcceptedTracksQCn;
980 //_____________________________________________________________________________
981 Double_t AliAnalysisTaskRhoVnModulation::QCnM11() {
982 // get multiplicity weights for the weighted two particle cumulant
983 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
984 return (QCnS(2,1) - QCnS(1,2));
986 //_____________________________________________________________________________
987 Double_t AliAnalysisTaskRhoVnModulation::QCnM1111() {
988 // get multiplicity weights for the weighted four particle cumulant
989 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
990 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));
992 //_____________________________________________________________________________
993 Bool_t AliAnalysisTaskRhoVnModulation::QCnRecovery(Double_t psi2, Double_t psi3) {
994 // decides how to deal with the situation where c2 or c3 is negative
995 // returns kTRUE depending on whether or not a modulated rho is used for the jet background
996 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
997 if(TMath::AreEqualAbs(fFitModulation->GetParameter(3), .0, 1e-10) && TMath::AreEqualAbs(fFitModulation->GetParameter(7), .0,1e-10)) {
998 fFitModulation->SetParameter(7, 0);
999 fFitModulation->SetParameter(3, 0);
1000 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1001 return kTRUE; // v2 and v3 have physical null values
1003 switch (fQCRecovery) {
1004 case kFixedRho : { // roll back to the original rho
1005 fFitModulation->SetParameter(7, 0);
1006 fFitModulation->SetParameter(3, 0);
1007 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1008 return kFALSE; // rho is forced to be fixed
1010 case kNegativeVn : {
1011 Double_t c2(fFitModulation->GetParameter(3));
1012 Double_t c3(fFitModulation->GetParameter(7));
1013 if( c2 < 0 ) c2 = -1.*TMath::Sqrt(-1.*c2);
1014 if( c3 < 0 ) c3 = -1.*TMath::Sqrt(-1.*c3);
1015 fFitModulation->SetParameter(3, c2);
1016 fFitModulation->SetParameter(7, c3);
1017 return kTRUE; // is this a physical quantity ?
1020 fitModulationType tempType(fFitModulationType); // store temporarily
1021 fFitModulationType = kCombined;
1022 fFitModulation->SetParameter(7, 0);
1023 fFitModulation->SetParameter(3, 0);
1024 Bool_t pass(CorrectRho(psi2, psi3)); // do the fit and all quality checks
1025 fFitModulationType = tempType; // roll back for next event
1028 default : return kFALSE;
1032 //_____________________________________________________________________________
1033 Bool_t AliAnalysisTaskRhoVnModulation::CorrectRho(Double_t psi2, Double_t psi3)
1035 // get rho' -> rho(phi)
1036 // two routines are available, both can be used with or without pt weights
1037 // [1] get vn from q-cumulants or as an integrated value from a user supplied histogram
1038 // in case of cumulants, both cumulants and vn values are stored. in both cases, v2 and v3
1039 // are expected. a check is performed to see if rho has no negative local minimum
1040 // for full description, see Phys. Rev. C 83, 044913
1041 // since the cn distribution has negative values, vn = sqrt(cn) can be imaginary sometimes
1042 // in this case one can either roll back to the 'original' rixed rho, do a fit for vn or take use
1043 // vn = - sqrt(|cn|)
1044 // [2] fitting a fourier expansion to the de/dphi distribution
1045 // the fit can be done with either v2, v3 or a combination.
1046 // in all cases, a cut can be made on the p-value of the chi-squared value of the fit
1047 // and a check can be performed to see if rho has no negative local minimum
1048 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1049 Int_t freeParams(2); // free parameters of the fit (for NDF)
1050 switch (fFitModulationType) { // for approaches where no fitting is required
1052 fFitModulation->FixParameter(4, psi2);
1053 fFitModulation->FixParameter(6, psi3);
1054 fFitModulation->FixParameter(3, CalculateQC2(2)); // set here with cn, vn = sqrt(cn)
1055 fFitModulation->FixParameter(7, CalculateQC2(3));
1056 // first fill the histos of the raw cumulant distribution
1057 if (fUsePtWeight) { // use weighted weights
1058 Double_t dQCnM11 = (fNoEventWeightsForQC) ? 1. : QCnM11();
1059 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), dQCnM11);
1060 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), dQCnM11);
1062 Double_t dQCnM = (fNoEventWeightsForQC) ? 2. : QCnM();
1063 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), dQCnM*(dQCnM-1));
1064 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), dQCnM*(dQCnM-1));
1066 // then see if one of the cn value is larger than zero and vn is readily available
1067 if(fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) {
1068 fFitModulation->FixParameter(3, TMath::Sqrt(fFitModulation->GetParameter(3)));
1069 fFitModulation->FixParameter(7, TMath::Sqrt(fFitModulation->GetParameter(7)));
1070 } else if (!QCnRecovery(psi2, psi3)) return kFALSE; // try to recover the cumulant, this will set v2 and v3
1071 if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { // general check
1072 fFitModulation->SetParameter(7, 0);
1073 fFitModulation->SetParameter(3, 0);
1074 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1080 fFitModulation->FixParameter(4, psi2);
1081 fFitModulation->FixParameter(6, psi3);
1082 fFitModulation->FixParameter(3, CalculateQC4(2)); // set here with cn, vn = sqrt(cn)
1083 fFitModulation->FixParameter(7, CalculateQC4(3));
1084 // first fill the histos of the raw cumulant distribution
1085 if (fUsePtWeight) { // use weighted weights
1086 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
1087 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
1089 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
1090 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
1092 // then see if one of the cn value is larger than zero and vn is readily available
1093 if(fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) {
1094 fFitModulation->FixParameter(3, TMath::Sqrt(fFitModulation->GetParameter(3)));
1095 fFitModulation->FixParameter(7, TMath::Sqrt(fFitModulation->GetParameter(7)));
1096 } else if (!QCnRecovery(psi2, psi3)) return kFALSE; // try to recover the cumulant, this will set v2 and v3
1097 if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { // general check
1098 fFitModulation->SetParameter(7, 0);
1099 fFitModulation->SetParameter(3, 0);
1100 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1104 case kIntegratedFlow : {
1105 // use v2 and v3 values from an earlier iteration over the data
1106 fFitModulation->FixParameter(3, fUserSuppliedV2->GetBinContent(fUserSuppliedV2->GetXaxis()->FindBin(fCent)));
1107 fFitModulation->FixParameter(4, psi2);
1108 fFitModulation->FixParameter(6, psi3);
1109 fFitModulation->FixParameter(7, fUserSuppliedV3->GetBinContent(fUserSuppliedV3->GetXaxis()->FindBin(fCent)));
1110 if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) {
1111 fFitModulation->SetParameter(7, 0);
1112 fFitModulation->SetParameter(3, 0);
1113 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1120 TString detector("");
1121 switch (fDetectorType) {
1122 case kTPC : detector+="TPC";
1124 case kVZEROA : detector+="VZEROA";
1126 case kVZEROC : detector+="VZEROC";
1128 case kVZEROComb : detector+="VZEROComb";
1132 Int_t iTracks(fTracks->GetEntriesFast());
1133 Double_t excludeInEta = -999;
1134 Double_t excludeInPhi = -999;
1135 Double_t excludeInPt = -999;
1136 if(iTracks <= 0 || fLocalRho->GetVal() <= 0 ) return kFALSE; // no use fitting an empty event ...
1137 if(fExcludeLeadingJetsFromFit > 0 ) {
1139 excludeInEta = fLeadingJet->Eta();
1140 excludeInPhi = fLeadingJet->Phi();
1141 excludeInPt = fLeadingJet->Pt();
1144 // check the acceptance of the track selection that will be used
1145 // if one uses e.g. semi-good tpc tracks, accepance in phi is reduced to 0 < phi < 4
1146 // the defaults (-10 < phi < 10) which accept all, are then overwritten
1147 Double_t lowBound(0.), upBound(TMath::TwoPi()); // bounds for fit
1148 if(GetParticleContainer()->GetParticlePhiMin() > lowBound) lowBound = GetParticleContainer()->GetParticlePhiMin();
1149 if(GetParticleContainer()->GetParticlePhiMax() < upBound) upBound = GetParticleContainer()->GetParticlePhiMax();
1151 fHistSwap->Reset(); // clear the histogram
1152 TH1F _tempSwap; // on stack for quick access
1153 TH1F _tempSwapN; // on stack for quick access, bookkeeping histogram
1154 if(fRebinSwapHistoOnTheFly) {
1155 if(fNAcceptedTracks < 49) fNAcceptedTracks = 49; // avoid aliasing effects
1156 _tempSwap = TH1F("_tempSwap", "_tempSwap", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), lowBound, upBound);
1157 if(fUsePtWeightErrorPropagation) _tempSwapN = TH1F("_tempSwapN", "_tempSwapN", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), lowBound, upBound);
1158 if(fUsePtWeight) _tempSwap.Sumw2();
1160 else _tempSwap = *fHistSwap; // now _tempSwap holds the desired histo
1161 // non poissonian error when using pt weights
1162 Double_t totalpts(0.), totalptsquares(0.), totalns(0.);
1163 for(Int_t i(0); i < iTracks; i++) {
1164 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1165 if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta) < GetJetContainer()->GetJetRadius()*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - GetJetContainer()->GetJetRadius() - GetJetContainer()->GetJetEtaMax() ) > 0 )) continue;
1166 if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
1168 _tempSwap.Fill(track->Phi(), track->Pt());
1169 if(fUsePtWeightErrorPropagation) {
1170 totalpts += track->Pt();
1171 totalptsquares += track->Pt()*track->Pt();
1173 _tempSwapN.Fill(track->Phi());
1176 else _tempSwap.Fill(track->Phi());
1178 if(fUsePtWeight && fUsePtWeightErrorPropagation) {
1179 // in the case of pt weights overwrite the poissonian error estimate which is assigned by root by a more sophisticated appraoch
1180 // 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
1181 // 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
1182 // event, for the latter use a poissonian estimate. the two contrubitions are assumed to be uncorrelated
1183 if(totalns < 1) return kFALSE; // not one track passes the cuts
1184 for(Int_t l = 0; l < _tempSwap.GetNbinsX(); l++) {
1185 if(_tempSwapN.GetBinContent(l+1) == 0) {
1186 _tempSwap.SetBinContent(l+1,0);
1187 _tempSwap.SetBinError(l+1,0);
1190 Double_t vartimesnsq = totalptsquares*totalns - totalpts*totalpts;
1191 Double_t variance = vartimesnsq/(totalns*(totalns-1.));
1192 Double_t SDOMSq = variance / _tempSwapN.GetBinContent(l+1);
1193 Double_t SDOMSqOverMeanSq = SDOMSq * _tempSwapN.GetBinContent(l+1) * _tempSwapN.GetBinContent(l+1) / (_tempSwapN.GetBinContent(l+1) * _tempSwapN.GetBinContent(l+1));
1194 Double_t poissonfrac = 1./_tempSwapN.GetBinContent(l+1);
1195 Double_t vartotalfrac = SDOMSqOverMeanSq + poissonfrac;
1196 Double_t vartotal = vartotalfrac * _tempSwap.GetBinContent(l+1) * _tempSwap.GetBinContent(l+1);
1197 if(vartotal > 0.0001) _tempSwap.SetBinError(l+1,TMath::Sqrt(vartotal));
1199 _tempSwap.SetBinContent(l+1,0);
1200 _tempSwap.SetBinError(l+1,0);
1206 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1207 switch (fFitModulationType) {
1209 fFitModulation->FixParameter(0, fLocalRho->GetVal() );
1213 fFitModulation->FixParameter(4, psi2);
1217 fFitModulation->FixParameter(4, psi3);
1221 fFitModulation->FixParameter(4, psi2);
1222 fFitModulation->FixParameter(6, psi3);
1225 case kFourierSeries : {
1226 // in this approach, an explicit calculation will be made of vn = sqrt(xn^2+yn^2)
1227 // where x[y] = Integrate[r(phi)cos[sin](n phi)dphi, 0, 2pi]
1228 Double_t cos2(0), sin2(0), cos3(0), sin3(0), sumPt(0);
1229 for(Int_t i(0); i < iTracks; i++) {
1230 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1231 if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
1232 sumPt += track->Pt();
1233 cos2 += track->Pt()*TMath::Cos(2*PhaseShift(track->Phi()-psi2));
1234 sin2 += track->Pt()*TMath::Sin(2*PhaseShift(track->Phi()-psi2));
1235 cos3 += track->Pt()*TMath::Cos(3*PhaseShift(track->Phi()-psi3));
1236 sin3 += track->Pt()*TMath::Sin(3*PhaseShift(track->Phi()-psi3));
1238 fFitModulation->SetParameter(3, TMath::Sqrt(cos2*cos2+sin2*sin2)/fLocalRho->GetVal());
1239 fFitModulation->SetParameter(4, psi2);
1240 fFitModulation->SetParameter(6, psi3);
1241 fFitModulation->SetParameter(7, TMath::Sqrt(cos3*cos3+sin3*sin3)/fLocalRho->GetVal());
1246 // toy mc, just here to check procedure, azimuthal profile is filled from hypothesis so p-value distribution should be flat
1247 Int_t _bins = _tempSwap.GetXaxis()->GetNbins();
1248 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());
1249 _tempFit->SetParameter(0, fFitModulation->GetParameter(0)); // normalization
1250 _tempFit->SetParameter(3, 0.1); // v2
1251 _tempFit->FixParameter(1, 1.); // constant
1252 _tempFit->FixParameter(2, 2.); // constant
1253 _tempFit->FixParameter(5, 3.); // constant
1254 _tempFit->FixParameter(4, fFitModulation->GetParameter(4));
1255 _tempFit->FixParameter(6, fFitModulation->GetParameter(6));
1256 _tempFit->SetParameter(7, 0.1); // v3
1257 _tempSwap.Reset(); // rese bin content
1258 for(int _binsI = 0; _binsI < _bins*_bins; _binsI++) _tempSwap.Fill(_tempFit->GetRandom());
1260 _tempSwap.Fit(fFitModulation, fFitModulationOptions.Data(), "", lowBound, upBound);
1261 // the quality of the fit is evaluated from 1 - the cdf of the chi square distribution
1262 // three methods are available, all with their drawbacks. all are stored, one is selected to do the cut
1263 Int_t NDF(_tempSwap.GetXaxis()->GetNbins()-freeParams);
1264 if(NDF == 0) return kFALSE;
1265 Double_t CDF(1.-ChiSquareCDF(NDF, ChiSquare(_tempSwap, fFitModulation)));
1266 Double_t CDFROOT(1.-ChiSquareCDF(NDF, fFitModulation->GetChisquare()));
1267 Double_t CDFKolmogorov(KolmogorovTest(_tempSwap, fFitModulation));
1268 // fill the values and centrality correlation (redundant but easy on the eyes)
1269 fHistPvalueCDF->Fill(CDF);
1270 fHistPvalueCDFCent->Fill(fCent, CDF);
1271 fHistPvalueCDFROOT->Fill(CDFROOT);
1272 fHistPvalueCDFROOTCent->Fill(fCent, CDFROOT);
1273 fHistKolmogorovTest->Fill(CDFKolmogorov);
1274 fHistChi2ROOTCent->Fill(fCent, fFitModulation->GetChisquare()/((float)NDF));
1275 fHistChi2Cent->Fill(fCent, ChiSquare(_tempSwap, fFitModulation)/((float)NDF));
1276 fHistKolmogorovTestCent->Fill(fCent, CDFKolmogorov);
1277 fHistPChi2Root->Fill(CDFROOT, fFitModulation->GetChisquare()/((float)NDF));
1278 fHistPChi2->Fill(CDF, ChiSquare(_tempSwap, fFitModulation)/((float)NDF));
1279 fHistPKolmogorov->Fill(CDF, CDFKolmogorov);
1281 // variable CDF is used for making cuts, so we fill it with the selected p-value
1282 switch (fFitGoodnessTest) {
1286 case kChi2Poisson : break; // CDF is already CDF
1287 case kKolmogorov : {
1288 CDF = CDFKolmogorov;
1294 // as an additional quality check, see if fitting a control fit has a higher significance
1295 _tempSwap.Fit(fFitControl, fFitModulationOptions.Data(), "", lowBound, upBound);
1296 Double_t CDFControl(-1.);
1297 switch (fFitGoodnessTest) {
1299 CDFControl = 1.-ChiSquareCDF(fFitControl->GetNDF(), fFitModulation->GetChisquare());
1301 case kChi2Poisson : {
1302 CDFControl = 1.-ChiSquareCDF(fFitControl->GetNDF(), ChiSquare(_tempSwap, fFitModulation));
1304 case kKolmogorov : {
1305 CDFControl = KolmogorovTest(_tempSwap, fFitControl);
1309 if(CDFControl > CDF) {
1310 CDF = -1.; // control fit is more significant, so throw out the 'old' fit
1311 fHistRhoStatusCent->Fill(fCent, -1);
1314 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
1315 // part of phase space, the requirement that energy desntiy is larger than zero is applied
1316 // to the FULL spectrum
1317 fHistRhoStatusCent->Fill(fCent, 0.);
1318 // for LOCAL didactic purposes, save the best and the worst fits
1319 // this routine can produce a lot of output histograms (it's not memory 'safe') and will not work on GRID
1320 // since the output will become unmergeable (i.e. different nodes may produce conflicting output)
1321 switch (fRunModeType) {
1323 if(fRandom->Uniform(0, 100) > fPercentageOfFits) break;
1324 static Int_t didacticCounterBest(0);
1325 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
1326 TF1* didacticFit = (TF1*)fFitModulation->Clone(Form("fit_%i_CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
1327 switch(fFitModulationType) {
1329 // to make a nice picture also plot the separate components (v2 and v3) of the fit
1330 // only done for cobined fit where there are actually components to split ...
1331 TF1* v2(new TF1("dfit_kV2", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
1332 v2->SetParameter(0, didacticFit->GetParameter(0)); // normalization
1333 v2->SetParameter(3, didacticFit->GetParameter(3)); // v2
1334 v2->FixParameter(1, 1.); // constant
1335 v2->FixParameter(2, 2.); // constant
1336 v2->FixParameter(4, didacticFit->GetParameter(4)); // psi2
1337 v2->SetLineColor(kGreen);
1338 didacticProfile->GetListOfFunctions()->Add(v2);
1339 TF1* v3(new TF1("dfit_kV3", "[0]*([1]+[2]*[3]*TMath::Cos([5]*(x-[4])))", 0, TMath::TwoPi()));
1340 v3->SetParameter(0, didacticFit->GetParameter(0)); // normalization
1341 v3->SetParameter(3, didacticFit->GetParameter(7)); // v3
1342 v3->FixParameter(1, 1.); // constant
1343 v3->FixParameter(2, 2.); // constant
1344 v3->FixParameter(4, didacticFit->GetParameter(6)); // psi3
1345 v3->FixParameter(5, 3.); // constant
1346 v3->SetLineColor(kYellow);
1347 didacticProfile->GetListOfFunctions()->Add(v3);
1351 didacticProfile->GetListOfFunctions()->Add(didacticFit);
1352 fOutputListGood->Add(didacticProfile);
1353 didacticCounterBest++;
1354 TH2F* didacticSurface = BookTH2F(Form("surface_%s", didacticProfile->GetName()), "#phi", "#eta", 50, 0, TMath::TwoPi(), 50, -1, 1, -1, kFALSE);
1355 for(Int_t i(0); i < iTracks; i++) {
1356 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1357 if(PassesCuts(track)) {
1358 if(fUsePtWeight) didacticSurface->Fill(track->Phi(), track->Eta(), track->Pt());
1359 else didacticSurface->Fill(track->Phi(), track->Eta());
1362 if(fExcludeLeadingJetsFromFit) { // visualize the excluded region
1363 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);
1364 f2->SetParameters(excludeInPt/3.,excludeInPhi,.1,excludeInEta,.1);
1365 didacticSurface->GetListOfFunctions()->Add(f2);
1367 fOutputListGood->Add(didacticSurface);
1371 } else { // if the fit is of poor quality revert to the original rho estimate
1372 switch (fRunModeType) { // again see if we want to save the fit
1374 static Int_t didacticCounterWorst(0);
1375 if(fRandom->Uniform(0, 100) > fPercentageOfFits) break;
1376 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data() ));
1377 TF1* didacticFit = (TF1*)fFitModulation->Clone(Form("fit_%i_p_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data()));
1378 didacticProfile->GetListOfFunctions()->Add(didacticFit);
1379 fOutputListBad->Add(didacticProfile);
1380 didacticCounterWorst++;
1384 switch (fFitModulationType) {
1385 case kNoFit : break; // nothing to do
1386 case kCombined : fFitModulation->SetParameter(7, 0); // no break
1387 case kFourierSeries : fFitModulation->SetParameter(7, 0); // no break
1388 default : { // needs to be done if there was a poor fit
1389 fFitModulation->SetParameter(3, 0);
1390 fFitModulation->SetParameter(0, fLocalRho->GetVal());
1393 if(CDF > -.5) fHistRhoStatusCent->Fill(fCent, 1.);
1394 return kFALSE; // return false if the fit is rejected
1398 //_____________________________________________________________________________
1399 Bool_t AliAnalysisTaskRhoVnModulation::PassesCuts(AliVEvent* event)
1402 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1403 // determine the run number to see if the track and jet cuts should be refreshed for semi-good TPC runs
1404 // only done if the runnumber changes, could be moved to a call to AliAnalysisTaskSE::Notify()
1405 if(fRunNumber != InputEvent()->GetRunNumber()) {
1406 fRunNumber = InputEvent()->GetRunNumber(); // set the current run number
1407 if(fDebug > 0) printf("__FUNC__ %s > NEW RUNNUMBER DETECTED \n ", __func__);
1408 // reset the cuts. should be a pointless operation except for the case where the run number changes
1409 // from semi-good back to good on one node, which is not a likely scenario
1410 AliAnalysisTaskEmcal::SetTrackPhiLimits(-10., 10.);
1411 AliAnalysisTaskEmcalJet::SetJetPhiLimits(-10., 10.);
1412 if(fCachedRho) { // if there's a cached rho, it's the default, so switch back
1413 if(fDebug > 0) printf("__FUNC__ %s > replacing rho with cached rho \n ", __func__);
1414 fRho = fCachedRho; // reset rho back to cached value. again, should be pointless
1416 Bool_t flaggedAsSemiGood(kFALSE); // not flagged as anything
1417 for(Int_t i(0); i < fExpectedSemiGoodRuns->GetSize(); i++) {
1418 if(fExpectedSemiGoodRuns->At(i) == fRunNumber) { // run is semi-good
1419 if(fDebug > 0) printf("__FUNC__ %s > semi-good tpc run detected, adjusting acceptance \n ", __func__);
1420 flaggedAsSemiGood = kTRUE;
1421 AliAnalysisTaskEmcalJet::SetJetPhiLimits(fSemiGoodJetMinPhi, fSemiGoodJetMaxPhi); // just an acceptance cut, jets are obtained from full azimuth, so no edge effects
1422 AliAnalysisTaskEmcal::SetTrackPhiLimits(fSemiGoodTrackMinPhi, fSemiGoodTrackMaxPhi); // only affects vn extraction, NOT jet finding
1423 // for semi-good runs, also try to get the 'small rho' estimate, if it is available
1424 AliRhoParameter* tempRho(dynamic_cast<AliRhoParameter*>(InputEvent()->FindListObject(fNameSmallRho.Data())));
1426 if(fDebug > 0) printf("__FUNC__ %s > switching to small rho, caching normal rho \n ", __func__);
1427 fHistAnalysisSummary->SetBinContent(54, 1.); // bookkeep the fact that small rho is used
1428 fCachedRho = fRho; // cache the original rho ...
1429 fRho = tempRho; // ... and use the small rho
1433 if(!flaggedAsSemiGood) {
1434 // in case the run is not a semi-good run, check if it is recognized as another run
1435 // only done to catch unexpected runs
1436 for(Int_t i(0); i < fExpectedRuns->GetSize(); i++) {
1437 if(fExpectedRuns->At(i) == fRunNumber) break; // run is known, break the loop else store the number in a random bin
1438 fHistUndeterminedRunQA->SetBinContent(TMath::Nint(10.*gRandom->Uniform(0.,.9))+1, fRunNumber);
1440 fHistAnalysisSummary->SetBinContent(53, 1.); // bookkeep which rho estimate is used
1443 // continue with event selection
1444 if(!event || !AliAnalysisTaskEmcal::IsEventSelected()) return kFALSE;
1445 if(fSemiCentralInclusive && ! (event->GetTriggerMask() & (ULong64_t(1)<<7))) return kFALSE;
1446 if(TMath::Abs(InputEvent()->GetPrimaryVertex()->GetZ()) > 10.) return kFALSE;
1447 // aod and esd specific checks
1448 switch (fDataType) {
1450 AliESDEvent* esdEvent = static_cast<AliESDEvent*>(InputEvent());
1451 if( (!esdEvent) || (TMath::Abs(esdEvent->GetPrimaryVertexSPD()->GetZ() - esdEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE;
1454 AliAODEvent* aodEvent = static_cast<AliAODEvent*>(InputEvent());
1455 if( (!aodEvent) || (TMath::Abs(aodEvent->GetPrimaryVertexSPD()->GetZ() - aodEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE;
1459 fCent = InputEvent()->GetCentrality()->GetCentralityPercentile("V0M");
1460 if(fCent <= fCentralityClasses->At(0) || fCent >= fCentralityClasses->At(fCentralityClasses->GetSize()-1) || TMath::Abs(fCent-InputEvent()->GetCentrality()->GetCentralityPercentile("TRK")) > 5.) return kFALSE;
1461 // determine centrality class
1462 fInCentralitySelection = -1;
1463 for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
1464 if(fCent >= fCentralityClasses->At(i) && fCent <= fCentralityClasses->At(1+i)) {
1465 fInCentralitySelection = i;
1468 if(fInCentralitySelection<0) return kFALSE; // should be null op
1469 if(fExplicitOutlierCut == 2010 || fExplicitOutlierCut == 2011) {
1470 if(!PassesCuts(fExplicitOutlierCut)) return kFALSE;
1472 if(fTracks->GetEntries() < 1) return kFALSE;
1473 if(fRho->GetVal() <= 0 ) return kFALSE;
1476 //_____________________________________________________________________________
1477 Bool_t AliAnalysisTaskRhoVnModulation::PassesCuts(Int_t year)
1479 // additional centrality cut based on relation between tpc and global multiplicity
1480 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1481 AliAODEvent* event(dynamic_cast<AliAODEvent*>(InputEvent()));
1482 if(!event) return kFALSE;
1483 Int_t multTPC(0), multGlob(0), nTracks(InputEvent()->GetNumberOfTracks());
1484 for(Int_t iTracks = 0; iTracks < nTracks; iTracks++) {
1485 AliAODTrack* track = event->GetTrack(iTracks);
1486 if(!track) continue;
1487 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
1488 if (track->TestFilterBit(1) && track->Chi2perNDF() > 0.2) multTPC++;
1489 if (!track->TestFilterBit(16) || track->Chi2perNDF() < 0.1) continue;
1490 Double_t b[2] = {-99., -99.};
1491 Double_t bCov[3] = {-99., -99., -99.};
1492 if (track->PropagateToDCA(event->GetPrimaryVertex(), event->GetMagneticField(), 100., b, bCov) && TMath::Abs(b[0]) < 0.3 && TMath::Abs(b[1]) < 0.3) multGlob++;
1494 if(year == 2010 && multTPC > (-40.3+1.22*multGlob) && multTPC < (32.1+1.59*multGlob)) return kTRUE;
1495 if(year == 2011 && multTPC > (-36.73 + 1.48*multGlob) && multTPC < (62.87 + 1.78*multGlob)) return kTRUE;
1498 //_____________________________________________________________________________
1499 Bool_t AliAnalysisTaskRhoVnModulation::PassesCuts(const AliVCluster* cluster) const
1502 if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1503 if(!cluster) return kFALSE;
1506 //_____________________________________________________________________________
1507 void AliAnalysisTaskRhoVnModulation::FillHistogramsAfterSubtraction(Double_t psi2, Double_t psi3, Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc)
1510 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1511 FillTrackHistograms();
1512 /* FillClusterHistograms(); */
1513 FillJetHistograms(psi2, psi3);
1514 /* FillCorrectedClusterHistograms(); */
1515 if(fFillQAHistograms) FillEventPlaneHistograms(vzero, vzeroComb, tpc);
1516 FillRhoHistograms();
1517 FillDeltaPtHistograms(psi2, psi3);
1519 //_____________________________________________________________________________
1520 void AliAnalysisTaskRhoVnModulation::FillTrackHistograms() const
1522 // fill track histograms
1523 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1524 Int_t iTracks(fTracks->GetEntriesFast()), iAcceptedTracks(0);
1525 for(Int_t i(0); i < iTracks; i++) {
1526 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1527 if(!PassesCuts(track)) continue;
1529 fHistPicoTrackPt[fInCentralitySelection]->Fill(track->Pt());
1530 if(fFillQAHistograms) FillQAHistograms(track);
1532 fHistPicoTrackMult[fInCentralitySelection]->Fill(iAcceptedTracks);
1534 //_____________________________________________________________________________
1535 void AliAnalysisTaskRhoVnModulation::FillClusterHistograms() const
1537 // fill cluster histograms
1538 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1539 Int_t iClusters(fCaloClusters->GetEntriesFast());
1540 for(Int_t i(0); i < iClusters; i++) {
1541 AliVCluster* cluster = static_cast<AliVCluster*>(fCaloClusters->At(iClusters));
1542 if (!PassesCuts(cluster)) continue;
1543 TLorentzVector clusterLorentzVector;
1544 cluster->GetMomentum(clusterLorentzVector, const_cast<Double_t*>(fVertex));
1545 //fHistClusterPt[fInCentralitySelection]->Fill(clusterLorentzVector.Pt());
1546 //fHistClusterEta[fInCentralitySelection]->Fill(clusterLorentzVector.Eta());
1547 //fHistClusterPhi[fInCentralitySelection]->Fill(clusterLorentzVector.Phi());
1551 //_____________________________________________________________________________
1552 void AliAnalysisTaskRhoVnModulation::FillCorrectedClusterHistograms() const
1554 // fill clusters after hadronic correction FIXME implement
1555 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1557 //_____________________________________________________________________________
1558 void AliAnalysisTaskRhoVnModulation::FillEventPlaneHistograms(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc) const
1560 // fill event plane histograms
1561 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1562 fHistPsiControl->Fill(0.5, vzero[0][0]); // vzero a psi2
1563 fHistPsiControl->Fill(1.5, vzero[1][0]); // vzero c psi2
1564 fHistPsiControl->Fill(2.5, tpc[0]); // tpc psi 2
1565 fHistPsiControl->Fill(5.5, vzero[0][1]); // vzero a psi3
1566 fHistPsiControl->Fill(6.5, vzero[1][1]); // vzero b psi3
1567 fHistPsiControl->Fill(7.5, tpc[1]); // tpc psi 3
1568 fHistPsiVZEROA->Fill(vzero[0][0]);
1569 fHistPsiVZEROC->Fill(vzero[1][0]);
1570 fHistPsiVZERO->Fill(vzeroComb[0]);
1571 fHistPsiTPC->Fill(tpc[0]);
1572 fHistPsiSpread->Fill(0.5, TMath::Abs(vzero[0][0]-vzero[1][0]));
1573 fHistPsiSpread->Fill(1.5, TMath::Abs(vzero[0][0]-tpc[0]));
1574 fHistPsiSpread->Fill(2.5, TMath::Abs(vzero[1][0]-tpc[0]));
1575 // event plane vs centrality QA histo's to check recentering
1576 Double_t TRK(InputEvent()->GetCentrality()->GetCentralityPercentile("TRK"));
1577 Double_t V0M(InputEvent()->GetCentrality()->GetCentralityPercentile("V0M"));
1578 fHistPsiVZEROAV0M->Fill(V0M, vzero[0][0]);
1579 fHistPsiVZEROCV0M->Fill(V0M, vzero[1][0]);
1580 fHistPsiVZEROVV0M->Fill(V0M, vzeroComb[0]);
1581 fHistPsiTPCiV0M->Fill(V0M, tpc[0]);
1582 fHistPsiVZEROATRK->Fill(TRK, vzero[0][0]);
1583 fHistPsiVZEROCTRK->Fill(TRK, vzero[1][0]);
1584 fHistPsiVZEROTRK->Fill(TRK, vzeroComb[0]);
1585 fHistPsiTPCTRK->Fill(TRK, tpc[0]);
1587 //_____________________________________________________________________________
1588 void AliAnalysisTaskRhoVnModulation::FillRhoHistograms()
1590 // fill rho histograms
1591 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1592 fHistRhoPackage[fInCentralitySelection]->Fill(fLocalRho->GetVal()); // save the rho estimate from the emcal jet package
1593 // get multiplicity FIXME inefficient
1594 Int_t iJets(fJets->GetEntriesFast());
1595 Double_t rho(fLocalRho->GetLocalVal(TMath::Pi(), TMath::Pi(), fLocalRho->GetVal()));
1596 fHistRho[fInCentralitySelection]->Fill(rho);
1597 fHistRhoVsMult->Fill(fTracks->GetEntries(), rho);
1598 fHistRhoVsCent->Fill(fCent, rho);
1599 for(Int_t i(0); i < iJets; i++) {
1600 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
1601 if(!PassesCuts(jet)) continue;
1602 fHistRhoAVsMult->Fill(fTracks->GetEntries(), rho * jet->Area());
1603 fHistRhoAVsCent->Fill(fCent, rho * jet->Area());
1606 //_____________________________________________________________________________
1607 void AliAnalysisTaskRhoVnModulation::FillDeltaPtHistograms(Double_t psi2, Double_t psi3) const
1609 // fill delta pt histograms
1610 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1612 const Float_t areaRC = fRandomConeRadius*fRandomConeRadius*TMath::Pi();
1613 // we're retrieved the leading jet, now get a random cone
1614 for(i = 0; i < fMaxCones; i++) {
1615 Float_t pt(0), eta(0), phi(0);
1616 // get a random cone without constraints on leading jet position
1617 CalculateRandomCone(pt, eta, phi, 0x0);
1619 if(fFillQAHistograms) fHistRCPhiEta[fInCentralitySelection]->Fill(phi, eta);
1620 fHistRhoVsRCPt[fInCentralitySelection]->Fill(pt, fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal())*areaRC);
1621 fHistRCPt[fInCentralitySelection]->Fill(pt);
1622 fHistDeltaPtDeltaPhi2[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1623 fHistDeltaPtDeltaPhi3[fInCentralitySelection]->Fill(PhaseShift(phi-psi3, 3.), pt - areaRC*fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1625 // get a random cone excluding leading jet area
1626 CalculateRandomCone(pt, eta, phi, fLeadingJet);
1628 if(fFillQAHistograms) fHistRCPhiEtaExLJ[fInCentralitySelection]->Fill(phi, eta);
1629 fHistRhoVsRCPtExLJ[fInCentralitySelection]->Fill(pt, fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal())*areaRC);
1630 fHistRCPtExLJ[fInCentralitySelection]->Fill(pt);
1631 fHistDeltaPtDeltaPhi2ExLJ[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1632 fHistDeltaPtDeltaPhi3ExLJ[fInCentralitySelection]->Fill(PhaseShift(phi-psi3, 3.), pt - areaRC*fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1636 //_____________________________________________________________________________
1637 void AliAnalysisTaskRhoVnModulation::FillJetHistograms(Double_t psi2, Double_t psi3)
1639 // fill jet histograms
1640 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1641 Int_t iJets(fJets->GetEntriesFast());
1642 for(Int_t i(0); i < iJets; i++) {
1643 AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
1644 if(PassesCuts(jet)) {
1645 Double_t pt(jet->Pt()), area(jet->Area()), eta(jet->Eta()), phi(jet->Phi());
1646 Double_t rho(fLocalRho->GetLocalVal(phi, GetJetContainer()->GetJetRadius(), fLocalRho->GetVal()));
1647 fHistJetPtRaw[fInCentralitySelection]->Fill(pt);
1648 fHistJetPt[fInCentralitySelection]->Fill(pt-area*rho);
1649 if(fFillQAHistograms) fHistJetEtaPhi[fInCentralitySelection]->Fill(eta, phi);
1650 fHistJetPtArea[fInCentralitySelection]->Fill(pt-area*rho, area);
1651 fHistJetPsi2Pt[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt-area*rho);
1652 fHistJetPsi3Pt[fInCentralitySelection]->Fill(PhaseShift(phi-psi3, 3.), pt-area*rho);
1653 fHistJetPtConstituents[fInCentralitySelection]->Fill(pt-area*rho, jet->Nch());
1654 fHistJetEtaRho[fInCentralitySelection]->Fill(eta, pt/area);
1655 if(fSubtractJetPt) jet->SetPtSub(pt-area*rho); // if requested, save the subtracted jet pt
1656 } else if(fSubtractJetPt) jet->SetPtSub(-999.);
1659 //_____________________________________________________________________________
1660 void AliAnalysisTaskRhoVnModulation::FillQAHistograms(AliVTrack* vtrack) const
1662 // fill qa histograms for pico tracks
1664 AliPicoTrack* track = static_cast<AliPicoTrack*>(vtrack);
1665 fHistRunnumbersPhi->Fill(fMappedRunNumber, track->Phi());
1666 fHistRunnumbersEta->Fill(fMappedRunNumber, track->Eta());
1667 Int_t type((int)(track->GetTrackType()));
1670 fHistPicoCat1[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1673 fHistPicoCat2[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1676 fHistPicoCat3[fInCentralitySelection]->Fill(track->Eta(), track->Phi());
1681 //_____________________________________________________________________________
1682 void AliAnalysisTaskRhoVnModulation::FillQAHistograms(AliVEvent* vevent)
1684 // fill qa histograms for events
1686 fHistVertexz->Fill(vevent->GetPrimaryVertex()->GetZ());
1687 fHistCentrality->Fill(fCent);
1688 Int_t runNumber(InputEvent()->GetRunNumber());
1689 for(fMappedRunNumber = 0; fExpectedRuns->GetSize()+1; fMappedRunNumber++) {
1690 if(fExpectedRuns->At(fMappedRunNumber) == runNumber) break;
1693 //_____________________________________________________________________________
1694 void AliAnalysisTaskRhoVnModulation::FillAnalysisSummaryHistogram() const
1696 // fill the analysis summary histrogram, saves all relevant analysis settigns
1697 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1698 fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fJetRadius");
1699 fHistAnalysisSummary->SetBinContent(2, GetJetContainer()->GetJetRadius());
1700 fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fJetEtaMin");
1701 fHistAnalysisSummary->SetBinContent(3, GetJetContainer()->GetJetEtaMin());
1702 fHistAnalysisSummary->GetXaxis()->SetBinLabel(4, "fJetEtaMax");
1703 fHistAnalysisSummary->SetBinContent(4, GetJetContainer()->GetJetEtaMax());
1704 fHistAnalysisSummary->GetXaxis()->SetBinLabel(5, "fJetPhiMin");
1705 fHistAnalysisSummary->SetBinContent(5, GetJetContainer()->GetJetPhiMin());
1706 fHistAnalysisSummary->GetXaxis()->SetBinLabel(6, "fJetPhiMax");
1707 fHistAnalysisSummary->SetBinContent(6, GetJetContainer()->GetJetPhiMin());
1708 fHistAnalysisSummary->GetXaxis()->SetBinLabel(16, "fForceBeamType");
1709 fHistAnalysisSummary->SetBinContent(16, fForceBeamType);
1710 fHistAnalysisSummary->GetXaxis()->SetBinLabel(17, "fMinCent");
1711 fHistAnalysisSummary->SetBinContent(17, fMinCent);
1712 fHistAnalysisSummary->GetXaxis()->SetBinLabel(18, "fMaxCent");
1713 fHistAnalysisSummary->SetBinContent(18, fMaxCent);
1714 fHistAnalysisSummary->GetXaxis()->SetBinLabel(19, "fMinVz");
1715 fHistAnalysisSummary->SetBinContent(19, fMinVz);
1716 fHistAnalysisSummary->GetXaxis()->SetBinLabel(20, "fMaxVz");
1717 fHistAnalysisSummary->SetBinContent(20, fMaxVz);
1718 fHistAnalysisSummary->GetXaxis()->SetBinLabel(21, "fOffTrigger");
1719 fHistAnalysisSummary->SetBinContent(21, fOffTrigger);
1720 fHistAnalysisSummary->GetXaxis()->SetBinLabel(33, "fRandomConeRadius");
1721 fHistAnalysisSummary->SetBinContent(33, fRandomConeRadius);
1722 fHistAnalysisSummary->GetXaxis()->SetBinLabel(34, "fitModulationType");
1723 fHistAnalysisSummary->SetBinContent(34, (int)fFitModulationType);
1724 fHistAnalysisSummary->GetXaxis()->SetBinLabel(35, "runModeType");
1725 fHistAnalysisSummary->SetBinContent(35, (int)fRunModeType);
1726 fHistAnalysisSummary->GetXaxis()->SetBinLabel(36, "data type");
1727 fHistAnalysisSummary->SetBinContent(36, (int)fDataType);
1728 fHistAnalysisSummary->GetXaxis()->SetBinLabel(37, "iterator");
1729 fHistAnalysisSummary->SetBinContent(37, 1.);
1730 fHistAnalysisSummary->GetXaxis()->SetBinLabel(38, "fMinPvalue");
1731 fHistAnalysisSummary->SetBinContent(38, fMinPvalue);
1732 fHistAnalysisSummary->GetXaxis()->SetBinLabel(39, "fMaxPvalue");
1733 fHistAnalysisSummary->SetBinContent(39, fMaxPvalue);
1734 fHistAnalysisSummary->GetXaxis()->SetBinLabel(40, "fExcludeLeadingJetsFromFit");
1735 fHistAnalysisSummary->SetBinContent(40, fExcludeLeadingJetsFromFit);
1736 fHistAnalysisSummary->GetXaxis()->SetBinLabel(41, "fRebinSwapHistoOnTheFly");
1737 fHistAnalysisSummary->SetBinContent(41, (int)fRebinSwapHistoOnTheFly);
1738 fHistAnalysisSummary->GetXaxis()->SetBinLabel(42, "fUsePtWeight");
1739 fHistAnalysisSummary->SetBinContent(42, (int)fUsePtWeight);
1740 fHistAnalysisSummary->GetXaxis()->SetBinLabel(43, "fMinLeadingHadronPt");
1741 fHistAnalysisSummary->SetBinContent(43, fMinLeadingHadronPt);
1742 fHistAnalysisSummary->GetXaxis()->SetBinLabel(44, "fExplicitOutlierCut");
1743 fHistAnalysisSummary->SetBinContent(44, fExplicitOutlierCut);
1744 fHistAnalysisSummary->GetXaxis()->SetBinLabel(45, "fLocalJetMinEta");
1745 fHistAnalysisSummary->SetBinContent(45,fLocalJetMinEta );
1746 fHistAnalysisSummary->GetXaxis()->SetBinLabel(46, "fLocalJetMaxEta");
1747 fHistAnalysisSummary->SetBinContent(46, fLocalJetMaxEta);
1748 fHistAnalysisSummary->GetXaxis()->SetBinLabel(47, "fLocalJetMinPhi");
1749 fHistAnalysisSummary->SetBinContent(47, fLocalJetMinPhi);
1750 fHistAnalysisSummary->GetXaxis()->SetBinLabel(48, "fLocalJetMaxPhi");
1751 fHistAnalysisSummary->SetBinContent(48, fLocalJetMaxPhi);
1752 fHistAnalysisSummary->GetXaxis()->SetBinLabel(49, "fSoftTrackMinPt");
1753 fHistAnalysisSummary->SetBinContent(49, fSoftTrackMinPt);
1754 fHistAnalysisSummary->GetXaxis()->SetBinLabel(50, "fSoftTrackMaxPt");
1755 fHistAnalysisSummary->SetBinContent(50, fSoftTrackMaxPt);
1756 fHistAnalysisSummary->GetXaxis()->SetBinLabel(51, "fMaxCones");
1757 fHistAnalysisSummary->SetBinContent(51, fMaxCones);
1758 fHistAnalysisSummary->GetXaxis()->SetBinLabel(52, "fUseScaledRho");
1759 fHistAnalysisSummary->SetBinContent(52, fUseScaledRho);
1760 fHistAnalysisSummary->GetXaxis()->SetBinLabel(53, "used rho");
1761 fHistAnalysisSummary->GetXaxis()->SetBinLabel(54, "used small rho");
1763 //_____________________________________________________________________________
1764 void AliAnalysisTaskRhoVnModulation::Terminate(Option_t *)
1767 switch (fRunModeType) {
1769 printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1770 AliAnalysisTaskRhoVnModulation::Dump();
1771 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));
1776 //_____________________________________________________________________________
1777 void AliAnalysisTaskRhoVnModulation::SetModulationFit(TF1* fit)
1779 // set modulation fit
1780 if (fFitModulation) delete fFitModulation;
1781 fFitModulation = fit;
1783 //_____________________________________________________________________________
1784 void AliAnalysisTaskRhoVnModulation::SetUseControlFit(Bool_t c)
1787 if (fFitControl) delete fFitControl;
1789 fFitControl = new TF1("controlFit", "pol0", 0, TMath::TwoPi());
1790 } else fFitControl = 0x0;
1792 //_____________________________________________________________________________
1793 TH1F* AliAnalysisTaskRhoVnModulation::GetResolutionFromOuptutFile(detectorType det, Int_t h, TArrayD* cen)
1795 // INTERFACE METHOD FOR OUTPUTFILE
1796 // get the detector resolution, user has ownership of the returned histogram
1798 printf(" > Please add fOutputList first < \n");
1802 (cen) ? r = new TH1F("R", "R", cen->GetSize()-1, cen->GetArray()) : r = new TH1F("R", "R", 10, 0, 10);
1803 if(!cen) r->GetXaxis()->SetTitle("number of centrality bin");
1804 r->GetYaxis()->SetTitle(Form("Resolution #Psi_{%i}", h));
1805 for(Int_t i(0); i < 10; i++) {
1806 TProfile* temp((TProfile*)fOutputList->FindObject(Form("fProfV%iResolution_%i", h, i)));
1808 Double_t a(temp->GetBinContent(3)), b(temp->GetBinContent(5)), c(temp->GetBinContent(7));
1809 Double_t d(temp->GetBinContent(9)), e(temp->GetBinContent(10)), f(temp->GetBinContent(11));
1810 Double_t _a(temp->GetBinError(3)), _b(temp->GetBinError(5)), _c(temp->GetBinError(7));
1811 Double_t _d(temp->GetBinError(9)), _e(temp->GetBinError(10)), _f(temp->GetBinError(11));
1812 if(a <= 0 || b <= 0 || c <= 0 || d <= 0 || e <= 0 || f <= 0) continue;
1815 r->SetBinContent(1+i, TMath::Sqrt((a*b)/c));
1816 if(i==0) r->SetNameTitle("VZEROA resolution", "VZEROA resolution");
1817 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
1820 r->SetBinContent(1+i, TMath::Sqrt((a*c)/b));
1821 if(i==0) r->SetNameTitle("VZEROC resolution", "VZEROC resolution");
1822 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
1825 r->SetBinContent(1+i, TMath::Sqrt((b*c)/a));
1826 if(i==0) r->SetNameTitle("TPC resolution", "TPC resolution");
1827 r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
1830 r->SetBinContent(1+i, TMath::Sqrt((d*e)/f));
1831 if(i==0) r->SetNameTitle("VZEROComb resolution", "VZEROComb resolution");
1832 r->SetBinError(1+i, TMath::Sqrt(_d*_d+_e*_e+_f*_f));
1839 //_____________________________________________________________________________
1840 TH1F* AliAnalysisTaskRhoVnModulation::CorrectForResolutionDiff(TH1F* v, detectorType det, TArrayD* cen, Int_t c, Int_t h)
1842 // INTERFACE METHOD FOR OUTPUT FILE
1843 // correct the supplied differential vn histogram v for detector resolution
1844 TH1F* r(GetResolutionFromOuptutFile(det, h, cen));
1846 printf(" > Couldn't find resolution < \n");
1849 Double_t res(1./r->GetBinContent(1+r->FindBin(c)));
1850 TF1* line = new TF1("line", "pol0", 0, 200);
1851 line->SetParameter(0, res);
1855 //_____________________________________________________________________________
1856 TH1F* AliAnalysisTaskRhoVnModulation::CorrectForResolutionInt(TH1F* v, detectorType det, TArrayD* cen, Int_t h)
1858 // INTERFACE METHOD FOR OUTPUT FILE
1859 // correct the supplied intetrated vn histogram v for detector resolution
1860 // integrated vn must have the same centrality binning as the resolotion correction
1861 TH1F* r(GetResolutionFromOuptutFile(det, h, cen));
1865 //_____________________________________________________________________________
1866 TH1F* AliAnalysisTaskRhoVnModulation::GetDifferentialQC(TProfile* refCumulants, TProfile* diffCumlants, TArrayD* ptBins, Int_t h)
1868 // get differential QC
1869 Double_t r(refCumulants->GetBinContent(h-1)); // v2 reference flow
1870 if(r > 0) r = TMath::Sqrt(r);
1871 TH1F* qc = new TH1F(Form("QC2v%i", h), Form("QC2v%i", h), ptBins->GetSize()-1, ptBins->GetArray());
1872 Double_t a(0), b(0), c(0); // dummy variables
1873 for(Int_t i(0); i < ptBins->GetSize(); i++) {
1875 a = diffCumlants->GetBinContent(1+i);
1876 b = diffCumlants->GetBinError(1+i);
1878 qc->SetBinContent(1+i, c);
1879 (a <= 0 || b <= 0) ? qc->SetBinError(1+i, b) : qc->SetBinError(1+i, TMath::Sqrt(c*c*b*b/(a*a)));
1885 //_____________________________________________________________________________