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