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