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