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