]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskRhoVnModulation.cxx
1) additional interface methods for AliAnalysisTaskjetFlow 2) extra histograms to...
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / UserTasks / AliAnalysisTaskRhoVnModulation.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* 
17  * analysis task for jet flow preparation
18  *
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:
21  *      - (anti-kt) jets
22  *      - background estimate rho
23  *      - pico tracks
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
27  *
28  * author: Redmer Alexander Bertens, Utrecht Univeristy, Utrecht, Netherlands
29  * rbertens@cern.ch, rbertens@nikhef.nl, r.a.bertens@uu.nl 
30  */
31
32 // root includes
33 #include <TStyle.h>
34 #include <TRandom3.h>
35 #include <TChain.h>
36 #include <TMath.h>
37 #include <TF1.h>
38 #include <TF2.h>
39 #include <TH1F.h>
40 #include <TH2F.h>
41 #include <TProfile.h>
42 // aliroot includes
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 // local includes
55 #include "AliAnalysisTaskRhoVnModulation.h"
56
57
58 class AliAnalysisTaskRhoVnModulation;
59 using namespace std;
60
61 ClassImp(AliAnalysisTaskRhoVnModulation)
62
63 AliAnalysisTaskRhoVnModulation::AliAnalysisTaskRhoVnModulation() : AliAnalysisTaskEmcalJet("AliAnalysisTaskRhoVnModulation", kTRUE), 
64     fDebug(0), fInitialized(0), fFillQAHistograms(kTRUE), fReduceBinsXByFactor(1), fReduceBinsYByFactor(1), fCentralityClasses(0), fPtBinsHybrids(0), fPtBinsJets(0), fUserSuppliedV2(0), fUserSuppliedV3(0), fUserSuppliedR2(0), fUserSuppliedR3(0), fNAcceptedTracks(0), fNAcceptedTracksQCn(0), fFitModulationType(kNoFit), fUsePtWeight(kTRUE), fDetectorType(kTPC), fFitModulationOptions("Q"), fRunModeType(kGrid), fDataType(kESD), fRandom(0), fMappedRunNumber(0), fInCentralitySelection(-1), fFitModulation(0), fMinPvalue(0), fMaxPvalue(1), fNameJetClones(0), fNamePicoTrackClones(0), fNameRho(0), fLocalJetMinEta(-10), fLocalJetMaxEta(-10), fLocalJetMinPhi(-10), fLocalJetMaxPhi(-10), fSoftTrackMinPt(0.15), fSoftTrackMaxPt(5.), fAbsVertexZ(10), fHistCentrality(0), fHistVertexz(0), fHistRunnumbersPhi(0), fHistRunnumbersEta(0), fHistPvaluePDF(0), fHistPvalueCDF(0), fMinDisanceRCtoLJ(0), fRandomConeRadius(-1.), fAbsVnHarmonics(kTRUE), fExcludeLeadingJetsFromFit(1.), fRebinSwapHistoOnTheFly(kTRUE), fPercentageOfFits(10.), fUseV0EventPlaneFromHeader(kFALSE), fSetPtSub(kFALSE), fExplicitOutlierCut(-1), fMinLeadingHadronPt(0), fOutputList(0), fOutputListGood(0), fOutputListBad(0), fHistAnalysisSummary(0), fHistSwap(0), fProfV2(0), fProfV2Cumulant(0), fProfV3(0), fProfV3Cumulant(0), fHistPsiControl(0), fHistPsiSpread(0), fHistPsiVZEROA(0), fHistPsiVZEROC(0), fHistPsiTPC(0), fHistRhoVsMult(0), fHistRhoVsCent(0), fHistRhoAVsMult(0), fHistRhoAVsCent(0) {
65     for(Int_t i(0); i < 10; i++) {
66         fProfV2Resolution[i] = 0;
67         fProfV3Resolution[i] = 0;
68         fHistPicoTrackPt[i] = 0;
69         fHistPicoCat1[i] = 0;
70         fHistPicoCat2[i] = 0;
71         fHistPicoCat3[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;
79         fHistRho[i] = 0;
80         fHistRCPhiEta[i] = 0;
81         fHistRhoVsRCPt[i] = 0;
82         fHistRCPt[i] = 0;
83         fHistDeltaPtDeltaPhi2TPC[i] = 0;
84         fHistDeltaPtDeltaPhi2V0A[i] = 0;
85         fHistDeltaPtDeltaPhi2V0C[i] = 0;
86         fHistDeltaPtDeltaPhi3TPC[i] = 0;
87         fHistDeltaPtDeltaPhi3V0A[i] = 0;
88         fHistDeltaPtDeltaPhi3V0C[i] = 0;
89         fHistRCPhiEtaExLJ[i] = 0;
90         fHistRhoVsRCPtExLJ[i] = 0;
91         fHistRCPtExLJ[i] = 0;
92         fHistDeltaPtDeltaPhi2ExLJTPC[i] = 0;
93         fHistDeltaPtDeltaPhi2ExLJV0A[i] = 0;
94         fHistDeltaPtDeltaPhi2ExLJV0C[i] = 0;
95         fHistDeltaPtDeltaPhi3ExLJTPC[i] = 0;
96         fHistDeltaPtDeltaPhi3ExLJV0A[i] = 0;
97         fHistDeltaPtDeltaPhi3ExLJV0C[i] = 0;
98         /* fHistRCPhiEtaRand[i] = 0; */
99         /* fHistRhoVsRCPtRand[i] = 0; */
100         /* fHistRCPtRand[i] = 0; */
101         /* fHistDeltaPtDeltaPhi2Rand[i] = 0; */
102         /* fHistDeltaPtDeltaPhi3Rand[i] = 0; */
103         fHistJetPtRaw[i] = 0;
104         fHistJetPt[i] = 0;
105         fHistJetEtaPhi[i] = 0;
106         fHistJetPtArea[i] = 0;
107         fHistJetPtConstituents[i] = 0;
108         fHistJetEtaRho[i] = 0;
109         fHistJetPsiTPCPt[i] = 0;
110         fHistJetPsiVZEROAPt[i] = 0;
111         fHistJetPsiVZEROCPt[i] = 0;
112         fHistDeltaPhi2VZEROA[i] = 0;
113         fHistDeltaPhi2VZEROC[i] = 0;
114         fHistDeltaPhi2TPC[i] = 0;
115         fHistDeltaPhi3VZEROA[i] = 0;
116         fHistDeltaPhi3VZEROC[i] = 0;
117         fHistDeltaPhi3TPC[i] = 0;
118    }
119     // default constructor
120 }
121 //_____________________________________________________________________________
122 AliAnalysisTaskRhoVnModulation::AliAnalysisTaskRhoVnModulation(const char* name, runModeType type) : AliAnalysisTaskEmcalJet(name, kTRUE),
123   fDebug(0), fInitialized(0), fFillQAHistograms(kTRUE), fReduceBinsXByFactor(1), fReduceBinsYByFactor(1), fCentralityClasses(0), fPtBinsHybrids(0), fPtBinsJets(0), fUserSuppliedV2(0), fUserSuppliedV3(0), fUserSuppliedR2(0), fUserSuppliedR3(0), fNAcceptedTracks(0), fNAcceptedTracksQCn(0), fFitModulationType(kNoFit), fUsePtWeight(kTRUE), fDetectorType(kTPC), fFitModulationOptions("Q"), fRunModeType(type), fDataType(kESD), fRandom(0), fMappedRunNumber(0), fInCentralitySelection(-1), fFitModulation(0), fMinPvalue(0), fMaxPvalue(1), fNameJetClones(0), fNamePicoTrackClones(0), fNameRho(0), fLocalJetMinEta(-10), fLocalJetMaxEta(-10), fLocalJetMinPhi(-10), fLocalJetMaxPhi(-10), fSoftTrackMinPt(0.15), fSoftTrackMaxPt(5.),  fAbsVertexZ(10), fHistCentrality(0), fHistVertexz(0), fHistRunnumbersPhi(0), fHistRunnumbersEta(0), fHistPvaluePDF(0), fHistPvalueCDF(0), fMinDisanceRCtoLJ(0), fRandomConeRadius(-1.), fAbsVnHarmonics(kTRUE), fExcludeLeadingJetsFromFit(1.), fRebinSwapHistoOnTheFly(kTRUE), fPercentageOfFits(10.), fUseV0EventPlaneFromHeader(kFALSE), fSetPtSub(kFALSE), fExplicitOutlierCut(-1), fMinLeadingHadronPt(0), fOutputList(0), fOutputListGood(0), fOutputListBad(0), fHistAnalysisSummary(0), fHistSwap(0), fProfV2(0), fProfV2Cumulant(0), fProfV3(0), fProfV3Cumulant(0), fHistPsiControl(0), fHistPsiSpread(0), fHistPsiVZEROA(0), fHistPsiVZEROC(0), fHistPsiTPC(0), fHistRhoVsMult(0), fHistRhoVsCent(0), fHistRhoAVsMult(0), fHistRhoAVsCent(0) {
124     for(Int_t i(0); i < 10; i++) {
125         fProfV2Resolution[i] = 0;
126         fProfV3Resolution[i] = 0;
127         fHistPicoTrackPt[i] = 0;
128         fHistPicoCat1[i] = 0;
129         fHistPicoCat2[i] = 0;
130         fHistPicoCat3[i] = 0;
131         /* fHistClusterPt[i] = 0; */
132         /* fHistClusterPhi[i] = 0; */
133         /* fHistClusterEta[i] = 0; */ 
134         /* fHistClusterCorrPt[i] = 0; */
135         /* fHistClusterCorrPhi[i] = 0; */
136         /* fHistClusterCorrEta[i] = 0; */
137         fHistRhoPackage[i] = 0;
138         fHistRho[i] = 0;
139         fHistRCPhiEta[i] = 0;
140         fHistRhoVsRCPt[i] = 0;
141         fHistRCPt[i] = 0;
142         fHistDeltaPtDeltaPhi2TPC[i] = 0;
143         fHistDeltaPtDeltaPhi2V0A[i] = 0;
144         fHistDeltaPtDeltaPhi2V0C[i] = 0;
145         fHistDeltaPtDeltaPhi3TPC[i] = 0;
146         fHistDeltaPtDeltaPhi3V0A[i] = 0;
147         fHistDeltaPtDeltaPhi3V0C[i] = 0;
148         fHistRCPhiEtaExLJ[i] = 0;
149         fHistRhoVsRCPtExLJ[i] = 0;
150         fHistRCPtExLJ[i] = 0;
151         fHistDeltaPtDeltaPhi2ExLJTPC[i] = 0;
152         fHistDeltaPtDeltaPhi2ExLJV0A[i] = 0;
153         fHistDeltaPtDeltaPhi2ExLJV0C[i] = 0;
154         fHistDeltaPtDeltaPhi3ExLJTPC[i] = 0;
155         fHistDeltaPtDeltaPhi3ExLJV0A[i] = 0;
156         fHistDeltaPtDeltaPhi3ExLJV0C[i] = 0;
157         /* fHistRCPhiEtaRand[i] = 0; */
158         /* fHistRhoVsRCPtRand[i] = 0; */
159         /* fHistRCPtRand[i] = 0; */
160         /* fHistDeltaPtDeltaPhi2Rand[i] = 0; */
161         /* fHistDeltaPtDeltaPhi3Rand[i] = 0; */
162         fHistJetPtRaw[i] = 0;
163         fHistJetPt[i] = 0;
164         fHistJetEtaPhi[i] = 0;
165         fHistJetPtArea[i] = 0;
166         fHistJetPtConstituents[i] = 0;
167         fHistJetEtaRho[i] = 0;
168         fHistJetPsiTPCPt[i] = 0;
169         fHistJetPsiVZEROAPt[i] = 0;
170         fHistJetPsiVZEROCPt[i] = 0;
171         fHistDeltaPhi2VZEROA[i] = 0;
172         fHistDeltaPhi2VZEROC[i] = 0;
173         fHistDeltaPhi2TPC[i] = 0;
174         fHistDeltaPhi3VZEROA[i] = 0;
175         fHistDeltaPhi3VZEROC[i] = 0;
176         fHistDeltaPhi3TPC[i] = 0;
177     }
178     // constructor
179     DefineInput(0, TChain::Class());
180     DefineOutput(1, TList::Class());
181     switch (fRunModeType) {
182         case kLocal : {
183             gStyle->SetOptFit(1);
184             DefineOutput(2, TList::Class());
185             DefineOutput(3, TList::Class());
186         } break;
187         default: fDebug = -1;   // suppress debug info explicitely when not running locally
188     }
189 }
190 //_____________________________________________________________________________
191 AliAnalysisTaskRhoVnModulation::~AliAnalysisTaskRhoVnModulation()
192 {
193     // destructor
194     if(fOutputList)             delete fOutputList;
195     if(fOutputListGood)         delete fOutputListGood;
196     if(fOutputListBad)          delete fOutputListBad;
197     if(fFitModulation)          delete fFitModulation;
198     if(fHistSwap)               delete fHistSwap;
199     if(fCentralityClasses)      delete fCentralityClasses;
200 }
201 //_____________________________________________________________________________
202 Bool_t AliAnalysisTaskRhoVnModulation::InitializeAnalysis() 
203 {
204     // initialize the anaysis
205     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
206     if(fRandomConeRadius <= 0) fRandomConeRadius = fJetRadius;
207     if(fLocalJetMinEta > -10 && fLocalJetMaxEta > -10) SetJetEtaLimits(fLocalJetMinEta, fLocalJetMaxEta);
208     if(fLocalJetMinPhi > -10 && fLocalJetMaxPhi > -10) SetJetPhiLimits(fLocalJetMinPhi, fLocalJetMaxPhi);
209     if(fMinDisanceRCtoLJ==0) fMinDisanceRCtoLJ = .5*fJetRadius;
210     if(dynamic_cast<AliAODEvent*>(InputEvent())) fDataType = kAOD; // determine the datatype
211     else if(dynamic_cast<AliESDEvent*>(InputEvent())) fDataType = kESD;
212     fHistAnalysisSummary->SetBinContent(36, (int)fDataType);
213     if(!fRandom) fRandom = new TRandom3(0);  // get a randomized if one hasn't been user-supplied
214     switch (fFitModulationType)  {
215         case kNoFit : { SetModulationFit(new TF1("fix_kNoFit", "[0]", 0, TMath::TwoPi())); } break;
216         case kV2 : {
217             SetModulationFit(new TF1("fit_kV2", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
218             fFitModulation->SetParameter(0, 0.);        // normalization
219             fFitModulation->SetParameter(3, 0.2);       // v2
220             fFitModulation->FixParameter(1, 1.);        // constant
221             fFitModulation->FixParameter(2, 2.);        // constant
222         } break;
223         case kV3: {
224             SetModulationFit(new TF1("fit_kV3", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
225             fFitModulation->SetParameter(0, 0.);        // normalization
226             fFitModulation->SetParameter(3, 0.2);       // v3
227             fFitModulation->FixParameter(1, 1.);        // constant
228             fFitModulation->FixParameter(2, 3.);        // constant
229         } break;
230         default : { // for the combined fit, the 'direct fourier series' or the user supplied vn values we use v2 and v3
231              SetModulationFit(new TF1("fit_kCombined", "[0]*([1]+[2]*([3]*TMath::Cos([2]*(x-[4]))+[7]*TMath::Cos([5]*(x-[6]))))", 0, TMath::TwoPi()));
232              fFitModulation->SetParameter(0, 0.);       // normalization
233              fFitModulation->SetParameter(3, 0.2);      // v2
234              fFitModulation->FixParameter(1, 1.);       // constant
235              fFitModulation->FixParameter(2, 2.);       // constant
236              fFitModulation->FixParameter(5, 3.);       // constant
237              fFitModulation->SetParameter(7, 0.2);      // v3
238         } break;
239     }
240     switch (fRunModeType) {
241         case kGrid : { fFitModulationOptions += "N0"; } break;
242         default : break;
243     }
244     FillAnalysisSummaryHistogram();
245     return kTRUE;
246 }
247 //_____________________________________________________________________________
248 TH1F* AliAnalysisTaskRhoVnModulation::BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Int_t c, Bool_t append)
249 {
250     // book a TH1F and connect it to the output container
251     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
252     if(fReduceBinsXByFactor > 0 ) bins = TMath::Nint(bins/(double)fReduceBinsXByFactor);
253     if(!fOutputList) return 0x0;
254     TString title(name);
255     if(c!=-1) { // format centrality dependent histograms accordingly
256         name = Form("%s_%i", name, c);
257         title += Form("_%i-%i", fCentralityClasses->At(c), fCentralityClasses->At(1+c));
258     }
259     title += Form(";%s;[counts]", x);
260     TH1F* histogram = new TH1F(name, title.Data(), bins, min, max);
261     histogram->Sumw2();
262     if(append) fOutputList->Add(histogram);
263     return histogram;   
264 }
265 //_____________________________________________________________________________
266 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)
267 {
268     // book a TH2F and connect it to the output container
269     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
270     if(fReduceBinsXByFactor > 0 ) binsx = TMath::Nint(binsx/(double)fReduceBinsXByFactor);
271     if(fReduceBinsYByFactor > 0 ) binsy = TMath::Nint(binsy/(double)fReduceBinsYByFactor);
272     if(!fOutputList) return 0x0;
273     TString title(name);
274     if(c!=-1) { // format centrality dependent histograms accordingly
275         name = Form("%s_%i", name, c);
276         title += Form("_%i-%i", fCentralityClasses->At(c), fCentralityClasses->At(1+c));
277     }
278     title += Form(";%s;%s", x, y);
279     TH2F* histogram = new TH2F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy);
280     histogram->Sumw2();
281     if(append) fOutputList->Add(histogram);
282     return histogram;   
283 }
284 //_____________________________________________________________________________
285 void AliAnalysisTaskRhoVnModulation::UserCreateOutputObjects()
286 {
287     // create output objects
288     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
289     fOutputList = new TList();
290     fOutputList->SetOwner(kTRUE);
291     if(!fCentralityClasses) {   // classes must be defined at this point
292         Int_t c[] = {0, 20, 40, 60, 80, 100};
293         fCentralityClasses = new TArrayI(sizeof(c)/sizeof(c[0]), c);
294     }
295     // global QA
296     fHistCentrality =           BookTH1F("fHistCentrality", "centrality", 102, -2, 100);
297     fHistVertexz =              BookTH1F("fHistVertexz", "vertex z (cm)", 100, -12, 12);
298
299     // pico track kinematics
300     for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) { 
301         fHistPicoTrackPt[i] =          BookTH1F("fHistPicoTrackPt", "p_{t} [GeV/c]", 100, 0, 50, i);
302         if(fFillQAHistograms) {
303             fHistPicoCat1[i] =             BookTH2F("fHistPicoCat1", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
304             fHistPicoCat2[i] =             BookTH2F("fHistPicoCat2", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
305             fHistPicoCat3[i] =             BookTH2F("fHistPicoCat3", "#eta", "#phi", 50, -1, 1, 50, 0, TMath::TwoPi(), i);
306         }
307         // emcal kinematics
308         /* fHistClusterPt[i] =            BookTH1F("fHistClusterPt", "p_{t} [GeV/c]", 100, 0, 100, i); */
309         /* fHistClusterPhi[i] =           BookTH1F("fHistClusterPhi", "#phi", 100, 0, TMath::TwoPi(), i); */
310         /* fHistClusterEta[i] =           BookTH1F("fHistClusterEta", "#eta", 100, -5, 5); */
311
312         // emcal kinematics after hadronic correction
313         /* fHistClusterCorrPt[i] =        BookTH1F("fHistClusterCorrPt", "p_{t} [GeV/c]", 100, 0, 100, i); */
314         /* fHistClusterCorrPhi[i] =       BookTH1F("fHistClusterCorrPhi", "#phi", 100, 0, TMath::TwoPi(), i); */
315         /* fHistClusterCorrEta[i] =       BookTH1F("fHistClusterCorrEta", "#eta", 100, -5, 5, i); */
316     }
317
318     // event plane estimates and quality
319     fHistPsiControl =           new TProfile("fHistPsiControl", "fHistPsiControl", 10, 0, 10);
320     fHistPsiControl->Sumw2();
321     fHistPsiSpread =            new TProfile("fHistPsiSpread", "fHistPsiSpread", 4, 0, 4);
322     fHistPsiSpread->Sumw2();
323     fHistPsiControl->GetXaxis()->SetBinLabel(1, "<#Psi_{2, VZEROA}>");
324     fHistPsiControl->GetXaxis()->SetBinLabel(2, "<#Psi_{2, VZEROC}>");
325     fHistPsiControl->GetXaxis()->SetBinLabel(3, "<#Psi_{2, TPC}>");
326     fHistPsiControl->GetXaxis()->SetBinLabel(4, "<#Psi_{2, TPC, #eta < 0}>");
327     fHistPsiControl->GetXaxis()->SetBinLabel(5, "<#Psi_{2, TPC, #eta > 0}>");
328     fHistPsiControl->GetXaxis()->SetBinLabel(6, "<#Psi_{3, VZEROA}>");
329     fHistPsiControl->GetXaxis()->SetBinLabel(7, "<#Psi_{3, VZEROC}>");
330     fHistPsiControl->GetXaxis()->SetBinLabel(8, "<#Psi_{3, TPC}>");
331     fHistPsiControl->GetXaxis()->SetBinLabel(9, "<#Psi_{3, TPC, #eta < 0}>");
332     fHistPsiControl->GetXaxis()->SetBinLabel(10, "<#Psi_{3, TPC, #eta > 0}>");
333     fHistPsiSpread->GetXaxis()->SetBinLabel(1, "<#Psi_{2, VZEROA} - #Psi_{2, VZEROC}>");
334     fHistPsiSpread->GetXaxis()->SetBinLabel(2, "<#Psi_{2, VZEROC} - #Psi_{2, TPC}>");
335     fHistPsiSpread->GetXaxis()->SetBinLabel(3, "<#Psi_{2, VZEROC} - #Psi_{2, TPC}>");
336     fHistPsiSpread->GetXaxis()->SetBinLabel(4, "<#Psi_{2, TPC, #eta < 0} - #Psi_{2, TPC, #eta > 0}>");
337     fOutputList->Add(fHistPsiControl);
338     fOutputList->Add(fHistPsiSpread);
339     fHistPsiVZEROA =            BookTH1F("fHistPsiVZEROA", "#Psi_{VZEROA}", 100, -.5*TMath::Pi(), .5*TMath::Pi());
340     fHistPsiVZEROC =            BookTH1F("fHistPsiVZEROC", "#Psi_{VZEROC}", 100, -.5*TMath::Pi(), .5*TMath::Pi());
341     fHistPsiTPC =               BookTH1F("fHistPsiTPC", "#Psi_{TPC}", 100, -.5*TMath::Pi(), .5*TMath::Pi());
342     // background
343     for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i ++) {
344         fHistRhoPackage[i] =           BookTH1F("fHistRhoPackage",  "#rho [GeV/c]", 100, 0, 150, i);
345         fHistRho[i] =                  BookTH1F("fHistRho", "#rho [GeV/c]", 100, 0, 150, i);
346     }
347     fHistRhoVsMult =            BookTH2F("fHistRhoVsMult", "multiplicity", "#rho [GeV/c]", 100, 0, 4000, 100, 0, 250);
348     fHistRhoVsCent =            BookTH2F("fHistRhoVsCent", "centrality", "#rho [GeV/c]", 100, 0, 100, 100, 0, 250);
349     fHistRhoAVsMult =           BookTH2F("fHistRhoAVsMult", "multiplicity", "#rho * A (jet) [GeV/c]", 100, 0, 4000, 100, 0, 50);
350     fHistRhoAVsCent =           BookTH2F("fHistRhoAVsCent", "centrality", "#rho * A (jet) [GeV/c]", 100, 0, 100, 100, 0, 50);
351
352     // delta pt distributions
353     for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i ++) {
354         if(fFillQAHistograms)   fHistRCPhiEta[i] = BookTH2F("fHistRCPhiEta", "#phi (RC)", "#eta (RC)", 100, 0, TMath::TwoPi(), 100, -1, 1, i);
355         fHistRhoVsRCPt[i] =            BookTH2F("fHistRhoVsRCPt", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i);
356         fHistRCPt[i] =                 BookTH1F("fHistRCPt", "p_{t} (RC) [GeV/c]", 130, -20, 150, i);
357         if(fFillQAHistograms)   fHistRCPhiEtaExLJ[i] = BookTH2F("fHistRCPhiEtaExLJ", "#phi (RC)", "#eta (RC)", 100, 0, TMath::TwoPi(), 100, -1, 1, i);
358         fHistDeltaPtDeltaPhi2TPC[i] =  BookTH2F("fHistDeltaPtDeltaPhi2TPC", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::Pi(), 100, -50, 100, i);
359         fHistDeltaPtDeltaPhi2V0A[i] =  BookTH2F("fHistDeltaPtDeltaPhi2V0A", "#phi - #Psi_{V0A}", "#delta p_{t} [GeV/c]", 50, 0, TMath::Pi(), 100, -50, 100, i);
360         fHistDeltaPtDeltaPhi2V0C[i] =  BookTH2F("fHistDeltaPtDeltaPhi2V0C", "#phi - #Psi_{V0C}", "#delta p_{t} [GeV/c]", 50, 0, TMath::Pi(), 100, -50, 100, i);
361         fHistDeltaPtDeltaPhi3TPC[i] =  BookTH2F("fHistDeltaPtDeltaPhi3TPC", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::TwoPi()/3., 100, -50, 100, i);
362         fHistDeltaPtDeltaPhi3V0A[i] =  BookTH2F("fHistDeltaPtDeltaPhi3V0A", "#phi - #Psi_{V0A}", "#delta p_{t} [GeV/c]", 50, 0, TMath::TwoPi()/3., 100, -50, 100, i);
363         fHistDeltaPtDeltaPhi3V0C[i] =  BookTH2F("fHistDeltaPtDeltaPhi3V0C", "#phi - #Psi_{V0C}", "#delta p_{t} [GeV/c]", 50, 0, TMath::TwoPi()/3., 100, -50, 100, i);
364         fHistRhoVsRCPtExLJ[i] =        BookTH2F("fHistRhoVsRCPtExLJ", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i);
365         fHistRCPtExLJ[i] =             BookTH1F("fHistRCPtExLJ", "p_{t} (RC) [GeV/c]", 130, -20, 150, i);
366         /* fHistRCPhiEtaRand[i] =         BookTH2F("fHistRCPhiEtaRand", "#phi (RC)", "#eta (RC)", 100, 0, TMath::TwoPi(), 100, -1, 1, i); */
367         fHistDeltaPtDeltaPhi2ExLJTPC[i] = BookTH2F("fHistDeltaPtDeltaPhi2ExLJTPC", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::Pi(), 100, -50, 100, i);
368         fHistDeltaPtDeltaPhi2ExLJV0A[i] = BookTH2F("fHistDeltaPtDeltaPhi2ExLJV0A", "#phi - #Psi_{V0A}", "#delta p_{t} [GeV/c]", 50, 0, TMath::Pi(), 100, -50, 100, i);
369         fHistDeltaPtDeltaPhi2ExLJV0C[i] = BookTH2F("fHistDeltaPtDeltaPhi2ExLJV0C", "#phi - #Psi_{V0C}", "#delta p_{t} [GeV/c]", 50, 0, TMath::Pi(), 100, -50, 100, i);
370         fHistDeltaPtDeltaPhi3ExLJTPC[i] = BookTH2F("fHistDeltaPtDeltaPhi3ExLJTPC", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::TwoPi()/3., 100, -50, 100, i);
371         fHistDeltaPtDeltaPhi3ExLJV0A[i] = BookTH2F("fHistDeltaPtDeltaPhi3ExLJV0A", "#phi - #Psi_{V0A}", "#delta p_{t} [GeV/c]", 50, 0, TMath::TwoPi()/3., 100, -50, 100, i);
372         fHistDeltaPtDeltaPhi3ExLJV0C[i] = BookTH2F("fHistDeltaPtDeltaPhi3ExLJV0C", "#phi - #Psi_{V0C}", "#delta p_{t} [GeV/c]", 50, 0, TMath::TwoPi()/3., 100, -50, 100, i);
373         /* fHistRhoVsRCPtRand[i] =        BookTH2F("fHistRhoVsRCPtRand", "p_{t} (RC) [GeV/c]", "#rho * A (RC) [GeV/c]", 100, 0, 300, 100, 0, 350, i); */
374         /* fHistRCPtRand[i] =             BookTH1F("fHistRCPtRand", "p_{t} (RC) [GeV/c]", 130, -20, 150, i); */
375         /* fHistDeltaPtDeltaPhi2Rand[i] =  BookTH2F("fHistDeltaPtDeltaPhi2Rand", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::Pi(), 100, -50, 100, i); */
376         /* fHistDeltaPtDeltaPhi3Rand[i] =  BookTH2F("fHistDeltaPtDeltaPhi3Rand", "#phi - #Psi_{TPC}", "#delta p_{t} [GeV/c]", 50, 0, TMath::TwoPi()/3., 100, -50, 100, i); */
377         // jet histograms (after kinematic cuts)
378         fHistJetPtRaw[i] =             BookTH1F("fHistJetPtRaw", "p_{t} RAW [GeV/c]", 200, -50, 150, i);
379         fHistJetPt[i] =                BookTH1F("fHistJetPt", "p_{t} [GeV/c]", 350, -100, 250, i);
380         if(fFillQAHistograms)   fHistJetEtaPhi[i] =            BookTH2F("fHistJetEtaPhi", "#eta", "#phi", 100, -1, 1, 100, 0, TMath::TwoPi(), i);
381         fHistJetPtArea[i] =            BookTH2F("fHistJetPtArea", "p_{t} [GeV/c]", "Area", 175, -100, 250, 30, 0, 0.9, i);
382         fHistJetPtConstituents[i] =    BookTH2F("fHistJetPtConstituents", "p_{t} [GeV/c]", "Area", 350, -100, 250, 60, 0, 150, i);
383         fHistJetEtaRho[i] =            BookTH2F("fHistJetEtaRho", "#eta", "#rho", 100, -1, 1, 100, 0, 300, i);
384         // in plane and out of plane spectra
385         fHistJetPsiTPCPt[i] =          BookTH2F("fHistJetPsiTPCPt", "#phi_{jet} - #Psi_{2, TPC}", "p_{t} [GeV/c]", 50, 0., TMath::Pi(), 700, -100, 250, i);
386         fHistJetPsiVZEROAPt[i] =       BookTH2F("fHistJetPsiVZEROAPt", "#phi_{jet} - #Psi_{2, VZEROA}", "p_{t} [GeV/c]", 50, 0., TMath::Pi(), 700, -100, 250, i);
387         fHistJetPsiVZEROCPt[i] =       BookTH2F("fHistJetPsiVZEROCPt", "#phi_{jet} - #Psi_{V2, ZEROC}", "p_{t} [GeV/c]", 50, 0., TMath::Pi(), 700, -100, 250, i);
388         // phi minus psi
389         fHistDeltaPhi2VZEROA[i] =       BookTH1F("fHistDeltaPhi2VZEROA", "#phi_{jet} - #Psi_{2, VZEROA}", 50, 0, TMath::Pi(), i);
390         fHistDeltaPhi2VZEROC[i] =       BookTH1F("fHistDeltaPhi2VZEROC", "#phi_{jet} - #Psi_{2, VZEROC}", 50, 0, TMath::Pi(), i);
391         fHistDeltaPhi2TPC[i] =          BookTH1F("fHistDeltaPhi2TPC", "#phi_{jet} - #Psi_{2, TPC}", 50, 0, TMath::Pi(), i);
392         fHistDeltaPhi3VZEROA[i] =       BookTH1F("fHistDeltaPhi3VZEROA", "#phi_{jet} - #Psi_{2, VZEROA}", 50, 0, TMath::TwoPi()/3., i);
393         fHistDeltaPhi3VZEROC[i] =       BookTH1F("fHistDeltaPhi3VZEROC", "#phi_{jet} - #Psi_{2, VZEROC}", 50, 0, TMath::TwoPi()/3., i);
394         fHistDeltaPhi3TPC[i] =          BookTH1F("fHistDeltaPhi3TPC", "#phi_{jet} - #Psi_{2, TPC}", 50, 0, TMath::TwoPi()/3., i);
395
396         fProfV2Resolution[i] = new TProfile(Form("fProfV2Resolution_%i", i), Form("fProfV2Resolution_%i", i), 8, -0.5, 7.5);
397         fProfV2Resolution[i]->GetXaxis()->SetBinLabel(3, "<cos(2(#Psi_{VZEROA} - #Psi_{VZEROC}))>");
398         fProfV2Resolution[i]->GetXaxis()->SetBinLabel(4, "<cos(2(#Psi_{VZEROC} - #Psi_{VZEROA}))>");
399         fProfV2Resolution[i]->GetXaxis()->SetBinLabel(5, "<cos(2(#Psi_{VZEROA} - #Psi_{TPC}))>");
400         fProfV2Resolution[i]->GetXaxis()->SetBinLabel(6, "<cos(2(#Psi_{TPC} - #Psi_{VZEROA}))>");
401         fProfV2Resolution[i]->GetXaxis()->SetBinLabel(7, "<cos(2(#Psi_{VZEROC} - #Psi_{TPC}))>");
402         fProfV2Resolution[i]->GetXaxis()->SetBinLabel(8, "<cos(2(#Psi_{TPC} - #Psi_{VZEROC}))>");
403         fOutputList->Add(fProfV2Resolution[i]); 
404         fProfV3Resolution[i] = new TProfile(Form("fProfV3Resolution_%i", i), Form("fProfV3Resolution_%i", i), 8, -0.5, 7.5);
405         fProfV3Resolution[i]->GetXaxis()->SetBinLabel(3, "<cos(3(#Psi_{VZEROA} - #Psi_{VZEROC}))>");
406         fProfV3Resolution[i]->GetXaxis()->SetBinLabel(4, "<cos(3(#Psi_{VZEROC} - #Psi_{VZEROA}))>");
407         fProfV3Resolution[i]->GetXaxis()->SetBinLabel(5, "<cos(3(#Psi_{VZEROA} - #Psi_{TPC}))>");
408         fProfV3Resolution[i]->GetXaxis()->SetBinLabel(6, "<cos(3(#Psi_{TPC} - #Psi_{VZEROA}))>");
409         fProfV3Resolution[i]->GetXaxis()->SetBinLabel(7, "<cos(3(#Psi_{VZEROC} - #Psi_{TPC}))>");
410         fProfV3Resolution[i]->GetXaxis()->SetBinLabel(8, "<cos(3(#Psi_{TPC} - #Psi_{VZEROC}))>");
411         fOutputList->Add(fProfV3Resolution[i]); 
412     }
413     // cdf and pdf of chisquare distribution
414     fHistPvaluePDF = BookTH1F("fHistPvaluePDF", "PDF #chi^{2}", 500, 0, 1);
415     fHistPvalueCDF = BookTH1F("fHistPvalueCDF", "CDF #chi^{2}", 500, 0, 1);
416     // vn profile
417     Float_t temp[fCentralityClasses->GetSize()];
418     for(Int_t i(0); i < fCentralityClasses->GetSize(); i++) temp[i] = fCentralityClasses->At(i);
419     fProfV2 = new TProfile("fProfV2", "fProfV2", fCentralityClasses->GetSize()-1, temp);
420     fProfV3 = new TProfile("fProfV3", "fProfV3", fCentralityClasses->GetSize()-1, temp);
421     fOutputList->Add(fProfV2);
422     fOutputList->Add(fProfV3);
423     switch (fFitModulationType) {
424         case kQC2 : {
425             fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
426             fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
427             fOutputList->Add(fProfV2Cumulant);
428             fOutputList->Add(fProfV3Cumulant);
429         } break;
430         case kQC4 : {
431             fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
432             fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
433             fOutputList->Add(fProfV2Cumulant);
434             fOutputList->Add(fProfV3Cumulant);
435         } break;
436         default : break;
437     }
438     // for the histograms initialized below, binning is fixed to runnumbers or flags
439     fReduceBinsXByFactor = 1;
440     fReduceBinsYByFactor = 1;
441     if(fFillQAHistograms) {
442         fHistRunnumbersEta = new TH2F("fHistRunnumbersEta", "fHistRunnumbersEta", 100, -.5, 99.5, 100, -1.1, 1.1);
443         fHistRunnumbersEta->Sumw2();
444         fOutputList->Add(fHistRunnumbersEta);
445         fHistRunnumbersPhi = new TH2F("fHistRunnumbersPhi", "fHistRunnumbersPhi", 100, -.5, 99.5, 100, -0.2, TMath::TwoPi()+0.2);
446         fHistRunnumbersPhi->Sumw2();
447         fOutputList->Add(fHistRunnumbersPhi);
448     }
449     fHistAnalysisSummary = BookTH1F("fHistAnalysisSummary", "flag", 50, -0.5, 50.5);
450     fHistSwap = new TH1F("fHistSwap", "fHistSwap", 20, 0, TMath::TwoPi());
451     if(fUsePtWeight) fHistSwap->Sumw2();
452
453     if(fUserSuppliedV2) fOutputList->Add(fUserSuppliedV2);
454     if(fUserSuppliedV3) fOutputList->Add(fUserSuppliedV3);
455     if(fUserSuppliedR2) fOutputList->Add(fUserSuppliedR2);
456     if(fUserSuppliedR3) fOutputList->Add(fUserSuppliedR3);
457     // increase readability of output list
458     fOutputList->Sort();
459     PostData(1, fOutputList);
460
461     switch (fRunModeType) {
462         case kLocal : {
463             fOutputListGood = new TList();
464             fOutputListGood->SetOwner(kTRUE);
465             fOutputListBad = new TList();
466             fOutputListBad->SetOwner(kTRUE);
467             PostData(2, fOutputListGood);
468             PostData(3, fOutputListBad);
469         } break;
470         default: break;
471     }
472 }
473 //_____________________________________________________________________________
474 Bool_t AliAnalysisTaskRhoVnModulation::Run()
475 {
476     // user exec: execute once for each event
477     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
478     if(!fInitialized) fInitialized = InitializeAnalysis();
479     // reject the event if expected data is missing
480     if(!PassesCuts(InputEvent())) return kFALSE;
481     if(!(fTracks||fJets||fRho)) return kFALSE;
482     if(!fCaloClusters && fDebug > 0) printf(" > Warning: couldn't retreive calo clusters! < \n");
483     // [0][0] psi2a     [1,0]   psi2c
484     // [0][1] psi3a     [1,1]   psi3c
485     Double_t vzero[2][2];
486     CalculateEventPlaneVZERO(vzero);
487     // [0] psi2         [1] psi3
488     Double_t tpc[2];
489     CalculateEventPlaneTPC(tpc);
490     Double_t psi2(-1), psi3(-1);
491     // arrays which will hold the fit parameters
492     switch (fDetectorType) {    // determine the detector type for the rho fit
493         case kTPC :     { psi2 = tpc[0];        psi3 = tpc[1]; }        break;
494         case kVZEROA :  { psi2 = vzero[0][0];   psi3 = vzero[0][1]; }   break;  
495         case kVZEROC :  { psi2 = vzero[1][0];   psi3 = vzero[1][1]; }   break;
496         default : break;
497     }
498     switch (fFitModulationType) { // do the fits
499         case kNoFit : { fFitModulation->FixParameter(0, RhoVal()); } break;
500         case kV2 : {
501             if(CorrectRho(psi2, psi3)) {
502                 fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
503                 if(fUserSuppliedR2) {
504                     Double_t r(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
505                     if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
506                 }
507                 CalculateEventPlaneResolution(vzero, tpc);
508             }
509         } break;
510         case kV3 : {
511             if(CorrectRho(psi2, psi3)) {
512                 if(fUserSuppliedR3) {
513                     Double_t r(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
514                     if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
515                 }
516                 fProfV3->Fill(fCent, fFitModulation->GetParameter(3));
517                 CalculateEventPlaneResolution(vzero, tpc);
518             }
519         } break;
520         case kQC2 : {
521             Bool_t QC2(CorrectRho(psi2, psi3));
522             if(fUserSuppliedR2 && fUserSuppliedR3) {
523                 // note for the qc method, resolution is REVERSED to go back to v2obs
524                 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
525                 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
526                 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
527                 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(3)*r3);
528             }
529             if (fUsePtWeight) { // use weighted weights
530                 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), QCnM11());
531                 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), QCnM11());
532                 if(QC2) {        // how to deal with negative results from the cumulants ?
533                     fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5), QCnM11());
534                     fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5), QCnM11()); 
535                 }
536             } else {
537                 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), QCnM()*(QCnM()-1));
538                 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), QCnM()*(QCnM()-1));
539                 if(QC2) {
540                     fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5), QCnM()*(QCnM()-1));
541                     fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5), QCnM()*(QCnM()-1));
542                 }
543             }
544             CalculateEventPlaneResolution(vzero, tpc);
545         } break;
546         case kQC4 : {
547             Bool_t QC4(CorrectRho(psi2, psi3));
548             if(fUserSuppliedR2 && fUserSuppliedR3) {
549                 // note for the qc method, resolution is REVERSED to go back to v2obs
550                 Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
551                 Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
552                 if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
553                 if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(3)*r3);
554             }
555             if (fUsePtWeight) { // use weighted weights
556                 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
557                 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
558                 if(QC4) {        // how to deal with negative values of cumulants ?
559                     fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM1111()*/);
560                     fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM1111()*/); 
561                 }
562             } else {
563                 fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
564                 fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*,  QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
565                 if(QC4) {
566                     fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.25)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
567                     fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.25)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
568                 }
569             }
570             CalculateEventPlaneResolution(vzero, tpc);
571         } break;
572         default : {
573             if(CorrectRho(psi2, psi3)) {
574                 if(fUserSuppliedR2 && fUserSuppliedR3) {
575                     Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
576                     Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
577                     if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r2);
578                     if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(3)/r3);
579                 }
580                 fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
581                 fProfV3->Fill(fCent, fFitModulation->GetParameter(7));
582                 CalculateEventPlaneResolution(vzero, tpc);
583             }
584         } break;
585     }
586     // fill a number of histograms 
587     FillHistogramsAfterSubtraction(vzero, tpc);
588     // send the output to the connected output container
589     PostData(1, fOutputList);
590     switch (fRunModeType) {
591         case kLocal : {
592             PostData(2, fOutputListGood);
593             PostData(3, fOutputListBad);
594         } break;
595         default: break;
596     }
597     return kTRUE;
598 }
599 //_____________________________________________________________________________
600 void AliAnalysisTaskRhoVnModulation::CalculateEventPlaneVZERO(Double_t vzero[2][2]) const 
601 {
602     // get the vzero event plane
603     if(fUseV0EventPlaneFromHeader) {    // use the vzero from the header
604         Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0), h(0);
605         vzero[0][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 2, a, b);
606         vzero[1][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 2, c, d);
607         vzero[0][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 3, e, f);
608         vzero[1][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 3, g, h);
609         return;
610     }
611     // grab the vzero event plane without recentering
612     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
613     Double_t qxa2(0), qya2(0), qxc2(0), qyc2(0);    // for psi2
614     Double_t qxa3(0), qya3(0), qxc3(0), qyc3(0);    // for psi3
615     for(Int_t iVZERO(0); iVZERO < 64; iVZERO++) {
616         Double_t phi(TMath::PiOver4()*(.5+iVZERO%8)), /* eta(0), */ weight(InputEvent()->GetVZEROEqMultiplicity(iVZERO));
617 //        (iVZERO<32) ? eta = -3.45+.5*(iVZERO/8) : eta = 4.8-.6*((iVZERO/8)-4);
618         if(iVZERO<32) {
619             qxa2 += weight*TMath::Cos(2.*phi);
620             qya2 += weight*TMath::Sin(2.*phi);
621             qxa3 += weight*TMath::Cos(3.*phi);
622             qya3 += weight*TMath::Sin(3.*phi);
623         }
624         else {
625             qxc2 += weight*TMath::Cos(2.*phi);
626             qyc2 += weight*TMath::Sin(2.*phi);
627             qxc3 += weight*TMath::Cos(3.*phi);
628             qyc3 += weight*TMath::Sin(3.*phi);
629        }
630     }
631     vzero[0][0] = .5*TMath::ATan2(qya2, qxa2);
632     vzero[1][0] = .5*TMath::ATan2(qyc2, qxc2);
633     vzero[0][1] = (1./3.)*TMath::ATan2(qya3, qxa3);
634     vzero[1][1] = (1./3.)*TMath::ATan2(qyc3, qxc3);
635 }
636 //_____________________________________________________________________________
637 void AliAnalysisTaskRhoVnModulation::CalculateEventPlaneTPC(Double_t* tpc)
638 {
639    // grab the TPC event plane
640    if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
641    fNAcceptedTracks = 0;                // reset the track counter
642    Double_t qx2(0), qy2(0);     // for psi2
643    Double_t qx3(0), qy3(0);     // for psi3
644    if(fTracks) {
645        Float_t excludeInEta[] = {-999, -999};
646        if(fExcludeLeadingJetsFromFit > 0 ) {    // remove the leading jet from ep estimate
647            AliEmcalJet* leadingJet[] = {0x0, 0x0};
648            static Int_t lJets[9999] = {-1};
649            GetSortedArray(lJets, fJets);
650            for(Int_t i(0); i < fJets->GetEntriesFast(); i++) {     // get the two leading jets
651                if (1 + i > fJets->GetEntriesFast()) break;
652                leadingJet[0] = static_cast<AliEmcalJet*>(fJets->At(lJets[i]));
653                leadingJet[1] = static_cast<AliEmcalJet*>(fJets->At(lJets[i+1]));
654                if(PassesCuts(leadingJet[0]) && PassesCuts(leadingJet[1])) break;
655            }
656            if(leadingJet[0] && leadingJet[1]) {
657                for(Int_t i(0); i < 2; i++) excludeInEta[i] = leadingJet[i]->Eta();
658            }
659        }
660        Int_t iTracks(fTracks->GetEntriesFast());
661        for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
662            AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
663            if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
664            if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta[0]) < fJetRadius*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - fJetRadius - fJetMaxEta ) > 0 )) continue;
665            fNAcceptedTracks++;
666            qx2+= TMath::Cos(2.*track->Phi());
667            qy2+= TMath::Sin(2.*track->Phi());
668            qx3+= TMath::Cos(3.*track->Phi());
669            qy3+= TMath::Sin(3.*track->Phi());
670        }
671    }
672    tpc[0] = .5*TMath::ATan2(qy2, qx2);
673    tpc[1] = (1./3.)*TMath::ATan2(qy3, qx3);
674
675 //_____________________________________________________________________________
676 void AliAnalysisTaskRhoVnModulation::CalculateEventPlaneResolution(Double_t vzero[2][2], Double_t* tpc) const
677 {
678     // fill the profiles for the resolution parameters
679     if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
680     fProfV2Resolution[fInCentralitySelection]->Fill(2., TMath::Cos(2.*(vzero[0][0] - vzero[1][0])));
681     fProfV2Resolution[fInCentralitySelection]->Fill(3., TMath::Cos(2.*(vzero[1][0] - vzero[0][0])));
682     fProfV2Resolution[fInCentralitySelection]->Fill(4., TMath::Cos(2.*(vzero[0][0] - tpc[0])));
683     fProfV2Resolution[fInCentralitySelection]->Fill(5., TMath::Cos(2.*(tpc[0] - vzero[0][0])));
684     fProfV2Resolution[fInCentralitySelection]->Fill(6., TMath::Cos(2.*(vzero[1][0] - tpc[0])));
685     fProfV2Resolution[fInCentralitySelection]->Fill(7., TMath::Cos(2.*(tpc[0] - vzero[1][0])));
686     fProfV3Resolution[fInCentralitySelection]->Fill(2., TMath::Cos(3.*(vzero[0][0] - vzero[1][0])));
687     fProfV3Resolution[fInCentralitySelection]->Fill(3., TMath::Cos(3.*(vzero[1][0] - vzero[0][0])));
688     fProfV3Resolution[fInCentralitySelection]->Fill(4., TMath::Cos(3.*(vzero[0][0] - tpc[0])));
689     fProfV3Resolution[fInCentralitySelection]->Fill(5., TMath::Cos(3.*(tpc[0] - vzero[0][0])));
690     fProfV3Resolution[fInCentralitySelection]->Fill(6., TMath::Cos(3.*(vzero[1][0] - tpc[0])));
691     fProfV3Resolution[fInCentralitySelection]->Fill(7., TMath::Cos(3.*(tpc[0] - vzero[1][0])));
692 }
693 //_____________________________________________________________________________
694 void AliAnalysisTaskRhoVnModulation::CalculateRandomCone(Float_t &pt, Float_t &eta, Float_t &phi, 
695         AliEmcalJet* jet, Bool_t randomize) const
696 {
697     // get a random cone
698     if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
699     pt = 0; eta = 0; phi = 0;
700     Float_t etaJet(999), phiJet(999), dJet(999);        // no jet: same as jet very far away
701     if(jet) { // if a leading jet is given, use its kinematic properties
702         etaJet = jet->Eta();
703         phiJet = jet->Phi();
704     }
705     // force the random cones to at least be within detector acceptance
706     Float_t minPhi(fJetMinPhi), maxPhi(fJetMaxPhi);
707     if(maxPhi > TMath::TwoPi()) maxPhi = TMath::TwoPi();
708     if(minPhi < 0 ) minPhi = 0;
709     Float_t diffRcRJR(TMath::Abs(fRandomConeRadius-fJetRadius));
710     // construct a random cone and see if it's far away enough from the leading jet
711     Int_t attempts(1000);
712     while(kTRUE) {
713         attempts--;
714         eta = gRandom->Uniform(fJetMinEta+diffRcRJR, fJetMaxEta-diffRcRJR);
715         phi = gRandom->Uniform(minPhi, maxPhi);
716
717         dJet = TMath::Sqrt((etaJet-eta)*(etaJet-eta)+(phiJet-phi)*(phiJet-phi));
718         if(dJet > fMinDisanceRCtoLJ) break;
719         else if (attempts == 0) {
720             printf(" > No random cone after 1000 tries, giving up ... !\n");
721             return;
722         }
723     }
724     if(fTracks) {
725         Int_t iTracks(fTracks->GetEntriesFast());
726         for(Int_t i(0); i < iTracks; i++) {
727             AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
728             if(!PassesCuts(track)) continue;
729             Float_t etaTrack(track->Eta()), phiTrack(track->Phi()), ptTrack(track->Pt());
730             // if requested, randomize eta and phi to destroy any correlated fluctuations
731             if(randomize) {
732                 etaTrack = gRandom->Uniform(fTrackMinEta, fTrackMaxEta);
733                 phiTrack = gRandom->Uniform(minPhi, maxPhi);
734             }
735             // get distance from cone
736             if(TMath::Abs(phiTrack-phi) > TMath::Abs(phiTrack - phi + TMath::TwoPi())) phiTrack+=TMath::TwoPi();
737             if(TMath::Abs(phiTrack-phi) > TMath::Abs(phiTrack - phi - TMath::TwoPi())) phiTrack-=TMath::TwoPi();
738             if(TMath::Sqrt(TMath::Abs((etaTrack-eta)*(etaTrack-eta)+(phiTrack-phi)*(phiTrack-phi))) <= fRandomConeRadius) pt+=ptTrack;
739         }
740     }
741 }
742 //_____________________________________________________________________________
743 Double_t AliAnalysisTaskRhoVnModulation::CalculateQC2(Int_t harm) {
744     // get the second order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho
745     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
746     Double_t reQ(0), imQ(0), modQ(0), M11(0), M(0);
747     if(fUsePtWeight) {  // for the weighted 2-nd order q-cumulant
748         QCnQnk(harm, 1, reQ, imQ);      // get the weighted 2-nd order q-vectors
749         modQ = reQ*reQ+imQ*imQ;         // get abs Q-squared
750         M11 = QCnM11();                 // equals S2,1 - S1,2
751         return (M11 > 0) ? ((modQ - QCnS(1,2))/M11) : -999;
752     } // else return the non-weighted 2-nd order q-cumulant
753     QCnQnk(harm, 0, reQ, imQ);          // get the non-weighted 2-nd order q-vectors
754     modQ = reQ*reQ+imQ*imQ;             // get abs Q-squared
755     M = QCnM();
756     return (M > 1) ? (modQ - M)/(M*(M-1)) : -999;
757 }
758 //_____________________________________________________________________________
759 Double_t AliAnalysisTaskRhoVnModulation::CalculateQC4(Int_t harm) {
760     // get the fourth order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho
761     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
762     Double_t reQn1(0), imQn1(0), reQ2n2(0), imQ2n2(0), reQn3(0), imQn3(0), M1111(0), M(0);
763     Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0);  // terms of the calculation
764     if(fUsePtWeight) {  // for the weighted 4-th order q-cumulant
765         QCnQnk(harm, 1, reQn1, imQn1);
766         QCnQnk(harm*2, 2, reQ2n2, imQ2n2);
767         QCnQnk(harm, 3, reQn3, imQn3);
768         // fill in the terms ...
769         a = (reQn1*reQn1+imQn1*imQn1)*(reQn1*reQn1+imQn1*imQn1);
770         b = reQ2n2*reQ2n2 + imQ2n2*imQ2n2;
771         c = -2.*(reQ2n2*reQn1*reQn1-reQ2n2*imQn1*imQn1+2.*imQ2n2*reQn1*imQn1);
772         d = 8.*(reQn3*reQn1+imQn3*imQn1);
773         e = -4.*QCnS(1,2)*(reQn1*reQn1+imQn1*imQn1);
774         f = -6.*QCnS(1,4);
775         g = 2.*QCnS(2,2);
776         M1111 = QCnM1111();
777         return (M1111 > 0) ? (a+b+c+d+e+f+g)/M1111 : -999;
778     }   // else return the unweighted case
779     Double_t reQn(0), imQn(0), reQ2n(0), imQ2n(0);
780     QCnQnk(harm, 0, reQn, imQn);
781     QCnQnk(harm*2, 0, reQ2n, imQ2n);
782     // fill in the terms ...
783     M = QCnM();
784     if(M < 4) return -999;
785     a = (reQn*reQn+imQn*imQn)*(reQn*reQn+imQn*imQn);
786     b = reQ2n*reQ2n + imQ2n*imQ2n;
787     c = -2.*(reQ2n*reQn*reQn-reQ2n*imQn*imQn+2.*imQ2n*reQn*imQn);
788     e = -4.*(M-2)*(reQn*reQn+imQn*imQn);
789     f = 2.*M*(M-3);
790     return (a+b+c+e+f)/(M*(M-1)*(M-2)*(M-3));
791 }
792 //_____________________________________________________________________________
793 void AliAnalysisTaskRhoVnModulation::QCnQnk(Int_t n, Int_t k, Double_t &reQ, Double_t &imQ) {
794     // get the weighted n-th order q-vector, pass real and imaginary part as reference
795     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
796     if(!fTracks) return;
797     fNAcceptedTracksQCn = 0;
798     Int_t iTracks(fTracks->GetEntriesFast());
799     for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
800         AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
801         if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
802         fNAcceptedTracksQCn++;
803         // for the unweighted case, k equals zero and the weight doesn't contribute to the equation below
804         reQ += TMath::Power(track->Pt(), k) * TMath::Cos(((double)n)*track->Phi());
805         imQ += TMath::Power(track->Pt(), k) * TMath::Sin(((double)n)*track->Phi());
806     }
807 }
808 //_____________________________________________________________________________
809 void AliAnalysisTaskRhoVnModulation::QCnDiffentialFlowVectors(
810         TClonesArray* pois, TArrayD* ptBins, Bool_t vpart, Double_t* repn, Double_t* impn, 
811         Double_t *mp, Double_t *reqn, Double_t *imqn, Double_t* mq, Int_t n) 
812 {
813     // get  unweighted differential flow vectors
814     Int_t iPois(pois->GetEntriesFast());
815     if(vpart) {
816         for(Int_t i(0); i < iPois; i++) {
817             for(Int_t ptBin(0); ptBin < ptBins->GetSize()-1; ptBin++) {
818                 AliVTrack* poi = static_cast<AliVTrack*>(pois->At(i));
819                 if(PassesCuts(poi)) {
820                     if(poi->Pt() >= ptBins->At(ptBin) && poi->Pt() < ptBins->At(ptBin+1)) {
821                             // fill the flow vectors assuming that all poi's are in the rp selection (true by design)
822                             repn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
823                             impn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
824                             mp[ptBin]++;
825                             reqn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
826                             imqn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
827                             mq[ptBin]++;
828                     }
829                 }
830             }
831         }
832     } else {
833         for(Int_t i(0); i < iPois; i++) {
834             for(Int_t ptBin(0); ptBin < ptBins->GetSize()-1; ptBin++) {
835                 AliEmcalJet* poi = static_cast<AliEmcalJet*>(pois->At(i));
836                 if(poi && poi->PtSub() > 0) {   // note here that no cuts are needed since only accepted jets have PtSub set !    
837                     if(poi->PtSub() >= ptBins->At(ptBin) && poi->PtSub() < ptBins->At(ptBin+1)) {    
838                             // fill the flow vectors assuming that all poi's are in the rp selection (true by design)  
839                             repn[ptBin]+=TMath::Cos(((double)n)*poi->Phi());
840                             impn[ptBin]+=TMath::Sin(((double)n)*poi->Phi());
841                             mp[ptBin]++;        // qn isn't filled, no overlap between poi's and rp's
842                     }
843                 }
844             }
845         }
846     }
847 }
848 //_____________________________________________________________________________
849 Double_t AliAnalysisTaskRhoVnModulation::QCnS(Int_t i, Int_t j) {
850     // get the weighted ij-th order autocorrelation correction
851     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
852     if(!fTracks || i <= 0 || j <= 0) return -999;
853     Int_t iTracks(fTracks->GetEntriesFast());
854     Double_t Sij(0);
855     for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
856         AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
857         if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
858         Sij+=TMath::Power(track->Pt(), j);
859     }
860     return TMath::Power(Sij, i);
861 }
862 //_____________________________________________________________________________
863 Double_t AliAnalysisTaskRhoVnModulation::QCnM() {
864     // get multiplicity for unweighted q-cumulants. function QCnQnk should be called first
865     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
866     return (Double_t) fNAcceptedTracksQCn;
867 }
868 //_____________________________________________________________________________
869 Double_t AliAnalysisTaskRhoVnModulation::QCnM11() {
870     // get multiplicity weights for the weighted two particle cumulant
871     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
872     return (QCnS(2,1) - QCnS(1,2));
873 }
874 //_____________________________________________________________________________
875 Double_t AliAnalysisTaskRhoVnModulation::QCnM1111() {
876     // get multiplicity weights for the weighted four particle cumulant
877     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
878     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));
879 }
880 //_____________________________________________________________________________
881 Bool_t AliAnalysisTaskRhoVnModulation::CorrectRho(Double_t psi2, Double_t psi3) 
882 {
883     // get rho' -> rho(phi)
884     // two routines are available, both can be used with or without pt weights
885     //  [1] get vn from q-cumulants or as an integrated value from a user supplied histogram
886     //      in case of cumulants, both cumulants and vn values are stored. in both cases, v2 and v3
887     //      are expected. a check is performed to see if rho has no negative local minimum
888     //      for full description, see Phys. Rev. C 83, 044913
889     //  [2] fitting a fourier expansion to the de/dphi distribution
890     //      the fit can be done with either v2, v3 or a combination.
891     //      in all cases, a cut can be made on the p-value of the chi-squared value of the fit
892     //      and a check can be performed to see if rho has no negative local minimum
893     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
894     switch (fFitModulationType) {       // for approaches where no fitting is required
895         case kQC2 : {
896             fFitModulation->FixParameter(4, psi2); 
897             fFitModulation->FixParameter(6, psi3);
898             fFitModulation->FixParameter(3, CalculateQC2(2));
899             fFitModulation->FixParameter(7, CalculateQC2(3));
900             if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { 
901                 fFitModulation->SetParameter(7, 0);
902                 fFitModulation->SetParameter(3, 0);
903                 fFitModulation->SetParameter(0, RhoVal());
904                 return kFALSE;
905             }
906             return (fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) ? kTRUE : kFALSE;
907         } break;
908         case kQC4 : {
909             fFitModulation->FixParameter(4, psi2); 
910             fFitModulation->FixParameter(6, psi3);
911             fFitModulation->FixParameter(3, CalculateQC4(2));
912             fFitModulation->FixParameter(7, CalculateQC4(3));
913             if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { 
914                 fFitModulation->SetParameter(7, 0);
915                 fFitModulation->SetParameter(3, 0);
916                 fFitModulation->SetParameter(0, RhoVal());
917                 return kFALSE;
918             }
919             return (fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) ? kTRUE : kFALSE;
920         } break;
921         case kIntegratedFlow : {
922             // use v2 and v3 values from an earlier iteration over the data
923             fFitModulation->FixParameter(3, fUserSuppliedV2->GetBinContent(fUserSuppliedV2->GetXaxis()->FindBin(fCent)));
924             fFitModulation->FixParameter(4, psi2);
925             fFitModulation->FixParameter(6, psi3);
926             fFitModulation->FixParameter(7, fUserSuppliedV3->GetBinContent(fUserSuppliedV3->GetXaxis()->FindBin(fCent)));
927             if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { 
928                 fFitModulation->SetParameter(7, 0);
929                 fFitModulation->SetParameter(3, 0);
930                 fFitModulation->SetParameter(0, RhoVal());
931                 return kFALSE;
932             }
933             return kTRUE;
934         }
935         default : break;
936     }
937     TString detector("");
938     switch (fDetectorType) {
939         case kTPC : detector+="TPC";
940             break;
941         case kVZEROA : detector+="VZEROA";
942             break;
943         case kVZEROC : detector+="VZEROC";
944             break;
945         default: break;
946     }
947     Int_t iTracks(fTracks->GetEntriesFast());
948     Double_t excludeInEta[] = {-999, -999};
949     Double_t excludeInPhi[] = {-999, -999};
950     Double_t excludeInPt[]  = {-999, -999};
951     if(iTracks <= 0 || RhoVal() <= 0 ) return kFALSE;   // no use fitting an empty event ...
952     if(fExcludeLeadingJetsFromFit > 0 ) {
953         AliEmcalJet* leadingJet[] = {0x0, 0x0};
954         static Int_t lJets[9999] = {-1};
955         GetSortedArray(lJets, fJets);
956         for(Int_t i(0); i < fJets->GetEntriesFast(); i++) {     // get the two leading jets
957             if (1 + i > fJets->GetEntriesFast()) break;
958             leadingJet[0] = static_cast<AliEmcalJet*>(fJets->At(lJets[i]));
959             leadingJet[1] = static_cast<AliEmcalJet*>(fJets->At(lJets[i+1]));
960             if(PassesCuts(leadingJet[0]) && PassesCuts(leadingJet[1])) break;
961         }
962         if(leadingJet[0] && leadingJet[1]) {
963             for(Int_t i(0); i < 2; i++) {
964                 excludeInEta[i] = leadingJet[i]->Eta();
965                 excludeInPhi[i] = leadingJet[i]->Phi();
966                 excludeInPt[i]  = leadingJet[i]->Pt();
967             }
968         }
969     }
970     fHistSwap->Reset();                 // clear the histogram
971     TH1F _tempSwap;
972     if(fRebinSwapHistoOnTheFly) {
973         if(fNAcceptedTracks < 49) fNAcceptedTracks = 49;       // avoid aliasing effects
974         _tempSwap = TH1F("_tempSwap", "_tempSwap", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), 0, TMath::TwoPi());
975     }
976     else _tempSwap = *fHistSwap;         // now _tempSwap holds the desired histo
977     for(Int_t i(0); i < iTracks; i++) {
978             AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
979             if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta[0]) < fJetRadius*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - fJetRadius - fJetMaxEta ) > 0 )) continue;
980             if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
981             if(fUsePtWeight) _tempSwap.Fill(track->Phi(), track->Pt());
982             else _tempSwap.Fill(track->Phi());
983     }
984 //    for(Int_t i(0); i < _tempSwap.GetXaxis()->GetNbins(); i++) _tempSwap.SetBinError(1+i, TMath::Sqrt(_tempSwap.GetBinContent(1+i)));
985     fFitModulation->SetParameter(0, RhoVal());
986     switch (fFitModulationType) {
987         case kNoFit : { fFitModulation->FixParameter(0, RhoVal() ); 
988         } break;
989         case kV2 : { 
990             fFitModulation->FixParameter(4, psi2); 
991         } break;
992         case kV3 : { 
993             fFitModulation->FixParameter(4, psi3); 
994         } break;
995         case kCombined : { 
996             fFitModulation->FixParameter(4, psi2); 
997             fFitModulation->FixParameter(6, psi3);
998         } break;
999         case kFourierSeries : {
1000             // in this approach, an explicit calculation will be made of vn = sqrt(xn^2+yn^2)
1001             // where x[y] = Integrate[r(phi)cos[sin](n phi)dphi, 0, 2pi]
1002             Double_t cos2(0), sin2(0), cos3(0), sin3(0), sumPt(0);
1003             for(Int_t i(0); i < iTracks; i++) {
1004                 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1005                 if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
1006                 sumPt += track->Pt();
1007                 cos2 += track->Pt()*TMath::Cos(2*PhaseShift(track->Phi()-psi2)); 
1008                 sin2 += track->Pt()*TMath::Sin(2*PhaseShift(track->Phi()-psi2));
1009                 cos3 += track->Pt()*TMath::Cos(3*PhaseShift(track->Phi()-psi3)); 
1010                 sin3 += track->Pt()*TMath::Sin(3*PhaseShift(track->Phi()-psi3));
1011             }
1012             fFitModulation->SetParameter(3, TMath::Sqrt(cos2*cos2+sin2*sin2)/RhoVal());
1013             fFitModulation->SetParameter(4, psi2);
1014             fFitModulation->SetParameter(6, psi3);
1015             fFitModulation->SetParameter(7, TMath::Sqrt(cos3*cos3+sin3*sin3)/RhoVal());
1016         } break;
1017         default : break;
1018     }
1019     _tempSwap.Fit(fFitModulation, fFitModulationOptions.Data(), "", 0, TMath::TwoPi());
1020     // the quality of the fit is evaluated from 1 - the cdf of the chi square distribution
1021     Double_t CDF(1.-ChiSquareCDF(fFitModulation->GetNDF(), fFitModulation->GetChisquare()));
1022 //    Double_t PDF(ChiSquarePDF(fFitModulation->GetNDF(), fFitModulation->GetChisquare()));
1023     fHistPvalueCDF->Fill(CDF);
1024 //    fHistPvaluePDF->Fill(PDF);
1025     if(CDF > fMinPvalue && CDF < fMaxPvalue && ( fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) > 0)) { // fit quality
1026         // for LOCAL didactic purposes, save the  best and the worst fits
1027         // this routine can produce a lot of output histograms (it's not memory 'safe') and will not work on GRID 
1028         // since the output will become unmergeable (i.e. different nodes may produce conflicting output)
1029         switch (fRunModeType) {
1030             case kLocal : {
1031                 if(fRandom->Uniform(0, 100) > fPercentageOfFits) break;
1032                 static Int_t didacticCounterBest(0);
1033                 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
1034                 TF1* didactifFit = (TF1*)fFitModulation->Clone(Form("fit_%i_CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
1035                 didacticProfile->GetListOfFunctions()->Add(didactifFit);
1036                 fOutputListGood->Add(didacticProfile);
1037                 didacticCounterBest++;
1038                 TH2F* didacticSurface = BookTH2F(Form("surface_%s", didacticProfile->GetName()), "#phi", "#eta", 50, 0, TMath::TwoPi(), 50, -1, 1, -1, kFALSE);
1039                 for(Int_t i(0); i < iTracks; i++) {
1040                     AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1041                     if(PassesCuts(track)) {
1042                         if(fUsePtWeight) didacticSurface->Fill(track->Phi(), track->Eta(), track->Pt());
1043                         else didacticSurface->Fill(track->Phi(), track->Eta());
1044                     }
1045                 }
1046                 if(fExcludeLeadingJetsFromFit) {       // visualize the excluded region
1047                     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);
1048                     f2->SetParameters(excludeInPt[0]/3.,excludeInPhi[0],.1,excludeInEta[0],.1);
1049                     didacticSurface->GetListOfFunctions()->Add(f2);
1050                     TF2 *f3 = new TF2(Form("%s_NLJ", didacticSurface->GetName()),"[0]*TMath::Gaus(x,[1],[2])*TMath::Gaus(y,[3],[4])", 0, TMath::TwoPi(), -1, 1);
1051                     f3->SetParameters(excludeInPt[1]/3.,excludeInPhi[1],.1,excludeInEta[1],.1);
1052                     f3->SetLineColor(kGreen);
1053                     didacticSurface->GetListOfFunctions()->Add(f3);
1054                 }
1055                 fOutputListGood->Add(didacticSurface);
1056             } break;
1057             default : break;
1058         }
1059     } else {    // if the fit is of poor quality revert to the original rho estimate
1060         switch (fRunModeType) { // again see if we want to save the fit
1061             case kLocal : {
1062                 static Int_t didacticCounterWorst(0);
1063                 if(fRandom->Uniform(0, 100) > fPercentageOfFits) break;
1064                 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data() ));
1065                 TF1* didactifFit = (TF1*)fFitModulation->Clone(Form("fit_%i_p_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data()));
1066                 didacticProfile->GetListOfFunctions()->Add(didactifFit);
1067                 fOutputListBad->Add(didacticProfile);
1068                 didacticCounterWorst++;
1069                 } break;
1070             default : break;
1071         }
1072         switch (fFitModulationType) {
1073             case kNoFit : break;        // nothing to do
1074             case kCombined : fFitModulation->SetParameter(7, 0);        // no break
1075             case kFourierSeries : fFitModulation->SetParameter(7, 0);   // no break
1076             default : { // needs to be done if there was a poor fit
1077                  fFitModulation->SetParameter(3, 0);
1078                  fFitModulation->SetParameter(0, RhoVal());
1079             } break;
1080         }
1081         return kFALSE;  // return false if the fit is rejected
1082     }
1083     return kTRUE;
1084 }
1085 //_____________________________________________________________________________
1086 Bool_t AliAnalysisTaskRhoVnModulation::PassesCuts(AliVEvent* event)
1087 {
1088     // event cuts
1089     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1090     if(!event) return kFALSE;
1091     if(TMath::Abs(InputEvent()->GetPrimaryVertex()->GetZ()) > 10.) return kFALSE;
1092     // aod and esd specific checks
1093     switch (fDataType) {
1094        case kESD: {
1095             AliESDEvent* esdEvent = static_cast<AliESDEvent*>(InputEvent());
1096             if( (!esdEvent) || (TMath::Abs(esdEvent->GetPrimaryVertexSPD()->GetZ() - esdEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE; 
1097        } break;
1098        case kAOD: {
1099             AliAODEvent* aodEvent = static_cast<AliAODEvent*>(InputEvent());
1100             if( (!aodEvent) || (TMath::Abs(aodEvent->GetPrimaryVertexSPD()->GetZ() - aodEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE; 
1101        } break;
1102        default: break;
1103     }
1104     fCent = InputEvent()->GetCentrality()->GetCentralityPercentile("V0M");
1105     if(fCent <= fCentralityClasses->At(0) || fCent >= fCentralityClasses->At(fCentralityClasses->GetSize()-1) || TMath::Abs(fCent-InputEvent()->GetCentrality()->GetCentralityPercentile("TRK")) > 5.) return kFALSE;
1106     // determine centrality class
1107     for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
1108         if(fCent >= fCentralityClasses->At(i) && fCent <= fCentralityClasses->At(1+i)) {
1109             fInCentralitySelection = i;
1110             break; }
1111     } 
1112     if(fExplicitOutlierCut == 2010 || fExplicitOutlierCut == 2011) {
1113        if(!PassesCuts(fExplicitOutlierCut)) return kFALSE;
1114     }
1115     if(fFillQAHistograms) FillQAHistograms(event);
1116     return kTRUE;
1117 }
1118 //_____________________________________________________________________________
1119 Bool_t AliAnalysisTaskRhoVnModulation::PassesCuts(Int_t year) 
1120 {
1121     // additional centrality cut based on relation between tpc and global multiplicity
1122     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1123     AliAODEvent* event(dynamic_cast<AliAODEvent*>(InputEvent()));
1124     if(!event) return kFALSE;
1125     Int_t multTPC(0), multGlob(0), nTracks(InputEvent()->GetNumberOfTracks());
1126     for(Int_t iTracks = 0; iTracks < nTracks; iTracks++) { 
1127         AliAODTrack* track = event->GetTrack(iTracks);
1128         if(!track) continue;
1129         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
1130         if (track->TestFilterBit(1) && track->Chi2perNDF() > 0.2) multTPC++;
1131         if (!track->TestFilterBit(16) || track->Chi2perNDF() < 0.1) continue;
1132         Double_t b[2] = {-99., -99.};
1133         Double_t bCov[3] = {-99., -99., -99.};
1134         if (track->PropagateToDCA(event->GetPrimaryVertex(), event->GetMagneticField(), 100., b, bCov) && TMath::Abs(b[0]) < 0.3 && TMath::Abs(b[1]) < 0.3) multGlob++;
1135     }
1136     if(year == 2010 && multTPC > (-40.3+1.22*multGlob) && multTPC < (32.1+1.59*multGlob)) return kTRUE;
1137     if(year == 2011  && multTPC > (-36.73 + 1.48*multGlob) && multTPC < (62.87 + 1.78*multGlob)) return kTRUE;
1138     return kFALSE;
1139 }
1140 //_____________________________________________________________________________
1141 Bool_t AliAnalysisTaskRhoVnModulation::PassesCuts(const AliVCluster* cluster) const
1142 {
1143     // cluster cuts
1144     if(fDebug > 1) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1145     if(!cluster) return kFALSE;
1146     return kTRUE;
1147 }
1148 //_____________________________________________________________________________
1149 void AliAnalysisTaskRhoVnModulation::FillHistogramsAfterSubtraction(Double_t vzero[2][2], Double_t* tpc) const
1150 {
1151     // fill histograms 
1152     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1153     FillTrackHistograms();
1154     /* FillClusterHistograms(); */
1155     FillJetHistograms(vzero, tpc); 
1156     /* FillCorrectedClusterHistograms(); */
1157     FillEventPlaneHistograms(vzero, tpc);
1158     FillRhoHistograms();
1159     FillDeltaPtHistograms(vzero, tpc);
1160     FillDeltaPhiHistograms(vzero, tpc);
1161 }
1162 //_____________________________________________________________________________
1163 void AliAnalysisTaskRhoVnModulation::FillTrackHistograms() const
1164 {
1165     // fill track histograms
1166     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1167     Int_t iTracks(fTracks->GetEntriesFast());
1168     for(Int_t i(0); i < iTracks; i++) {
1169         AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
1170         if(!PassesCuts(track)) continue;
1171         fHistPicoTrackPt[fInCentralitySelection]->Fill(track->Pt());
1172         if(fFillQAHistograms) FillQAHistograms(track);
1173     }
1174 }
1175 //_____________________________________________________________________________
1176 void AliAnalysisTaskRhoVnModulation::FillClusterHistograms() const
1177 {
1178     // fill cluster histograms
1179     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1180     /* Int_t iClusters(fCaloClusters->GetEntriesFast());
1181     for(Int_t i(0); i < iClusters; i++) {
1182         AliVCluster* cluster = static_cast<AliVCluster*>(fCaloClusters->At(iClusters));
1183         if (!PassesCuts(cluster)) continue;
1184         TLorentzVector clusterLorentzVector;
1185         cluster->GetMomentum(clusterLorentzVector, const_cast<Double_t*>(fVertex));
1186         fHistClusterPt[fInCentralitySelection]->Fill(clusterLorentzVector.Pt());
1187         fHistClusterEta[fInCentralitySelection]->Fill(clusterLorentzVector.Eta());
1188         fHistClusterPhi[fInCentralitySelection]->Fill(clusterLorentzVector.Phi());
1189     }
1190     return; */
1191 }
1192 //_____________________________________________________________________________
1193 void AliAnalysisTaskRhoVnModulation::FillCorrectedClusterHistograms() const
1194 {
1195     // fill clusters after hadronic correction FIXME implement
1196     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1197 }
1198 //_____________________________________________________________________________
1199 void AliAnalysisTaskRhoVnModulation::FillEventPlaneHistograms(Double_t vzero[2][2], Double_t* tpc) const
1200 {
1201     // fill event plane histograms
1202     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1203     fHistPsiControl->Fill(0.5, vzero[0][0]);    // vzero a psi2
1204     fHistPsiControl->Fill(1.5, vzero[1][0]);    // vzero c psi2
1205     fHistPsiControl->Fill(2.5, tpc[0]);         // tpc psi 2
1206     fHistPsiControl->Fill(5.5, vzero[0][1]);    // vzero a psi3
1207     fHistPsiControl->Fill(6.5, vzero[1][1]);    // vzero b psi3
1208     fHistPsiControl->Fill(7.5, tpc[1]);         // tpc psi 3
1209     fHistPsiVZEROA->Fill(vzero[0][0]);
1210     fHistPsiVZEROC->Fill(vzero[1][0]);
1211     fHistPsiTPC->Fill(tpc[0]);
1212     fHistPsiSpread->Fill(0.5, TMath::Abs(vzero[0][0]-vzero[1][0]));
1213     fHistPsiSpread->Fill(1.5, TMath::Abs(vzero[0][0]-tpc[0]));
1214     fHistPsiSpread->Fill(2.5, TMath::Abs(vzero[1][0]-tpc[0]));
1215 }
1216 //_____________________________________________________________________________
1217 void AliAnalysisTaskRhoVnModulation::FillRhoHistograms() const
1218 {
1219     // fill rho histograms
1220     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1221     fHistRhoPackage[fInCentralitySelection]->Fill(RhoVal());    // save the rho estimate from the emcal jet package
1222     // get multiplicity FIXME inefficient
1223     Int_t iTracks(fTracks->GetEntriesFast()), mult(0), iJets(fJets->GetEntriesFast());
1224     for(Int_t i(0); i < iTracks; i ++) { if(PassesCuts(static_cast<AliVTrack*>(fTracks->At(i)))) mult++; }
1225     Double_t rho(RhoVal(TMath::Pi(), TMath::Pi(), fRho->GetVal()));
1226     fHistRho[fInCentralitySelection]->Fill(rho);
1227     fHistRhoVsMult->Fill(mult, rho);
1228     fHistRhoVsCent->Fill(fCent, rho);
1229     for(Int_t i(0); i < iJets; i++) {
1230         AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
1231         if(!PassesCuts(jet)) continue;
1232         fHistRhoAVsMult->Fill(mult, rho * jet->Area());
1233         fHistRhoAVsCent->Fill(fCent, rho * jet->Area());
1234     }
1235 }
1236 //_____________________________________________________________________________
1237 void AliAnalysisTaskRhoVnModulation::FillDeltaPtHistograms(Double_t vzero[2][2], Double_t* tpc) const
1238 {
1239     // fill delta pt histograms
1240     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1241     Int_t i(0), maxCones(20);
1242     AliEmcalJet* leadingJet(0x0);
1243     static Int_t sJets[9999] = {-1};
1244     GetSortedArray(sJets, fJets);
1245     do { // get the leading jet 
1246         leadingJet = static_cast<AliEmcalJet*>(fJets->At(sJets[i]));
1247         i++;
1248     }
1249     while (!PassesCuts(leadingJet)&&i<fJets->GetEntriesFast()); 
1250     if(!leadingJet && fDebug > 0) printf(" > failed to retrieve leading jet ! < \n");
1251     const Float_t areaRC = fRandomConeRadius*fRandomConeRadius*TMath::Pi();
1252     // we're retrieved the leading jet, now get a random cone
1253     for(i = 0; i < maxCones; i++) {
1254        Float_t pt(0), eta(0), phi(0);
1255        // get a random cone without constraints on leading jet position
1256        CalculateRandomCone(pt, eta, phi, 0x0);
1257        if(pt > 0) {
1258            if(fFillQAHistograms) fHistRCPhiEta[fInCentralitySelection]->Fill(phi, eta);
1259            fHistRhoVsRCPt[fInCentralitySelection]->Fill(pt, RhoVal(phi, fJetRadius, fRho->GetVal())*areaRC);
1260            fHistRCPt[fInCentralitySelection]->Fill(pt);
1261            fHistDeltaPtDeltaPhi2TPC[fInCentralitySelection]->Fill(PhaseShift(phi-tpc[0], 2.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1262            fHistDeltaPtDeltaPhi2V0A[fInCentralitySelection]->Fill(PhaseShift(phi-vzero[0][0], 2.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1263            fHistDeltaPtDeltaPhi2V0C[fInCentralitySelection]->Fill(PhaseShift(phi-vzero[1][0], 2.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1264            fHistDeltaPtDeltaPhi3TPC[fInCentralitySelection]->Fill(PhaseShift(phi-tpc[1], 3.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1265            fHistDeltaPtDeltaPhi3V0A[fInCentralitySelection]->Fill(PhaseShift(phi-vzero[0][1], 3.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1266            fHistDeltaPtDeltaPhi3V0C[fInCentralitySelection]->Fill(PhaseShift(phi-vzero[1][1], 3.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1267        }
1268        // get a random cone excluding leading jet area
1269        CalculateRandomCone(pt, eta, phi, leadingJet);
1270        if(pt > 0) {
1271            if(fFillQAHistograms) fHistRCPhiEtaExLJ[fInCentralitySelection]->Fill(phi, eta);
1272            fHistRhoVsRCPtExLJ[fInCentralitySelection]->Fill(pt, RhoVal(phi, fJetRadius, fRho->GetVal())*areaRC);
1273            fHistRCPtExLJ[fInCentralitySelection]->Fill(pt);
1274            fHistDeltaPtDeltaPhi2ExLJTPC[fInCentralitySelection]->Fill(PhaseShift(phi-tpc[0], 2.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1275            fHistDeltaPtDeltaPhi2ExLJV0A[fInCentralitySelection]->Fill(PhaseShift(phi-vzero[0][0], 2.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1276            fHistDeltaPtDeltaPhi2ExLJV0C[fInCentralitySelection]->Fill(PhaseShift(phi-vzero[1][0], 2.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1277            fHistDeltaPtDeltaPhi3ExLJTPC[fInCentralitySelection]->Fill(PhaseShift(phi-tpc[1], 3.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1278            fHistDeltaPtDeltaPhi3ExLJV0A[fInCentralitySelection]->Fill(PhaseShift(phi-vzero[0][1], 3.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1279            fHistDeltaPtDeltaPhi3ExLJV0C[fInCentralitySelection]->Fill(PhaseShift(phi-vzero[1][1], 3.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1280        }
1281        // get a random cone in an event with randomized phi and eta
1282        /* CalculateRandomCone(pt, eta, phi, 0x0, kTRUE);
1283        if( pt > 0) {
1284            fHistRCPhiEtaRand[fInCentralitySelection]->Fill(phi, eta);
1285            fHistRhoVsRCPtRand[fInCentralitySelection]->Fill(pt, RhoVal(phi, fJetRadius, fRho->GetVal())*areaRC);
1286            fHistRCPtRand[fInCentralitySelection]->Fill(pt);
1287            fHistDeltaPtDeltaPhi2Rand[fInCentralitySelection]->Fill(PhaseShift(phi-psi2, 2.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1288            fHistDeltaPtDeltaPhi3Rand[fInCentralitySelection]->Fill(PhaseShift(phi-psi3, 3.), pt - areaRC*RhoVal(phi, fJetRadius, fRho->GetVal()));
1289        } */
1290     } 
1291 }
1292 //_____________________________________________________________________________
1293 void AliAnalysisTaskRhoVnModulation::FillJetHistograms(Double_t vzero[2][2], Double_t* tpc) const
1294 {
1295     // fill jet histograms
1296     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1297     Int_t iJets(fJets->GetEntriesFast());
1298     for(Int_t i(0); i < iJets; i++) {
1299         AliEmcalJet* jet = static_cast<AliEmcalJet*>(fJets->At(i));
1300         if(PassesCuts(jet)) {
1301             Double_t pt(jet->Pt()), area(jet->Area()), eta(jet->Eta()), phi(jet->Phi());
1302             Double_t rho(RhoVal(phi, fJetRadius, fRho->GetVal()));
1303             fHistJetPtRaw[fInCentralitySelection]->Fill(pt);
1304             fHistJetPt[fInCentralitySelection]->Fill(pt-area*rho);
1305             if(fFillQAHistograms) fHistJetEtaPhi[fInCentralitySelection]->Fill(eta, phi);
1306             fHistJetPtArea[fInCentralitySelection]->Fill(pt-area*rho, area);
1307             fHistJetPsiTPCPt[fInCentralitySelection]->Fill(PhaseShift(phi-tpc[0], 2.), pt-area*rho);
1308             fHistJetPsiVZEROAPt[fInCentralitySelection]->Fill(PhaseShift(phi-vzero[0][0], 2.), pt-area*rho);
1309             fHistJetPsiVZEROCPt[fInCentralitySelection]->Fill(PhaseShift(phi-vzero[1][0], 2.), pt-area*rho);
1310             fHistJetPtConstituents[fInCentralitySelection]->Fill(pt-area*rho, jet->Nch());
1311             fHistJetEtaRho[fInCentralitySelection]->Fill(eta, pt/area);
1312             if(fSetPtSub) jet->SetPtSub(pt-area*rho);
1313         }
1314         else { // if the jet is rejected, excluded it for the flow analysis
1315             if(fSetPtSub) jet->SetPtSub(-999.);
1316         }
1317     }
1318 }
1319 //_____________________________________________________________________________
1320 void AliAnalysisTaskRhoVnModulation::FillDeltaPhiHistograms(Double_t vzero[2][2], Double_t* tpc) const
1321 {
1322    // fill phi minus psi histograms
1323    if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1324    if(fTracks) {
1325        Int_t iTracks(fTracks->GetEntriesFast());
1326        for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
1327            AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
1328            if(!PassesCuts(track)) continue;
1329            fHistDeltaPhi2VZEROA[fInCentralitySelection]->Fill(PhaseShift(track->Phi()-vzero[0][0], 2.));
1330            fHistDeltaPhi2VZEROC[fInCentralitySelection]->Fill(PhaseShift(track->Phi()-vzero[1][0], 2.));
1331            fHistDeltaPhi2TPC[fInCentralitySelection]->Fill(PhaseShift(track->Phi()-tpc[0], 2.));
1332            fHistDeltaPhi3VZEROA[fInCentralitySelection]->Fill(PhaseShift(track->Phi()-vzero[0][1], 3.));
1333            fHistDeltaPhi3VZEROC[fInCentralitySelection]->Fill(PhaseShift(track->Phi()-vzero[1][1], 3.));
1334            fHistDeltaPhi3TPC[fInCentralitySelection]->Fill(PhaseShift(track->Phi()-tpc[1], 3.));
1335        }
1336    }
1337 }
1338 //_____________________________________________________________________________
1339 void AliAnalysisTaskRhoVnModulation::FillQAHistograms(AliVTrack* vtrack) const
1340 {
1341     // fill qa histograms for pico tracks
1342     if(!vtrack) return;
1343     AliPicoTrack* track = static_cast<AliPicoTrack*>(vtrack);
1344     fHistRunnumbersPhi->Fill(fMappedRunNumber, track->Phi());
1345     fHistRunnumbersEta->Fill(fMappedRunNumber, track->Eta());
1346     Int_t type((int)(track->GetTrackType()));
1347     switch (type) {
1348         case 0:
1349            fHistPicoCat1[fInCentralitySelection]->Fill(track->Eta(), track->Phi()); 
1350            break;
1351         case 1:
1352            fHistPicoCat2[fInCentralitySelection]->Fill(track->Eta(), track->Phi()); 
1353            break;
1354         case 2:
1355            fHistPicoCat3[fInCentralitySelection]->Fill(track->Eta(), track->Phi()); 
1356            break;
1357         default: break;
1358     }
1359 }
1360 //_____________________________________________________________________________
1361 void AliAnalysisTaskRhoVnModulation::FillQAHistograms(AliVEvent* vevent) 
1362 {
1363     // fill qa histograms for events
1364     if(!vevent) return;
1365     fHistVertexz->Fill(vevent->GetPrimaryVertex()->GetZ());
1366     fHistCentrality->Fill(fCent);
1367     Int_t runNumber(InputEvent()->GetRunNumber());
1368     Int_t runs[] = {167813, 167988, 168066, 168068, 168069, 168076, 168104, 168212, 168311, 168322, 168325, 168341, 168361, 168362, 168458, 168460, 168461, 168992, 169091, 169094, 169138, 169143, 169167, 169417, 169835, 169837, 169838, 169846, 169855, 169858, 169859, 169923, 169956, 170027, 170036, 170081, 169975, 169981, 170038, 170040, 170083, 170084, 170085, 170088, 170089, 170091, 170152, 170155, 170159, 170163, 170193, 170195, 170203, 170204, 170205, 170228, 170230, 170264, 170268, 170269, 170270, 170306, 170308, 170309};
1369     for(fMappedRunNumber = 0; fMappedRunNumber < 64; fMappedRunNumber++) {
1370         if(runs[fMappedRunNumber]==runNumber) break;
1371     }
1372 }
1373 //_____________________________________________________________________________
1374 void AliAnalysisTaskRhoVnModulation::FillAnalysisSummaryHistogram() const
1375 {
1376     // fill the analysis summary histrogram, saves all relevant analysis settigns
1377     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1378     fHistAnalysisSummary->GetXaxis()->SetBinLabel(1, "fJetRadius"); 
1379     fHistAnalysisSummary->SetBinContent(1, fJetRadius);
1380     fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fPtBiasJetTrack");
1381     fHistAnalysisSummary->SetBinContent(2, fPtBiasJetTrack);
1382     fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fPtBiasJetClus");
1383     fHistAnalysisSummary->SetBinContent(3, fPtBiasJetClus);
1384     fHistAnalysisSummary->GetXaxis()->SetBinLabel(4, "fJetPtCut");
1385     fHistAnalysisSummary->SetBinContent(4, fJetPtCut);
1386     fHistAnalysisSummary->GetXaxis()->SetBinLabel(5, "fJetAreaCut");
1387     fHistAnalysisSummary->SetBinContent(5, fJetAreaCut);
1388     fHistAnalysisSummary->GetXaxis()->SetBinLabel(6, "fPercAreaCut");
1389     fHistAnalysisSummary->SetBinContent(6, fPercAreaCut);
1390     fHistAnalysisSummary->GetXaxis()->SetBinLabel(7, "fAreaEmcCut");
1391     fHistAnalysisSummary->SetBinContent(7, fAreaEmcCut);
1392     fHistAnalysisSummary->GetXaxis()->SetBinLabel(8, "fJetMinEta");
1393     fHistAnalysisSummary->SetBinContent(8, fJetMinEta);
1394     fHistAnalysisSummary->GetXaxis()->SetBinLabel(9, "fJetMaxEta");
1395     fHistAnalysisSummary->SetBinContent(9, fJetMaxEta);
1396     fHistAnalysisSummary->GetXaxis()->SetBinLabel(10, "fJetMinPhi");
1397     fHistAnalysisSummary->SetBinContent(10, fJetMinPhi);
1398     fHistAnalysisSummary->GetXaxis()->SetBinLabel(11, "fJetMaxPhi");
1399     fHistAnalysisSummary->SetBinContent(11, fJetMaxPhi);
1400     fHistAnalysisSummary->GetXaxis()->SetBinLabel(12, "fMaxClusterPt");
1401     fHistAnalysisSummary->SetBinContent(12, fMaxClusterPt);
1402     fHistAnalysisSummary->GetXaxis()->SetBinLabel(13, "fMaxTrackPt");
1403     fHistAnalysisSummary->SetBinContent(13, fMaxTrackPt);
1404     fHistAnalysisSummary->GetXaxis()->SetBinLabel(14, "fLeadingHadronType");
1405     fHistAnalysisSummary->SetBinContent(14, fLeadingHadronType);
1406     fHistAnalysisSummary->GetXaxis()->SetBinLabel(15, "fAnaType");
1407     fHistAnalysisSummary->SetBinContent(15, fAnaType);
1408     fHistAnalysisSummary->GetXaxis()->SetBinLabel(16, "fForceBeamType");
1409     fHistAnalysisSummary->SetBinContent(16, fForceBeamType);
1410     fHistAnalysisSummary->GetXaxis()->SetBinLabel(17, "fMinCent");
1411     fHistAnalysisSummary->SetBinContent(17, fMinCent);
1412     fHistAnalysisSummary->GetXaxis()->SetBinLabel(18, "fMaxCent");
1413     fHistAnalysisSummary->SetBinContent(18, fMaxCent);
1414     fHistAnalysisSummary->GetXaxis()->SetBinLabel(19, "fMinVz");
1415     fHistAnalysisSummary->SetBinContent(19, fMinVz);
1416     fHistAnalysisSummary->GetXaxis()->SetBinLabel(20, "fMaxVz");
1417     fHistAnalysisSummary->SetBinContent(20, fMaxVz);
1418     fHistAnalysisSummary->GetXaxis()->SetBinLabel(21, "fOffTrigger");
1419     fHistAnalysisSummary->SetBinContent(21, fOffTrigger);
1420     fHistAnalysisSummary->GetXaxis()->SetBinLabel(22, "fClusPtCut");
1421     fHistAnalysisSummary->SetBinContent(22, fClusPtCut);
1422     fHistAnalysisSummary->GetXaxis()->SetBinLabel(23, "fTrackPtCut");
1423     fHistAnalysisSummary->SetBinContent(23, fTrackPtCut);
1424     fHistAnalysisSummary->GetXaxis()->SetBinLabel(24, "fTrackMinEta");
1425     fHistAnalysisSummary->SetBinContent(24, fTrackMinEta);
1426     fHistAnalysisSummary->GetXaxis()->SetBinLabel(25, "fTrackMaxEta");
1427     fHistAnalysisSummary->SetBinContent(25, fTrackMaxEta);
1428     fHistAnalysisSummary->GetXaxis()->SetBinLabel(26, "fTrackMinPhi");
1429     fHistAnalysisSummary->SetBinContent(26, fTrackMinPhi);
1430     fHistAnalysisSummary->GetXaxis()->SetBinLabel(27, "fTrackMaxPhi");
1431     fHistAnalysisSummary->SetBinContent(27, fTrackMaxPhi);
1432     fHistAnalysisSummary->GetXaxis()->SetBinLabel(28, "fClusTimeCutLow");
1433     fHistAnalysisSummary->SetBinContent(28, fClusTimeCutLow);
1434     fHistAnalysisSummary->GetXaxis()->SetBinLabel(29, "fClusTimeCutUp");
1435     fHistAnalysisSummary->SetBinContent(29, fClusTimeCutUp);
1436     fHistAnalysisSummary->GetXaxis()->SetBinLabel(30, "fMinPtTrackInEmcal");
1437     fHistAnalysisSummary->SetBinContent(30, fMinPtTrackInEmcal);
1438     fHistAnalysisSummary->GetXaxis()->SetBinLabel(31, "fEventPlaneVsEmcal");
1439     fHistAnalysisSummary->SetBinContent(31, fEventPlaneVsEmcal);
1440     fHistAnalysisSummary->GetXaxis()->SetBinLabel(32, "fMinEventPlane");
1441     fHistAnalysisSummary->SetBinContent(32, fMaxEventPlane);
1442     fHistAnalysisSummary->GetXaxis()->SetBinLabel(33, "fRandomConeRadius");
1443     fHistAnalysisSummary->SetBinContent(33, fRandomConeRadius);
1444     fHistAnalysisSummary->GetXaxis()->SetBinLabel(34, "fitModulationType");
1445     fHistAnalysisSummary->SetBinContent(34, (int)fFitModulationType);
1446     fHistAnalysisSummary->GetXaxis()->SetBinLabel(35, "runModeType");
1447     fHistAnalysisSummary->SetBinContent(35, (int)fRunModeType);
1448     fHistAnalysisSummary->GetXaxis()->SetBinLabel(36, "data type");
1449     fHistAnalysisSummary->SetBinContent(36, (int)fDataType);
1450     fHistAnalysisSummary->GetXaxis()->SetBinLabel(37, "iterator");
1451     fHistAnalysisSummary->SetBinContent(37, 1.);
1452     fHistAnalysisSummary->GetXaxis()->SetBinLabel(38, "fMinPvalue");
1453     fHistAnalysisSummary->SetBinContent(38, fMinPvalue);
1454     fHistAnalysisSummary->GetXaxis()->SetBinLabel(39, "fMaxPvalue");
1455     fHistAnalysisSummary->SetBinContent(39, fMaxPvalue);
1456     fHistAnalysisSummary->GetXaxis()->SetBinLabel(40, "fExcludeLeadingJetsFromFit");
1457     fHistAnalysisSummary->SetBinContent(40, fExcludeLeadingJetsFromFit);
1458     fHistAnalysisSummary->GetXaxis()->SetBinLabel(41, "fRebinSwapHistoOnTheFly");
1459     fHistAnalysisSummary->SetBinContent(41, (int)fRebinSwapHistoOnTheFly);
1460     fHistAnalysisSummary->GetXaxis()->SetBinLabel(42, "fUsePtWeight");
1461     fHistAnalysisSummary->SetBinContent(42, (int)fUsePtWeight);
1462     fHistAnalysisSummary->GetXaxis()->SetBinLabel(43, "fMinLeadingHadronPt");
1463     fHistAnalysisSummary->SetBinContent(43, fMinLeadingHadronPt);
1464     fHistAnalysisSummary->GetXaxis()->SetBinLabel(44, "fExplicitOutlierCut");
1465     fHistAnalysisSummary->SetBinContent(44, fExplicitOutlierCut);
1466     fHistAnalysisSummary->GetXaxis()->SetBinLabel(45, "fLocalJetMinEta");
1467     fHistAnalysisSummary->SetBinContent(45,fLocalJetMinEta );
1468     fHistAnalysisSummary->GetXaxis()->SetBinLabel(46, "fLocalJetMaxEta");
1469     fHistAnalysisSummary->SetBinContent(46, fLocalJetMaxEta);
1470     fHistAnalysisSummary->GetXaxis()->SetBinLabel(47, "fLocalJetMinPhi");
1471     fHistAnalysisSummary->SetBinContent(47, fLocalJetMinPhi);
1472     fHistAnalysisSummary->GetXaxis()->SetBinLabel(48, "fLocalJetMaxPhi");
1473     fHistAnalysisSummary->SetBinContent(48, fLocalJetMaxPhi);
1474     fHistAnalysisSummary->GetXaxis()->SetBinLabel(49, "fSoftTrackMinPt");
1475     fHistAnalysisSummary->SetBinContent(49, fSoftTrackMinPt);
1476     fHistAnalysisSummary->GetXaxis()->SetBinLabel(50, "fSoftTrackMaxPt");
1477     fHistAnalysisSummary->SetBinContent(50, fSoftTrackMaxPt);
1478 }
1479 //_____________________________________________________________________________
1480 void AliAnalysisTaskRhoVnModulation::Terminate(Option_t *)
1481 {
1482     // terminate
1483     switch (fRunModeType) {
1484         case kLocal : {
1485         printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
1486         if(fFillQAHistograms) {
1487             Int_t runs[] = {167813, 167988, 168066, 168068, 168069, 168076, 168104, 168212, 168311, 168322, 168325, 168341, 168361, 168362, 168458, 168460, 168461, 168992, 169091, 169094, 169138, 169143, 169167, 169417, 169835, 169837, 169838, 169846, 169855, 169858, 169859, 169923, 169956, 170027, 170036, 170081, 169975, 169981, 170038, 170040, 170083, 170084, 170085, 170088, 170089, 170091, 170152, 170155, 170159, 170163, 170193, 170195, 170203, 170204, 170205, 170228, 170230, 170264, 170268, 170269, 170270, 170306, 170308, 170309};
1488             for(Int_t i(0); i < 64; i++) { 
1489                 fHistRunnumbersPhi->GetXaxis()->SetBinLabel(i+1, Form("%i", runs[i]));
1490                 fHistRunnumbersEta->GetXaxis()->SetBinLabel(i+1, Form("%i", runs[i]));
1491             }
1492             fHistRunnumbersPhi->GetXaxis()->SetBinLabel(65, "undetermined");
1493             fHistRunnumbersEta->GetXaxis()->SetBinLabel(65, "undetermined");
1494         }
1495         AliAnalysisTaskRhoVnModulation::Dump();
1496         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));
1497         } break;
1498         default : break;
1499     }
1500 }
1501 //_____________________________________________________________________________
1502 TH1F* AliAnalysisTaskRhoVnModulation::GetResolutionFromOuptutFile(detectorType det, Int_t h, TArrayD* cen)
1503 {
1504     // INTERFACE METHOD FOR OUTPUTFILE
1505     // get the detector resolution, user has ownership of the returned histogram
1506     if(!fOutputList) {
1507         printf(" > Please add fOutputList first < \n");
1508         return 0x0;
1509     }
1510     TH1F* r(0x0);
1511     (cen) ? r = new TH1F("R", "R", cen->GetSize()-1, cen->GetArray()) : r = new TH1F("R", "R", 10, 0, 10);
1512     if(!cen) r->GetXaxis()->SetTitle("number of centrality bin");
1513     r->GetYaxis()->SetTitle(Form("Resolution #Psi_{%i}", h));
1514     for(Int_t i(0); i < 10; i++) {
1515         TProfile* temp((TProfile*)fOutputList->FindObject(Form("fProfV%iResolution_%i", h, i)));
1516         if(!temp) break;
1517         Double_t a(temp->GetBinContent(3)), b(temp->GetBinContent(5)), c(temp->GetBinContent(7));
1518         Double_t _a(temp->GetBinError(3)), _b(temp->GetBinError(5)), _c(temp->GetBinError(7));
1519         if(a <= 0 || b <= 0 || c <= 0) continue;
1520         switch (det) {
1521             case kVZEROA : {
1522                 r->SetBinContent(1+i, TMath::Sqrt((a*b)/c));
1523                 if(i==0) r->SetNameTitle("VZEROA resolution", "VZEROA resolution");
1524             } break;
1525             case kVZEROC : {
1526                 r->SetBinContent(1+i, TMath::Sqrt((a*c)/b));
1527                 if(i==0) r->SetNameTitle("VZEROC resolution", "VZEROC resolution");
1528             } break;
1529             case kTPC : {
1530                 r->SetBinContent(1+i, TMath::Sqrt((b*c)/a));
1531                 if(i==0) r->SetNameTitle("TPC resolution", "TPC resolution");
1532             } break;
1533             default : break;
1534         }
1535         r->SetBinError(1+i, TMath::Sqrt(_a*_a+_b*_b+_c*_c));
1536     }
1537     return r;
1538 }
1539 //_____________________________________________________________________________
1540 TH1F* AliAnalysisTaskRhoVnModulation::CorrectForResolutionDiff(TH1F* v, detectorType det, TArrayD* cen, Int_t c, Int_t h)
1541 {
1542     // INTERFACE METHOD FOR OUTPUT FILE
1543     // correct the supplied differential vn histogram v for detector resolution
1544     TH1F* r(GetResolutionFromOuptutFile(det, h, cen));
1545     if(!r) {
1546         printf(" > Couldn't find resolution < \n");
1547         return 0x0;
1548     }
1549     Double_t res(1./r->GetBinContent(1+r->FindBin(c)));
1550     TF1* line = new TF1("line", "pol0", 0, 200);
1551     line->SetParameter(0, res);
1552     return (v->Multiply(line)) ? v : 0x0;
1553 }
1554 //_____________________________________________________________________________
1555 TH1F* AliAnalysisTaskRhoVnModulation::CorrectForResolutionInt(TH1F* v, detectorType det, TArrayD* cen, Int_t h)
1556 {
1557     // INTERFACE METHOD FOR OUTPUT FILE
1558     // correct the supplied intetrated vn histogram v for detector resolution
1559     // integrated vn must have the same centrality binning as the resolotion correction
1560     TH1F* r(GetResolutionFromOuptutFile(det, h, cen));
1561     return (v->Divide(v, r)) ? v : 0x0;
1562 }
1563 //_____________________________________________________________________________