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