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