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