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