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