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