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