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