]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/AliAnalysisTaskLocalRho.cxx
Minor fix for warning about destructor call for TF1 in SetModulationFit
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / AliAnalysisTaskLocalRho.cxx
1 // $Id$
2 // 
3 // analysis task to estimate an event's local energy density
4 //
5 // This task is part of the emcal jet framework and should be run in the emcaljet train
6 // The following extensions to an accepted AliVEvent are expected:
7 //      - (anti-kt) jets -> necessary if one wants to exclude leading jet contribution to the event plane
8 //      - background estimate of rho -> this task estimates modulation, not rho itself
9 //      - pico tracks -> a uniform track selection is necessary to estimate the contribution of v_n harmonics
10 //      aod's and esd's are handled transparently
11 // The task will estimates a phi-dependent background density rho 
12 // which is added to the event as a AliLocalRhoParamter object
13 //
14 // Author: Redmer Alexander Bertens, Utrecht Univeristy, Utrecht, Netherlands
15 // rbertens@cern.ch, rbertens@nikhef.nl, r.a.bertens@uu.nl 
16
17 // root includes
18 #include <TStyle.h>
19 #include <TRandom3.h>
20 #include <TChain.h>
21 #include <TMath.h>
22 #include <TF1.h>
23 #include <TF2.h>
24 #include <TH1F.h>
25 #include <TH2F.h>
26 #include <TProfile.h>
27 // aliroot includes
28 #include <AliAnalysisTask.h>
29 #include <AliAnalysisManager.h>
30 #include <AliCentrality.h>
31 #include <AliVVertex.h>
32 #include <AliESDEvent.h>
33 #include <AliAODEvent.h>
34 #include <AliAODTrack.h>
35 // emcal jet framework includes
36 #include <AliPicoTrack.h>
37 #include <AliEmcalJet.h>
38 #include <AliRhoParameter.h>
39 #include <AliLocalRhoParameter.h>
40 #include <AliAnalysisTaskLocalRho.h>
41
42 class AliAnalysisTaskLocalRho;
43 using namespace std;
44
45 ClassImp(AliAnalysisTaskLocalRho)
46
47 AliAnalysisTaskLocalRho::AliAnalysisTaskLocalRho() : AliAnalysisTaskEmcalJet("AliAnalysisTaskLocalRho", kTRUE), 
48     fDebug(0), fInitialized(0), fAttachToEvent(kTRUE), fFillHistograms(kFALSE), fNoEventWeightsForQC(kTRUE), fLocalRhoName(GetName()), fCentralityClasses(0), fUserSuppliedV2(0), fUserSuppliedV3(0), fUserSuppliedR2(0), fUserSuppliedR3(0), fNAcceptedTracks(0), fNAcceptedTracksQCn(0), fInCentralitySelection(-1), fFitModulationType(kNoFit), fQCRecovery(kTryFit), fUsePtWeight(kTRUE), fDetectorType(kTPC), fFitModulationOptions("WLQI"), fRunModeType(kGrid), fFitModulation(0), fMinPvalue(0.01), fMaxPvalue(1), fLocalRho(0), fLocalJetMinEta(-10), fLocalJetMaxEta(-10), fLocalJetMinPhi(-10), fLocalJetMaxPhi(-10), fSoftTrackMinPt(0.15), fSoftTrackMaxPt(5.), fHistPvalueCDF(0), fAbsVnHarmonics(kTRUE), fExcludeLeadingJetsFromFit(1.), fRebinSwapHistoOnTheFly(kTRUE), fPercentageOfFits(10.), fUseV0EventPlaneFromHeader(kTRUE), fOutputList(0), fOutputListGood(0), fOutputListBad(0), fHistSwap(0), fHistAnalysisSummary(0), fProfV2(0), fProfV2Cumulant(0), fProfV3(0), fProfV3Cumulant(0) {
49     for(Int_t i(0); i < 10; i++) {
50         fHistPsi2[i] = 0; 
51         fHistPsi3[i] = 0;
52     }
53     // default constructor
54 }
55 //_____________________________________________________________________________
56 AliAnalysisTaskLocalRho::AliAnalysisTaskLocalRho(const char* name, runModeType type) : AliAnalysisTaskEmcalJet(name, kTRUE),
57     fDebug(0), fInitialized(0), fAttachToEvent(kTRUE), fFillHistograms(kFALSE), fNoEventWeightsForQC(kTRUE), fLocalRhoName(GetName()), fCentralityClasses(0), fUserSuppliedV2(0), fUserSuppliedV3(0), fUserSuppliedR2(0), fUserSuppliedR3(0), fNAcceptedTracks(0), fNAcceptedTracksQCn(0), fInCentralitySelection(-1), fFitModulationType(kNoFit), fQCRecovery(kTryFit), fUsePtWeight(kTRUE), fDetectorType(kTPC), fFitModulationOptions("WLQI"), fRunModeType(type), fFitModulation(0), fMinPvalue(0.01), fMaxPvalue(1), fLocalRho(0), fLocalJetMinEta(-10), fLocalJetMaxEta(-10), fLocalJetMinPhi(-10), fLocalJetMaxPhi(-10), fSoftTrackMinPt(0.15), fSoftTrackMaxPt(5.), fHistPvalueCDF(0), fAbsVnHarmonics(kTRUE), fExcludeLeadingJetsFromFit(1.), fRebinSwapHistoOnTheFly(kTRUE), fPercentageOfFits(10.), fUseV0EventPlaneFromHeader(kTRUE), fOutputList(0), fOutputListGood(0), fOutputListBad(0), fHistSwap(0), fHistAnalysisSummary(0), fProfV2(0), fProfV2Cumulant(0), fProfV3(0), fProfV3Cumulant(0) {
58     for(Int_t i(0); i < 10; i++) {
59         fHistPsi2[i] = 0; 
60         fHistPsi3[i] = 0;
61     }
62     // constructor
63     DefineInput(0, TChain::Class());
64     DefineOutput(1, TList::Class());
65     switch (fRunModeType) {
66         case kLocal : {
67             gStyle->SetOptFit(1);
68             DefineOutput(2, TList::Class());
69             DefineOutput(3, TList::Class());
70         } break;
71         default: fDebug = -1;   // suppress debug info explicitely when not running locally
72     }
73 }
74 //_____________________________________________________________________________
75 AliAnalysisTaskLocalRho::~AliAnalysisTaskLocalRho()
76 {
77     // destructor
78     if(fOutputList)             delete fOutputList;
79     if(fOutputListGood)         delete fOutputListGood;
80     if(fOutputListBad)          delete fOutputListBad;
81     if(fFitModulation)          delete fFitModulation;
82     if(fHistSwap)               delete fHistSwap;
83 }
84 //_____________________________________________________________________________
85 Bool_t AliAnalysisTaskLocalRho::InitializeAnalysis() 
86 {
87     // initialize the anaysis
88     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
89     if(fLocalJetMinEta > -10 && fLocalJetMaxEta > -10) SetJetEtaLimits(fLocalJetMinEta, fLocalJetMaxEta);
90     if(fLocalJetMinPhi > -10 && fLocalJetMaxPhi > -10) SetJetPhiLimits(fLocalJetMinPhi, fLocalJetMaxPhi);
91     switch (fFitModulationType)  {
92         case kNoFit : { SetModulationFit(new TF1("fit_kNoFit", "[0]", 0, TMath::TwoPi())); } break;
93         case kV2 : {
94             SetModulationFit(new TF1("fit_kV2", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
95             fFitModulation->SetParameter(0, 0.);        // normalization
96             fFitModulation->SetParameter(3, 0.2);       // v2
97             fFitModulation->FixParameter(1, 1.);        // constant
98             fFitModulation->FixParameter(2, 2.);        // constant
99         } break;
100         case kV3: {
101             SetModulationFit(new TF1("fit_kV3", "[0]*([1]+[2]*[3]*TMath::Cos([2]*(x-[4])))", 0, TMath::TwoPi()));
102             fFitModulation->SetParameter(0, 0.);        // normalization
103             fFitModulation->SetParameter(3, 0.2);       // v3
104             fFitModulation->FixParameter(1, 1.);        // constant
105             fFitModulation->FixParameter(2, 3.);        // constant
106         } break;
107         default : { // for the combined fit, the 'direct fourier series' or the user supplied vn values we use v2 and v3
108              SetModulationFit(new TF1("fit_kCombined", "[0]*([1]+[2]*([3]*TMath::Cos([2]*(x-[4]))+[7]*TMath::Cos([5]*(x-[6]))))", 0, TMath::TwoPi()));
109              fFitModulation->SetParameter(0, 0.);       // normalization
110              fFitModulation->SetParameter(3, 0.2);      // v2
111              fFitModulation->FixParameter(1, 1.);       // constant
112              fFitModulation->FixParameter(2, 2.);       // constant
113              fFitModulation->FixParameter(5, 3.);       // constant
114              fFitModulation->SetParameter(7, 0.2);      // v3
115         } break;
116     }
117     switch (fRunModeType) {
118         case kGrid : { fFitModulationOptions += "N0"; } break;
119         default : break;
120     }
121     fLocalRho = new AliLocalRhoParameter(fLocalRhoName.Data(), 0); 
122     // add the local rho to the event if necessary
123     if(fAttachToEvent) {
124         if(!(InputEvent()->FindListObject(fLocalRho->GetName()))) {
125             InputEvent()->AddObject(fLocalRho);
126         } else {
127             AliFatal(Form("%s: Container with same name %s already present. Aborting", GetName(), fLocalRho->GetName()));
128         }
129     }
130     FillAnalysisSummaryHistogram();
131     return kTRUE;
132 }
133 //_____________________________________________________________________________
134 void AliAnalysisTaskLocalRho::UserCreateOutputObjects()
135 {
136     // create output objects
137     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
138     fHistSwap = new TH1F("fHistSwap", "fHistSwap", 20, 0, TMath::TwoPi());
139     if(!fCentralityClasses) {   // classes must be defined at this point
140         Int_t c[] = {0, 10, 20, 30, 40, 50, 60, 70, 80, 90, 100};
141         fCentralityClasses = new TArrayI(sizeof(c)/sizeof(c[0]), c);
142     }
143     fOutputList = new TList();
144     fOutputList->SetOwner(kTRUE);
145     // the analysis summary histo which stores all the analysis flags is always written to file
146     fHistAnalysisSummary = BookTH1F("fHistAnalysisSummary", "flag", 50, -0.5, 50.5);
147     if(!fFillHistograms) {
148         PostData(1, fOutputList);
149         return;
150     }
151     for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {     
152         fHistPsi2[i] = BookTH1F("fHistPsi2", "#Psi_{2}", 100, -.5*TMath::Pi(), .5*TMath::Pi(), i);
153         fHistPsi3[i] = BookTH1F("fHistPsi3", "#Psi_{3}", 100, -1.*TMath::Pi()/3., TMath::Pi()/3., i);
154     }
155     // cdf of chisquare distribution
156     fHistPvalueCDF = BookTH1F("fHistPvalueCDF", "CDF #chi^{2}", 500, 0, 1);
157     fOutputList->Add(fHistPvalueCDF);
158     // vn profiles
159     Float_t temp[fCentralityClasses->GetSize()];
160     for(Int_t i(0); i < fCentralityClasses->GetSize(); i++) temp[i] = fCentralityClasses->At(i);
161     fProfV2 = new TProfile("fProfV2", "fProfV2", fCentralityClasses->GetSize()-1, temp);
162     fProfV3 = new TProfile("fProfV3", "fProfV3", fCentralityClasses->GetSize()-1, temp);
163     fOutputList->Add(fProfV2);
164     fOutputList->Add(fProfV3);
165     switch (fFitModulationType) {
166         case kQC2 : {
167             fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
168             fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
169             fOutputList->Add(fProfV2Cumulant);
170             fOutputList->Add(fProfV3Cumulant);
171         } break;
172         case kQC4 : {
173             fProfV2Cumulant = new TProfile("fProfV2Cumulant", "fProfV2Cumulant", fCentralityClasses->GetSize()-1, temp);
174             fProfV3Cumulant = new TProfile("fProfV3Cumulant", "fProfV3Cumulant", fCentralityClasses->GetSize()-1, temp);
175             fOutputList->Add(fProfV2Cumulant);
176             fOutputList->Add(fProfV3Cumulant);
177         } break;
178         default : break;
179     }
180     if(fUsePtWeight) fHistSwap->Sumw2();
181     if(fUserSuppliedV2) fOutputList->Add(fUserSuppliedV2);
182     if(fUserSuppliedV3) fOutputList->Add(fUserSuppliedV3);
183     if(fUserSuppliedR2) fOutputList->Add(fUserSuppliedR2);
184     if(fUserSuppliedR3) fOutputList->Add(fUserSuppliedR3);
185         // increase readability of output list
186     fOutputList->Sort();
187     PostData(1, fOutputList);
188     switch (fRunModeType) {
189         case kLocal : {
190             fOutputListGood = new TList();
191             fOutputListGood->SetOwner(kTRUE);
192             fOutputListBad = new TList();
193             fOutputListBad->SetOwner(kTRUE);
194             PostData(2, fOutputListGood);
195             PostData(3, fOutputListBad);
196         } break;
197         default: break;
198     }
199 }
200 //_____________________________________________________________________________
201 TH1F* AliAnalysisTaskLocalRho::BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Int_t c, Bool_t append)
202 {
203     // book a TH1F and connect it to the output container
204     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
205     if(!fOutputList) return 0x0;
206     TString title(name);
207     if(c!=-1) { // format centrality dependent histograms accordingly
208         name = Form("%s_%i", name, c);
209         title += Form("_%i-%i", fCentralityClasses->At(c), fCentralityClasses->At(1+c));
210     }
211     title += Form(";%s;[counts]", x);
212     TH1F* histogram = new TH1F(name, title.Data(), bins, min, max);
213     histogram->Sumw2();
214     if(append) fOutputList->Add(histogram);
215     return histogram;   
216 }
217 //_____________________________________________________________________________
218 TH2F* AliAnalysisTaskLocalRho::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)
219 {
220     // book a TH2F and connect it to the output container
221     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
222     if(!fOutputList) return 0x0;
223     TString title(name);
224     if(c!=-1) { // format centrality dependent histograms accordingly
225         name = Form("%s_%i", name, c);
226         title += Form("_%i-%i", fCentralityClasses->At(c), fCentralityClasses->At(1+c));
227     }
228     title += Form(";%s;%s", x, y);
229     TH2F* histogram = new TH2F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy);
230     histogram->Sumw2();
231     if(append) fOutputList->Add(histogram);
232     return histogram;   
233 }
234 //_____________________________________________________________________________
235 Bool_t AliAnalysisTaskLocalRho::Run()
236 {
237     // execute once for each event
238     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
239     if(!(InputEvent()||fTracks||fJets||fRho)) return kFALSE;
240     if(!fInitialized) fInitialized = InitializeAnalysis();
241     // get the centrality bin (necessary for some control histograms
242     Double_t cent(InputEvent()->GetCentrality()->GetCentralityPercentile("V0M"));
243     for(Int_t i(0); i < fCentralityClasses->GetSize()-1; i++) {
244         if(cent >= fCentralityClasses->At(i) && cent <= fCentralityClasses->At(1+i)) {
245             fInCentralitySelection = i;
246             break; }
247     }
248     // set the rho value 
249     fLocalRho->SetVal(fRho->GetVal());
250     // set the correct event plane accordign to the requested reference detector
251     Double_t psi2(-1), psi3(-1);
252     switch (fDetectorType) {    // determine the detector type for the rho fit
253         case kTPC :     { 
254             // [0] psi2         [1] psi3
255             Double_t tpc[2];
256             CalculateEventPlaneTPC(tpc);
257             psi2 = tpc[0];         psi3 = tpc[1]; 
258         } break;
259         case kVZEROA :  { 
260             // [0][0] psi2a     [1,0]   psi2c
261             // [0][1] psi3a     [1,1]   psi3c
262             Double_t vzero[2][2];
263             CalculateEventPlaneVZERO(vzero);
264             psi2 = vzero[0][0];    psi3 = vzero[0][1]; 
265         }   break;  
266         case kVZEROC :  { 
267             // [0][0] psi2a     [1,0]   psi2c
268             // [0][1] psi3a     [1,1]   psi3c
269             Double_t vzero[2][2];
270             CalculateEventPlaneVZERO(vzero);
271             psi2 = vzero[1][0];    psi3 = vzero[1][1]; 
272         }   break;
273         case kVZEROComb : { 
274              /* for the combined vzero event plane
275              * [0] psi2         [1] psi3
276              * not fully implmemented yet, use with caution ! */
277              Double_t vzeroComb[2];
278              CalculateEventPlaneCombinedVZERO(vzeroComb);
279              psi2 = vzeroComb[0]; psi3 = vzeroComb[1];
280         } break;
281         default : break;
282     }
283     if(fFillHistograms) FillEventPlaneHistograms(psi2, psi3);
284     switch (fFitModulationType) { // do the fits
285         case kNoFit : { fFitModulation->FixParameter(0, fLocalRho->GetVal()); } break;
286         case kV2 : {    // only v2
287             if(CorrectRho(psi2, psi3)) {
288                 if(fFillHistograms) fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
289                 if(fUserSuppliedR2) {
290                     Double_t r(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
291                     if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
292                 }
293             }
294         } break;
295         case kV3 : {    // only v3
296             if(CorrectRho(psi2, psi3)) {
297                 if(fUserSuppliedR3) {
298                     Double_t r(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
299                     if(r > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r);
300                 }
301                 if(fFillHistograms) fProfV3->Fill(fCent, fFitModulation->GetParameter(3));
302             }
303         } break;
304         case kQC2 : {   // qc2 analysis - NOTE: not a wise idea to use this !
305             if(CorrectRho(psi2, psi3)) {
306                 if(fUserSuppliedR2 && fUserSuppliedR3) {
307                     // note for the qc method, resolution is REVERSED to go back to v2obs
308                     Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
309                     Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
310                     if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
311                     if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(3)*r3);
312                 }
313                 if (fUsePtWeight) { // use weighted weights
314                     Double_t dQCnM11 = (fNoEventWeightsForQC) ? 1. : QCnM11();
315                     if(fFillHistograms) {
316                         fProfV2->Fill(fCent, fFitModulation->GetParameter(3), dQCnM11);
317                         fProfV3->Fill(fCent, fFitModulation->GetParameter(7), dQCnM11); 
318                     }
319                 } else {
320                     Double_t dQCnM = (fNoEventWeightsForQC) ? 2. : QCnM();
321                     if(fFillHistograms) {
322                         fProfV2->Fill(fCent, fFitModulation->GetParameter(3), dQCnM*(dQCnM-1));
323                         fProfV3->Fill(fCent, fFitModulation->GetParameter(7), dQCnM*(dQCnM-1));
324                     }
325                 }
326             }
327         } break;
328         case kQC4 : {   // NOTE: see comment at kQC2
329             if(CorrectRho(psi2, psi3)) {
330                 if(fUserSuppliedR2 && fUserSuppliedR3) {
331                     // note for the qc method, resolution is REVERSED to go back to v2obs   
332                     Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
333                     Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
334                     if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)*r2);
335                     if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(3)*r3);
336                 }
337                 if (fUsePtWeight) { // use weighted weights
338                     if(fFillHistograms) {
339                         fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM1111()*/);
340                         fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM1111()*/); 
341                     }
342                 } else {
343                     if(fFillHistograms) {
344                         fProfV2->Fill(fCent, TMath::Power(fFitModulation->GetParameter(3),0.5)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
345                     fProfV3->Fill(fCent, TMath::Power(fFitModulation->GetParameter(7),0.5)/*, QCnM()*(QCnM()-1)*(QCnM()-2)*(QCnM()-3)*/);
346                     }
347                 }
348             }
349         } break;
350         default : {
351             if(CorrectRho(psi2, psi3)) {
352                 if(fUserSuppliedR2 && fUserSuppliedR3) {
353                     Double_t r2(fUserSuppliedR2->GetBinContent(fUserSuppliedR2->GetXaxis()->FindBin(fCent)));
354                     Double_t r3(fUserSuppliedR3->GetBinContent(fUserSuppliedR3->GetXaxis()->FindBin(fCent)));
355                     if(r2 > 0) fFitModulation->SetParameter(3, fFitModulation->GetParameter(3)/r2);
356                     if(r3 > 0) fFitModulation->SetParameter(7, fFitModulation->GetParameter(7)/r3);
357                 }
358                 if(fFillHistograms) {
359                     fProfV2->Fill(fCent, fFitModulation->GetParameter(3));
360                     fProfV3->Fill(fCent, fFitModulation->GetParameter(7));
361                 }
362             }
363         } break;
364     }
365     // if all went well, add local rho
366     fLocalRho->SetLocalRho(fFitModulation);
367     PostData(1, fOutputList);
368     return kTRUE;
369 }
370 //_____________________________________________________________________________
371 void AliAnalysisTaskLocalRho::CalculateEventPlaneVZERO(Double_t vzero[2][2]) const 
372 {
373     // get the vzero event plane
374     if(fUseV0EventPlaneFromHeader) {
375         // use the vzero event plane from the event header
376         // note: to use the calibrated vzero event plane, run 
377         // $ALICE_ROOT/ANALYSIS/macros/AddTaskVZEROEPSelection.C
378         // prior to this task (make sure the calibration is available for the dataset
379         // you want to use)
380         Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0), h(0);
381         vzero[0][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 2, a, b);
382         vzero[1][0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 2, c, d);
383         vzero[0][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 3, e, f);
384         vzero[1][1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 3, g, h);
385         return;
386     }
387     // grab the vzero event plane without recentering
388     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
389     Double_t qxa2(0), qya2(0), qxc2(0), qyc2(0);    // for psi2
390     Double_t qxa3(0), qya3(0), qxc3(0), qyc3(0);    // for psi3
391     for(Int_t iVZERO(0); iVZERO < 64; iVZERO++) {
392         Double_t phi(TMath::PiOver4()*(.5+iVZERO%8)), /* eta(0), */ weight(InputEvent()->GetVZEROEqMultiplicity(iVZERO));
393 //        (iVZERO<32) ? eta = -3.45+.5*(iVZERO/8) : eta = 4.8-.6*((iVZERO/8)-4);
394         if(iVZERO<32) {
395             qxa2 += weight*TMath::Cos(2.*phi);
396             qya2 += weight*TMath::Sin(2.*phi);
397             qxa3 += weight*TMath::Cos(3.*phi);
398             qya3 += weight*TMath::Sin(3.*phi);
399         }
400         else {
401             qxc2 += weight*TMath::Cos(2.*phi);
402             qyc2 += weight*TMath::Sin(2.*phi);
403             qxc3 += weight*TMath::Cos(3.*phi);
404             qyc3 += weight*TMath::Sin(3.*phi);
405        }
406     }
407     vzero[0][0] = .5*TMath::ATan2(qya2, qxa2);
408     vzero[1][0] = .5*TMath::ATan2(qyc2, qxc2);
409     vzero[0][1] = (1./3.)*TMath::ATan2(qya3, qxa3);
410     vzero[1][1] = (1./3.)*TMath::ATan2(qyc3, qxc3);
411 }
412 //_____________________________________________________________________________
413 void AliAnalysisTaskLocalRho::CalculateEventPlaneTPC(Double_t* tpc)
414 {
415    // grab the TPC event plane. if parameter fExcludeLeadingJetsFromFit is larger than 0, 
416    // strip in eta of width fExcludeLeadingJetsFromFit * fJetRadius around the leading jet (before
417    // subtraction of rho) will be exluded from the event plane estimate
418    if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
419    fNAcceptedTracks = 0;                // reset the track counter
420    Double_t qx2(0), qy2(0);     // for psi2
421    Double_t qx3(0), qy3(0);     // for psi3
422    if(fTracks) {
423        Float_t excludeInEta[] = {-999, -999};
424        if(fExcludeLeadingJetsFromFit > 0 ) {    // remove the leading jet from ep estimate
425            AliEmcalJet* leadingJet[] = {0x0, 0x0};
426            static Int_t lJets[9999] = {-1};
427            GetSortedArray(lJets, fJets);
428            for(Int_t i(0); i < fJets->GetEntriesFast(); i++) {     // get the two leading jets
429                if (1 + i > fJets->GetEntriesFast()) break;
430                leadingJet[0] = static_cast<AliEmcalJet*>(fJets->At(lJets[i]));
431                leadingJet[1] = static_cast<AliEmcalJet*>(fJets->At(lJets[i+1]));
432                if(PassesCuts(leadingJet[0]) && PassesCuts(leadingJet[1])) break;
433            }
434            if(leadingJet[0] && leadingJet[1]) {
435                for(Int_t i(0); i < 2; i++) excludeInEta[i] = leadingJet[i]->Eta();
436            }
437        }
438        Int_t iTracks(fTracks->GetEntriesFast());
439        for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
440            AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
441            if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
442            if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta[0]) < fJetRadius*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - fJetRadius - fJetMaxEta ) > 0 )) continue;
443            fNAcceptedTracks++;
444            qx2+= TMath::Cos(2.*track->Phi());
445            qy2+= TMath::Sin(2.*track->Phi());
446            qx3+= TMath::Cos(3.*track->Phi());
447            qy3+= TMath::Sin(3.*track->Phi());
448        }
449    }
450    tpc[0] = .5*TMath::ATan2(qy2, qx2);
451    tpc[1] = (1./3.)*TMath::ATan2(qy3, qx3);
452
453 //_____________________________________________________________________________
454 void AliAnalysisTaskLocalRho::CalculateEventPlaneCombinedVZERO(Double_t* comb) const
455 {
456     // grab the combined vzero event plane
457 //    if(fUseV0EventPlaneFromHeader) {    // use the vzero from the header
458         Double_t a(0), b(0), c(0), d(0);
459         comb[0] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 10, 2, a, b);
460         comb[1] = InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 10, 3, c, d);
461 // FIXME the rest of this function isn't impelmented yet (as of 01-07-2013)
462 // this means a default the combined vzero event plane from the header is used
463 // to get this value 'by hand', vzeroa and vzeroc event planes have to be combined
464 // according to their resolution - this will be added ...
465 //
466 //    } else {
467 //        Double_t qx2a(0), qy2a(0), qx2c(0), qy2c(0), qx3a(0), qy3a(0), qx3c(0), qy3c(0);
468 //        InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 2, qx2a, qy2a);
469 //        InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 2, qx2c, qy2c);
470 //        InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 8, 3, qx3a, qy3a);
471 //        InputEvent()->GetEventplane()->CalculateVZEROEventPlane(InputEvent(), 9, 3, qx3c, qy3c);
472 //        Double_t chi2A(-1), chi2C(-1), chi3A(-1), chi3C(-1);     // get chi from the resolution
473 //        Double_t qx2(chi2A*chi2A*qx2a+chi2C*chi2C*qx2c);
474 //        Double_t qy2(chi2A*chi2A*qy2a+chi2C*chi2C*qy2c);
475 //        Double_t qx3(chi3A*chi3A*qx3a+chi3C*chi3C*qx3c);
476 //        Double_t qy3(chi3A*chi3A*qy3a+chi3C*chi3C*qy3c);
477 //        comb[0] = .5*TMath::ATan2(qy2, qx2);
478 //        comb[1] = (1./3.)*TMath::ATan2(qy3, qx3);
479 //    }
480 }
481 //_____________________________________________________________________________
482 Double_t AliAnalysisTaskLocalRho::CalculateQC2(Int_t harm) {
483     // get the second order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho
484     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
485     Double_t reQ(0), imQ(0), modQ(0), M11(0), M(0);
486     if(fUsePtWeight) {  // for the weighted 2-nd order q-cumulant
487         QCnQnk(harm, 1, reQ, imQ);      // get the weighted 2-nd order q-vectors
488         modQ = reQ*reQ+imQ*imQ;         // get abs Q-squared
489         M11 = QCnM11();                 // equals S2,1 - S1,2
490         return (M11 > 0) ? ((modQ - QCnS(1,2))/M11) : -999;
491     } // else return the non-weighted 2-nd order q-cumulant
492     QCnQnk(harm, 0, reQ, imQ);          // get the non-weighted 2-nd order q-vectors
493     modQ = reQ*reQ+imQ*imQ;             // get abs Q-squared
494     M = QCnM();
495     return (M > 1) ? (modQ - M)/(M*(M-1)) : -999;
496 }
497 //_____________________________________________________________________________
498 Double_t AliAnalysisTaskLocalRho::CalculateQC4(Int_t harm) {
499     // get the fourth order q-cumulant, a -999 return will be caught in the qa routine of CorrectRho
500     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
501     Double_t reQn1(0), imQn1(0), reQ2n2(0), imQ2n2(0), reQn3(0), imQn3(0), M1111(0), M(0);
502     Double_t a(0), b(0), c(0), d(0), e(0), f(0), g(0);  // terms of the calculation
503     if(fUsePtWeight) {  // for the weighted 4-th order q-cumulant
504         QCnQnk(harm, 1, reQn1, imQn1);
505         QCnQnk(harm*2, 2, reQ2n2, imQ2n2);
506         QCnQnk(harm, 3, reQn3, imQn3);
507         // fill in the terms ...
508         a = (reQn1*reQn1+imQn1*imQn1)*(reQn1*reQn1+imQn1*imQn1);
509         b = reQ2n2*reQ2n2 + imQ2n2*imQ2n2;
510         c = -2.*(reQ2n2*reQn1*reQn1-reQ2n2*imQn1*imQn1+2.*imQ2n2*reQn1*imQn1);
511         d = 8.*(reQn3*reQn1+imQn3*imQn1);
512         e = -4.*QCnS(1,2)*(reQn1*reQn1+imQn1*imQn1);
513         f = -6.*QCnS(1,4);
514         g = 2.*QCnS(2,2);
515         M1111 = QCnM1111();
516         return (M1111 > 0) ? (a+b+c+d+e+f+g)/M1111 : -999;
517     }   // else return the unweighted case
518     Double_t reQn(0), imQn(0), reQ2n(0), imQ2n(0);
519     QCnQnk(harm, 0, reQn, imQn);
520     QCnQnk(harm*2, 0, reQ2n, imQ2n);
521     // fill in the terms ...
522     M = QCnM();
523     if(M < 4) return -999;
524     a = (reQn*reQn+imQn*imQn)*(reQn*reQn+imQn*imQn);
525     b = reQ2n*reQ2n + imQ2n*imQ2n;
526     c = -2.*(reQ2n*reQn*reQn-reQ2n*imQn*imQn+2.*imQ2n*reQn*imQn);
527     e = -4.*(M-2)*(reQn*reQn+imQn*imQn);
528     f = 2.*M*(M-3);
529     return (a+b+c+e+f)/(M*(M-1)*(M-2)*(M-3));
530 }
531 //_____________________________________________________________________________
532 void AliAnalysisTaskLocalRho::QCnQnk(Int_t n, Int_t k, Double_t &reQ, Double_t &imQ) {
533     // get the weighted n-th order q-vector, pass real and imaginary part as reference
534     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
535     if(!fTracks) return;
536     fNAcceptedTracksQCn = 0;
537     Int_t iTracks(fTracks->GetEntriesFast());
538     for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
539         AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
540         if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
541         fNAcceptedTracksQCn++;
542         // for the unweighted case, k equals zero and the weight doesn't contribute to the equation below
543         reQ += TMath::Power(track->Pt(), k) * TMath::Cos(((double)n)*track->Phi());
544         imQ += TMath::Power(track->Pt(), k) * TMath::Sin(((double)n)*track->Phi());
545     }
546 }
547 //_____________________________________________________________________________
548 Double_t AliAnalysisTaskLocalRho::QCnS(Int_t i, Int_t j) {
549     // get the weighted ij-th order autocorrelation correction
550     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
551     if(!fTracks || i <= 0 || j <= 0) return -999;
552     Int_t iTracks(fTracks->GetEntriesFast());
553     Double_t Sij(0);
554     for(Int_t iTPC(0); iTPC < iTracks; iTPC++) {
555         AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(iTPC));
556         if(!PassesCuts(track) || track->Pt() < fSoftTrackMinPt || track->Pt() > fSoftTrackMaxPt) continue;
557         Sij+=TMath::Power(track->Pt(), j);
558     }
559     return TMath::Power(Sij, i);
560 }
561 //_____________________________________________________________________________
562 Double_t AliAnalysisTaskLocalRho::QCnM() {
563     // get multiplicity for unweighted q-cumulants. function QCnQnk should be called first
564     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
565     return (Double_t) fNAcceptedTracksQCn;
566 }
567 //_____________________________________________________________________________
568 Double_t AliAnalysisTaskLocalRho::QCnM11() {
569     // get multiplicity weights for the weighted two particle cumulant
570     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
571     return (QCnS(2,1) - QCnS(1,2));
572 }
573 //_____________________________________________________________________________
574 Double_t AliAnalysisTaskLocalRho::QCnM1111() {
575     // get multiplicity weights for the weighted four particle cumulant
576     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
577     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));
578 }
579 //_____________________________________________________________________________
580 Bool_t AliAnalysisTaskLocalRho::QCnRecovery(Double_t psi2, Double_t psi3) {
581     // decides how to deal with the situation where c2 or c3 is negative 
582     // returns kTRUE depending on whether or not a modulated rho is used for the jet background
583     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
584     if(TMath::AreEqualAbs(fFitModulation->GetParameter(3), .0, 1e-10) && TMath::AreEqualAbs(fFitModulation->GetParameter(7), .0,1e-10)) {
585         fFitModulation->SetParameter(7, 0);
586         fFitModulation->SetParameter(3, 0);
587         fFitModulation->SetParameter(0, fLocalRho->GetVal());
588         return kTRUE;   // v2 and v3 have physical null values
589     }
590     switch (fQCRecovery) {
591         case kFixedRho : {      // roll back to the original rho
592            fFitModulation->SetParameter(7, 0);
593            fFitModulation->SetParameter(3, 0);
594            fFitModulation->SetParameter(0, fLocalRho->GetVal());
595            return kFALSE;       // rho is forced to be fixed
596         }
597         case kNegativeVn : {
598            Double_t c2(fFitModulation->GetParameter(3));
599            Double_t c3(fFitModulation->GetParameter(7));
600            if( c2 < 0 ) c2 = -1.*TMath::Sqrt(-1.*c2);
601            if( c3 < 0 ) c3 = -1.*TMath::Sqrt(-1.*c3);
602            fFitModulation->SetParameter(3, c2);
603            fFitModulation->SetParameter(7, c3);
604            return kTRUE;        // is this a physical quantity ?
605         }
606         case kTryFit : {
607            fitModulationType tempType(fFitModulationType);  // store temporarily
608            fFitModulationType = kCombined;
609            fFitModulation->SetParameter(7, 0);
610            fFitModulation->SetParameter(3, 0);
611            Bool_t pass(CorrectRho(psi2, psi3));         // do the fit and all quality checks
612            fFitModulationType = tempType;               // roll back for next event
613            return pass;
614         }
615         default : return kFALSE;
616     }
617     return kFALSE;
618 }
619 //_____________________________________________________________________________
620 Bool_t AliAnalysisTaskLocalRho::CorrectRho(Double_t psi2, Double_t psi3) 
621 {
622     // get rho' -> rho(phi)
623     // three routines are available, 1 and 2 can be used with or without pt weights
624     //  [1] get vn from q-cumulants
625     //      in case of cumulants, both cumulants and vn values are stored. in both cases, v2 and v3
626     //      are expected. a check is performed to see if rho has no negative local minimum
627     //      for full description, see Phys. Rev. C 83, 044913
628     //      since the cn distribution has negative values, vn = sqrt(cn) can be imaginary sometimes
629     //      in this case one can either roll back to the 'original' fixed rho, do a fit for vn or take use
630     //      vn = - sqrt(|cn|) note that because of this, use of q-cumulants is not safe ! 
631     //  [2] fitting a fourier expansion to the de/dphi distribution
632     //      the fit can be done with either v2, v3 or a combination.
633     //      in all cases, a cut can be made on the p-value of the chi-squared value of the fit
634     //      and a check can be performed to see if rho has no negative local minimum
635     //  [3] get v2 and v3 from user supplied histograms
636     //      in this way, a fixed value of v2 and v3 is subtracted w.r.t. whichever event plane is requested
637     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
638     switch (fFitModulationType) {       // for approaches where no fitting is required
639         case kQC2 : {
640             fFitModulation->FixParameter(4, psi2); 
641             fFitModulation->FixParameter(6, psi3);
642             fFitModulation->FixParameter(3, CalculateQC2(2));   // set here with cn, vn = sqrt(cn)
643             fFitModulation->FixParameter(7, CalculateQC2(3));
644             // first fill the histos of the raw cumulant distribution
645             if (fUsePtWeight) { // use weighted weights
646                 Double_t dQCnM11 = (fNoEventWeightsForQC) ? 1. : QCnM11();
647                 if(fFillHistograms) {
648                     fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), dQCnM11);
649                     fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), dQCnM11);
650                 }
651             } else {
652                 Double_t dQCnM = (fNoEventWeightsForQC) ? 2. : QCnM();
653                 if(fFillHistograms) {
654                     fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3), dQCnM*(dQCnM-1));
655                     fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7), dQCnM*(dQCnM-1));
656                 }
657             }
658             // then see if one of the cn value is larger than zero and vn is readily available
659             if(fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) {
660                 fFitModulation->FixParameter(3, TMath::Sqrt(fFitModulation->GetParameter(3)));
661                 fFitModulation->FixParameter(7, TMath::Sqrt(fFitModulation->GetParameter(7)));
662             } else if (!QCnRecovery(psi2, psi3)) return kFALSE;  // try to recover the cumulant, this will set v2 and v3
663             if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) {  // general check 
664                 fFitModulation->SetParameter(7, 0);
665                 fFitModulation->SetParameter(3, 0);
666                 fFitModulation->SetParameter(0, fLocalRho->GetVal());
667                 return kFALSE;
668             }
669             return kTRUE;
670         } break;
671         case kQC4 : {
672             fFitModulation->FixParameter(4, psi2); 
673             fFitModulation->FixParameter(6, psi3);
674             fFitModulation->FixParameter(3, CalculateQC4(2));   // set here with cn, vn = sqrt(cn)
675             fFitModulation->FixParameter(7, CalculateQC4(3));
676             // first fill the histos of the raw cumulant distribution
677             if (fUsePtWeight) { // use weighted weights
678                 if(fFillHistograms) {
679                     fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
680                     fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
681                 }
682             } else {
683                 if(fFillHistograms) {
684                     fProfV2Cumulant->Fill(fCent, fFitModulation->GetParameter(3)/*, QCnM1111()*/);
685                     fProfV3Cumulant->Fill(fCent, fFitModulation->GetParameter(7)/*, QCnM1111()*/);
686                 }
687             }
688             // then see if one of the cn value is larger than zero and vn is readily available
689             if(fFitModulation->GetParameter(3) > 0 && fFitModulation->GetParameter(7) > 0) {
690                 fFitModulation->FixParameter(3, TMath::Sqrt(fFitModulation->GetParameter(3)));
691                 fFitModulation->FixParameter(7, TMath::Sqrt(fFitModulation->GetParameter(7)));
692             } else if (!QCnRecovery(psi2, psi3)) return kFALSE;  // try to recover the cumulant, this will set v2 and v3
693             if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) {  // general check 
694                 fFitModulation->SetParameter(7, 0);
695                 fFitModulation->SetParameter(3, 0);
696                 fFitModulation->SetParameter(0, fLocalRho->GetVal());
697                 return kFALSE;
698             }
699         } break;
700         case kIntegratedFlow : {
701             // use v2 and v3 values from an earlier iteration over the data
702             fFitModulation->FixParameter(3, fUserSuppliedV2->GetBinContent(fUserSuppliedV2->GetXaxis()->FindBin(fCent)));
703             fFitModulation->FixParameter(4, psi2);
704             fFitModulation->FixParameter(6, psi3);
705             fFitModulation->FixParameter(7, fUserSuppliedV3->GetBinContent(fUserSuppliedV3->GetXaxis()->FindBin(fCent)));
706             if(fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) < 0) { 
707                 fFitModulation->SetParameter(7, 0);
708                 fFitModulation->SetParameter(3, 0);
709                 fFitModulation->SetParameter(0, fLocalRho->GetVal());
710                 return kFALSE;
711             }
712             return kTRUE;
713         }
714         default : break;
715     }
716     TString detector("");
717     switch (fDetectorType) {
718         case kTPC : detector+="TPC";
719             break;
720         case kVZEROA : detector+="VZEROA";
721             break;
722         case kVZEROC : detector+="VZEROC";
723             break;
724         case kVZEROComb : detector+="VZEROComb";
725             break; 
726         default: break;
727     }
728     Int_t iTracks(fTracks->GetEntriesFast());
729     Double_t excludeInEta[] = {-999, -999};
730     Double_t excludeInPhi[] = {-999, -999};
731     Double_t excludeInPt[]  = {-999, -999};
732     if(iTracks <= 0 || fLocalRho->GetVal() <= 0 ) return kFALSE;   // no use fitting an empty event ...
733     if(fExcludeLeadingJetsFromFit > 0 ) {
734         AliEmcalJet* leadingJet[] = {0x0, 0x0};
735         static Int_t lJets[9999] = {-1};
736         GetSortedArray(lJets, fJets);
737         for(Int_t i(0); i < fJets->GetEntriesFast(); i++) {     // get the two leading jets
738             if (1 + i > fJets->GetEntriesFast()) break;
739             leadingJet[0] = static_cast<AliEmcalJet*>(fJets->At(lJets[i]));
740             leadingJet[1] = static_cast<AliEmcalJet*>(fJets->At(lJets[i+1]));
741             if(PassesCuts(leadingJet[0]) && PassesCuts(leadingJet[1])) break;
742         }
743         if(leadingJet[0] && leadingJet[1]) {
744             for(Int_t i(0); i < 2; i++) {
745                 excludeInEta[i] = leadingJet[i]->Eta();
746                 excludeInPhi[i] = leadingJet[i]->Phi();
747                 excludeInPt[i]  = leadingJet[i]->Pt();
748             }
749         }
750     }
751     fHistSwap->Reset();                 // clear the histogram
752     TH1F _tempSwap;
753     if(fRebinSwapHistoOnTheFly) {
754         if(fNAcceptedTracks < 49) fNAcceptedTracks = 49;       // avoid aliasing effects
755         _tempSwap = TH1F("_tempSwap", "_tempSwap", TMath::CeilNint(TMath::Sqrt(fNAcceptedTracks)), 0, TMath::TwoPi());
756         if(fUsePtWeight) _tempSwap.Sumw2();
757     }
758     else _tempSwap = *fHistSwap;         // now _tempSwap holds the desired histo
759     for(Int_t i(0); i < iTracks; i++) {
760             AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
761             if(fExcludeLeadingJetsFromFit > 0 &&( (TMath::Abs(track->Eta() - excludeInEta[0]) < fJetRadius*fExcludeLeadingJetsFromFit ) || (TMath::Abs(track->Eta()) - fJetRadius - fJetMaxEta ) > 0 )) continue;
762             if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
763             if(fUsePtWeight) _tempSwap.Fill(track->Phi(), track->Pt());
764             else _tempSwap.Fill(track->Phi());
765     }
766 //    for(Int_t i(0); i < _tempSwap.GetXaxis()->GetNbins(); i++) _tempSwap.SetBinError(1+i, TMath::Sqrt(_tempSwap.GetBinContent(1+i)));
767     fFitModulation->SetParameter(0, fLocalRho->GetVal());
768     switch (fFitModulationType) {
769         case kNoFit : { fFitModulation->FixParameter(0, fLocalRho->GetVal() ); 
770         } break;
771         case kV2 : { 
772             fFitModulation->FixParameter(4, psi2); 
773         } break;
774         case kV3 : { 
775             fFitModulation->FixParameter(4, psi3); 
776         } break;
777         case kCombined : {
778             fFitModulation->FixParameter(4, psi2); 
779             fFitModulation->FixParameter(6, psi3);
780         } break;
781         case kFourierSeries : {
782             // in this approach, an explicit calculation will be made of vn = sqrt(xn^2+yn^2)
783             // where x[y] = Integrate[r(phi)cos[sin](n phi)dphi, 0, 2pi]
784             Double_t cos2(0), sin2(0), cos3(0), sin3(0), sumPt(0);
785             for(Int_t i(0); i < iTracks; i++) {
786                 AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
787                 if(!PassesCuts(track) || track->Pt() > fSoftTrackMaxPt || track->Pt() < fSoftTrackMinPt) continue;
788                 sumPt += track->Pt();
789                 cos2 += track->Pt()*TMath::Cos(2*PhaseShift(track->Phi()-psi2)); 
790                 sin2 += track->Pt()*TMath::Sin(2*PhaseShift(track->Phi()-psi2));
791                 cos3 += track->Pt()*TMath::Cos(3*PhaseShift(track->Phi()-psi3)); 
792                 sin3 += track->Pt()*TMath::Sin(3*PhaseShift(track->Phi()-psi3));
793             }
794             fFitModulation->SetParameter(3, TMath::Sqrt(cos2*cos2+sin2*sin2)/fLocalRho->GetVal());
795             fFitModulation->SetParameter(4, psi2);
796             fFitModulation->SetParameter(6, psi3);
797             fFitModulation->SetParameter(7, TMath::Sqrt(cos3*cos3+sin3*sin3)/fLocalRho->GetVal());
798         } break;
799         default : break;
800     }
801     _tempSwap.Fit(fFitModulation, fFitModulationOptions.Data(), "", 0, TMath::TwoPi());
802     // the quality of the fit is evaluated from 1 - the cdf of the chi square distribution
803     Double_t CDF(1.-ChiSquareCDF(fFitModulation->GetNDF(), fFitModulation->GetChisquare()));
804     if(fFillHistograms) fHistPvalueCDF->Fill(CDF);
805     if(CDF > fMinPvalue && CDF < fMaxPvalue && ( fAbsVnHarmonics && fFitModulation->GetMinimum(0, TMath::TwoPi()) > 0)) { // fit quality
806         // for LOCAL didactic purposes, save the  best and the worst fits
807         // this routine can produce a lot of output histograms (it's not memory 'safe') and will not work on GRID 
808         // since the output will become unmergeable (i.e. different nodes may produce conflicting output)
809         switch (fRunModeType) {
810             case kLocal : {
811                 if(gRandom->Uniform(0, 100) > fPercentageOfFits) break;
812                 static Int_t didacticCounterBest(0);
813                 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
814                 TF1* didactifFit = (TF1*)fFitModulation->Clone(Form("fit_%i_CDF_%.3f_cen_%i_%s", didacticCounterBest, CDF, fInCentralitySelection, detector.Data()));
815                 didacticProfile->GetListOfFunctions()->Add(didactifFit);
816                 fOutputListGood->Add(didacticProfile);
817                 didacticCounterBest++;
818                 TH2F* didacticSurface = BookTH2F(Form("surface_%s", didacticProfile->GetName()), "#phi", "#eta", 50, 0, TMath::TwoPi(), 50, -1, 1, -1, kFALSE);
819                 for(Int_t i(0); i < iTracks; i++) {
820                     AliVTrack* track = static_cast<AliVTrack*>(fTracks->At(i));
821                     if(PassesCuts(track)) {
822                         if(fUsePtWeight) didacticSurface->Fill(track->Phi(), track->Eta(), track->Pt());
823                         else didacticSurface->Fill(track->Phi(), track->Eta());
824                     }
825                 }
826                 if(fExcludeLeadingJetsFromFit) {       // visualize the excluded region
827                     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);
828                     f2->SetParameters(excludeInPt[0]/3.,excludeInPhi[0],.1,excludeInEta[0],.1);
829                     didacticSurface->GetListOfFunctions()->Add(f2);
830                     TF2 *f3 = new TF2(Form("%s_NLJ", didacticSurface->GetName()),"[0]*TMath::Gaus(x,[1],[2])*TMath::Gaus(y,[3],[4])", 0, TMath::TwoPi(), -1, 1);
831                     f3->SetParameters(excludeInPt[1]/3.,excludeInPhi[1],.1,excludeInEta[1],.1);
832                     f3->SetLineColor(kGreen);
833                     didacticSurface->GetListOfFunctions()->Add(f3);
834                 }
835                 fOutputListGood->Add(didacticSurface);
836             } break;
837             default : break;
838         }
839     } else {    // if the fit is of poor quality revert to the original rho estimate
840         switch (fRunModeType) { // again see if we want to save the fit
841             case kLocal : {
842                 static Int_t didacticCounterWorst(0);
843                 if(gRandom->Uniform(0, 100) > fPercentageOfFits) break;
844                 TProfile* didacticProfile = (TProfile*)_tempSwap.Clone(Form("Fit_%i_1-CDF_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data() ));
845                 TF1* didactifFit = (TF1*)fFitModulation->Clone(Form("fit_%i_p_%.3f_cen_%i_%s", didacticCounterWorst, CDF, fInCentralitySelection, detector.Data()));
846                 didacticProfile->GetListOfFunctions()->Add(didactifFit);
847                 fOutputListBad->Add(didacticProfile);
848                 didacticCounterWorst++;
849                 } break;
850             default : break;
851         }
852         switch (fFitModulationType) {
853             case kNoFit : break;        // nothing to do
854             case kCombined : fFitModulation->SetParameter(7, 0);        // no break
855             case kFourierSeries : fFitModulation->SetParameter(7, 0);   // no break
856             default : { // needs to be done if there was a poor fit
857                  fFitModulation->SetParameter(3, 0);
858                  fFitModulation->SetParameter(0, fLocalRho->GetVal());
859             } break;
860         }
861         return kFALSE;  // return false if the fit is rejected
862     }
863     return kTRUE;
864 }
865 //_____________________________________________________________________________
866 void AliAnalysisTaskLocalRho::FillAnalysisSummaryHistogram() const
867 {
868     // fill the analysis summary histrogram, saves all relevant analysis settigns
869     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
870     fHistAnalysisSummary->GetXaxis()->SetBinLabel(1, "fJetRadius"); 
871     fHistAnalysisSummary->SetBinContent(1, fJetRadius);
872     fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fPtBiasJetTrack");
873     fHistAnalysisSummary->SetBinContent(2, fPtBiasJetTrack);
874     fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fPtBiasJetClus");
875     fHistAnalysisSummary->SetBinContent(3, fPtBiasJetClus);
876     fHistAnalysisSummary->GetXaxis()->SetBinLabel(4, "fJetPtCut");
877     fHistAnalysisSummary->SetBinContent(4, fJetPtCut);
878     fHistAnalysisSummary->GetXaxis()->SetBinLabel(5, "fJetAreaCut");
879     fHistAnalysisSummary->SetBinContent(5, fJetAreaCut);
880     fHistAnalysisSummary->GetXaxis()->SetBinLabel(6, "fPercAreaCut");
881     fHistAnalysisSummary->SetBinContent(6, fPercAreaCut);
882     fHistAnalysisSummary->GetXaxis()->SetBinLabel(7, "fAreaEmcCut");
883     fHistAnalysisSummary->SetBinContent(7, fAreaEmcCut);
884     fHistAnalysisSummary->GetXaxis()->SetBinLabel(8, "fJetMinEta");
885     fHistAnalysisSummary->SetBinContent(8, fJetMinEta);
886     fHistAnalysisSummary->GetXaxis()->SetBinLabel(9, "fJetMaxEta");
887     fHistAnalysisSummary->SetBinContent(9, fJetMaxEta);
888     fHistAnalysisSummary->GetXaxis()->SetBinLabel(10, "fJetMinPhi");
889     fHistAnalysisSummary->SetBinContent(10, fJetMinPhi);
890     fHistAnalysisSummary->GetXaxis()->SetBinLabel(11, "fJetMaxPhi");
891     fHistAnalysisSummary->SetBinContent(11, fJetMaxPhi);
892     fHistAnalysisSummary->GetXaxis()->SetBinLabel(12, "fMaxClusterPt");
893     fHistAnalysisSummary->SetBinContent(12, fMaxClusterPt);
894     fHistAnalysisSummary->GetXaxis()->SetBinLabel(13, "fMaxTrackPt");
895     fHistAnalysisSummary->SetBinContent(13, fMaxTrackPt);
896     fHistAnalysisSummary->GetXaxis()->SetBinLabel(14, "fLeadingHadronType");
897     fHistAnalysisSummary->SetBinContent(14, fLeadingHadronType);
898     fHistAnalysisSummary->GetXaxis()->SetBinLabel(15, "fAnaType");
899     fHistAnalysisSummary->SetBinContent(15, fAnaType);
900     fHistAnalysisSummary->GetXaxis()->SetBinLabel(16, "fForceBeamType");
901     fHistAnalysisSummary->SetBinContent(16, fForceBeamType);
902     fHistAnalysisSummary->GetXaxis()->SetBinLabel(19, "fMinVz");
903     fHistAnalysisSummary->SetBinContent(19, fMinVz);
904     fHistAnalysisSummary->GetXaxis()->SetBinLabel(20, "fMaxVz");
905     fHistAnalysisSummary->SetBinContent(20, fMaxVz);
906     fHistAnalysisSummary->GetXaxis()->SetBinLabel(21, "fOffTrigger");
907     fHistAnalysisSummary->SetBinContent(21, fOffTrigger);
908     fHistAnalysisSummary->GetXaxis()->SetBinLabel(22, "fClusPtCut");
909     fHistAnalysisSummary->SetBinContent(22, fClusPtCut);
910     fHistAnalysisSummary->GetXaxis()->SetBinLabel(23, "fTrackPtCut");
911     fHistAnalysisSummary->SetBinContent(23, fTrackPtCut);
912     fHistAnalysisSummary->GetXaxis()->SetBinLabel(24, "fTrackMinEta");
913     fHistAnalysisSummary->SetBinContent(24, fTrackMinEta);
914     fHistAnalysisSummary->GetXaxis()->SetBinLabel(25, "fTrackMaxEta");
915     fHistAnalysisSummary->SetBinContent(25, fTrackMaxEta);
916     fHistAnalysisSummary->GetXaxis()->SetBinLabel(26, "fTrackMinPhi");
917     fHistAnalysisSummary->SetBinContent(26, fTrackMinPhi);
918     fHistAnalysisSummary->GetXaxis()->SetBinLabel(27, "fTrackMaxPhi");
919     fHistAnalysisSummary->SetBinContent(27, fTrackMaxPhi);
920     fHistAnalysisSummary->GetXaxis()->SetBinLabel(28, "fClusTimeCutLow");
921     fHistAnalysisSummary->SetBinContent(28, fClusTimeCutLow);
922     fHistAnalysisSummary->GetXaxis()->SetBinLabel(29, "fClusTimeCutUp");
923     fHistAnalysisSummary->SetBinContent(29, fClusTimeCutUp);
924     fHistAnalysisSummary->GetXaxis()->SetBinLabel(30, "fMinPtTrackInEmcal");
925     fHistAnalysisSummary->SetBinContent(30, fMinPtTrackInEmcal);
926     fHistAnalysisSummary->GetXaxis()->SetBinLabel(31, "fEventPlaneVsEmcal");
927     fHistAnalysisSummary->SetBinContent(31, fEventPlaneVsEmcal);
928     fHistAnalysisSummary->GetXaxis()->SetBinLabel(32, "fMinEventPlane");
929     fHistAnalysisSummary->SetBinContent(32, fMaxEventPlane);
930     fHistAnalysisSummary->GetXaxis()->SetBinLabel(34, "fitModulationType");
931     fHistAnalysisSummary->SetBinContent(34, (int)fFitModulationType);
932     fHistAnalysisSummary->GetXaxis()->SetBinLabel(35, "runModeType");
933     fHistAnalysisSummary->SetBinContent(35, (int)fRunModeType);
934     fHistAnalysisSummary->GetXaxis()->SetBinLabel(37, "iterator");
935     fHistAnalysisSummary->SetBinContent(37, 1.);
936     fHistAnalysisSummary->GetXaxis()->SetBinLabel(38, "fMinPvalue");
937     fHistAnalysisSummary->SetBinContent(38, fMinPvalue);
938     fHistAnalysisSummary->GetXaxis()->SetBinLabel(39, "fMaxPvalue");
939     fHistAnalysisSummary->SetBinContent(39, fMaxPvalue);
940     fHistAnalysisSummary->GetXaxis()->SetBinLabel(40, "fExcludeLeadingJetsFromFit");
941     fHistAnalysisSummary->SetBinContent(40, fExcludeLeadingJetsFromFit);
942     fHistAnalysisSummary->GetXaxis()->SetBinLabel(41, "fRebinSwapHistoOnTheFly");
943     fHistAnalysisSummary->SetBinContent(41, (int)fRebinSwapHistoOnTheFly);
944     fHistAnalysisSummary->GetXaxis()->SetBinLabel(42, "fUsePtWeight");
945     fHistAnalysisSummary->SetBinContent(42, (int)fUsePtWeight);
946     fHistAnalysisSummary->GetXaxis()->SetBinLabel(45, "fLocalJetMinEta");
947     fHistAnalysisSummary->SetBinContent(45,fLocalJetMinEta );
948     fHistAnalysisSummary->GetXaxis()->SetBinLabel(46, "fLocalJetMaxEta");
949     fHistAnalysisSummary->SetBinContent(46, fLocalJetMaxEta);
950     fHistAnalysisSummary->GetXaxis()->SetBinLabel(47, "fLocalJetMinPhi");
951     fHistAnalysisSummary->SetBinContent(47, fLocalJetMinPhi);
952     fHistAnalysisSummary->GetXaxis()->SetBinLabel(48, "fLocalJetMaxPhi");
953     fHistAnalysisSummary->SetBinContent(48, fLocalJetMaxPhi);
954     fHistAnalysisSummary->GetXaxis()->SetBinLabel(49, "fSoftTrackMinPt");
955     fHistAnalysisSummary->SetBinContent(49, fSoftTrackMinPt);
956     fHistAnalysisSummary->GetXaxis()->SetBinLabel(50, "fSoftTrackMaxPt");
957     fHistAnalysisSummary->SetBinContent(50, fSoftTrackMaxPt);
958 }
959 //_____________________________________________________________________________
960 void AliAnalysisTaskLocalRho::FillEventPlaneHistograms(Double_t psi2, Double_t psi3) const
961 {
962     // fill event plane histograms
963     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
964     fHistPsi2[fInCentralitySelection]->Fill(psi2);
965     fHistPsi3[fInCentralitySelection]->Fill(psi3);    
966 }
967 //_____________________________________________________________________________
968 void AliAnalysisTaskLocalRho::Terminate(Option_t *)
969 {
970     // terminate
971 }
972 //_____________________________________________________________________________
973 void AliAnalysisTaskLocalRho::SetModulationFit(TF1* fit) {
974     // Set function to fit modulation
975     if (fFitModulation) delete fFitModulation;
976     fFitModulation = fit; 
977 }