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