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