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