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