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