]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/FLOW/Tasks/AliAnalysisTaskPhiFlow.cxx
Carlos Perez: fix typo
[u/mrichter/AliRoot.git] / PWG / FLOW / Tasks / AliAnalysisTaskPhiFlow.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 // AliAnalysisTaskPhiFlow:
17 // author: Redmer Alexander Bertens (rbertens@cern.ch)
18 // analyis task for phi-meson reconstruction and determination of v_n
19 // AliPhiMesonHelperTrack provides a lightweight helper track for reconstruction
20 // new in this version (use with caution): vzero event plane, event mixing
21
22 #include "TChain.h"
23 #include "TTree.h"
24 #include "TH1.h"
25 #include "TH1F.h"
26 #include "TH2F.h"
27 #include "TCanvas.h"
28 #include "TMath.h"
29 #include "TObjArray.h"
30 #include "AliAnalysisTask.h"
31 #include "AliAnalysisTaskSE.h"
32 #include "AliAnalysisManager.h"
33 #include "AliESDEvent.h"
34 #include "AliAODEvent.h"
35 #include "AliESDInputHandler.h"
36 #include "AliAODInputHandler.h"
37 #include "AliCentrality.h"
38 #include "AliVEvent.h"
39 #include "AliAnalysisTaskPhiFlow.h"
40 #include "AliFlowBayesianPID.h"
41 #include "AliStack.h"
42 #include "AliESDtrackCuts.h"
43 #include "AliMCEvent.h"
44 #include "TProfile.h"
45 #include "AliFlowCandidateTrack.h"
46 #include "AliFlowTrackCuts.h"
47 #include "AliFlowEventSimple.h"
48 #include "AliFlowTrackSimple.h"
49 #include "AliFlowCommonConstants.h"
50 #include "AliFlowEvent.h"
51 #include "TVector3.h"
52 #include "AliESDVZERO.h"
53 #include "AliAODVZERO.h"
54 #include "AliPID.h"
55 #include "AliPIDResponse.h"
56 #include "AliAODMCParticle.h"
57 #include "AliAnalysisTaskVnV0.h"
58 #include "AliEventPoolManager.h"
59
60 class AliFlowTrackCuts;
61
62 using std::cout;
63 using std::endl;
64
65 ClassImp(AliAnalysisTaskPhiFlow)
66 ClassImp(AliPhiMesonHelperTrack)
67
68 AliAnalysisTaskPhiFlow::AliAnalysisTaskPhiFlow() : AliAnalysisTaskSE(),
69    fDebug(0), fIsMC(0), fEventMixing(0), fTypeMixing(0), fQA(0), fV0(0), fAODAnalysis(0), fMassBins(1), fMinMass(-1.), fMaxMass(0.), fCutsRP(NULL), fNullCuts(0), fPIDResponse(0), fFlowEvent(0), fBayesianResponse(0), fCandidates(0), fOldTrackParam(0), fRequireTPCStandAlone(0), fCandidateEtaPtCut(0), fCandidateMinEta(0), fCandidateMaxEta(0), fCandidateMinPt(0), fCandidateMaxPt(0), fNPtBins(18), fCentrality(999), fVertex(999), fESD(0), fAOD(0), fPoolManager(0), fOutputList(0), fEventStats(0), fCentralityPass(0), fCentralityNoPass(0), fNOPID(0), fPIDk(0),fNOPIDTOF(0), fPIDTOF(0), fPtP(0), fPtN(0), fPtKP(0), fPtKN(0), fMultCorAfterCuts(0), fMultvsCentr(0), fCentralityMin(0), fCentralityMax(100), fkCentralityMethodA(0), fkCentralityMethodB(0), fCentralityCut2010(0), fCentralityCut2011(0), fPOICuts(0), fVertexRange(0), fPhi(0), fPt(0), fEta(0), fVZEROA(0), fVZEROC(0), fTPCM(0), fDeltaDipAngle(0), fDeltaDipPt(0), fApplyDeltaDipCut(0), fDCAAll(0), fDCAXYQA(0), fDCAZQA(0), fDCAPrim(0), fDCASecondaryWeak(0), fDCAMaterial(0), fSubEventDPhiv2(0)
70 {
71    // Default constructor
72    for(Int_t i(0); i < 7; i++) fPIDConfig[i] = 0.;
73    for(Int_t i(0); i < 5; i++) fDCAConfig[i] = 0.;
74    for(Int_t i(0); i < 20; i++) {
75        fVertexMixingBins[i] = 0;
76        fCentralityMixingBins[i] = 0;
77    }
78    fMixingParameters[0] = 1000; fMixingParameters[1] = 50000; fMixingParameters[2] = 5;
79    for(Int_t i(0); i < 18; i++) {
80        for(Int_t j(0); j < 2; j++) fV0Data[i][j] = 0;
81        fInvMNP[i] = 0; fInvMNN[i] = 0; fInvMPP[i] = 0; fPtSpectra[i] = 0; fPtBins[i] = 0.;
82    }
83 }
84 //_____________________________________________________________________________
85 AliAnalysisTaskPhiFlow::AliAnalysisTaskPhiFlow(const char *name) : AliAnalysisTaskSE(name),
86    fDebug(0), fIsMC(0), fEventMixing(0), fTypeMixing(0), fQA(0), fV0(0), fAODAnalysis(0), fMassBins(1), fMinMass(-1.), fMaxMass(0.), fCutsRP(NULL), fNullCuts(0), fPIDResponse(0), fFlowEvent(0), fBayesianResponse(0), fCandidates(0), fOldTrackParam(0), fRequireTPCStandAlone(0), fCandidateEtaPtCut(0), fCandidateMinEta(0), fCandidateMaxEta(0), fCandidateMinPt(0), fCandidateMaxPt(0), fNPtBins(18), fCentrality(999), fVertex(999), fESD(0), fAOD(0), fPoolManager(0), fOutputList(0), fEventStats(0), fCentralityPass(0), fCentralityNoPass(0), fNOPID(0), fPIDk(0), fNOPIDTOF(0), fPIDTOF(0), fPtP(0), fPtN(0), fPtKP(0), fPtKN(0), fMultCorAfterCuts(0), fMultvsCentr(0), fCentralityMin(0), fCentralityMax(100), fkCentralityMethodA(0), fkCentralityMethodB(0), fCentralityCut2010(0), fCentralityCut2011(0), fPOICuts(0), fVertexRange(0), fPhi(0), fPt(0), fEta(0), fVZEROA(0), fVZEROC(0), fTPCM(0), fDeltaDipAngle(0), fDeltaDipPt(0), fApplyDeltaDipCut(0), fDCAAll(0), fDCAXYQA(0), fDCAZQA(0), fDCAPrim(0), fDCASecondaryWeak(0), fDCAMaterial(0), fSubEventDPhiv2(0)
87 {
88    // Constructor
89    for(Int_t i(0); i < 7; i++) fPIDConfig[i] = 0.;
90    for(Int_t i(0); i < 5; i++) fDCAConfig[i] = 0.;
91    for(Int_t i(0); i < 20; i++) {
92        fVertexMixingBins[i] = 0;
93        fCentralityMixingBins[i] = 0;
94    }
95    fMixingParameters[0] = 1000; fMixingParameters[1] = 50000; fMixingParameters[2] = 5;
96    for(Int_t i(0); i < 18; i++) {
97        for(Int_t j(0); j < 2; j++) fV0Data[i][j] = 0;
98        fInvMNP[i] = 0; fInvMNN[i] = 0; fInvMPP[i] = 0; fPtSpectra[i] = 0; fPtBins[i] = 0.;
99    }
100    DefineInput(0, TChain::Class());
101    DefineOutput(1, TList::Class());
102    DefineOutput(2, AliFlowEventSimple::Class());
103    if(fDebug > 0) cout << " === Created instance of AliAnaysisTaskPhiFlow === " << endl;
104 }
105 //_____________________________________________________________________________
106 AliAnalysisTaskPhiFlow::~AliAnalysisTaskPhiFlow()
107 {
108    // Destructor
109    if (fNullCuts) delete fNullCuts;
110    if (fOutputList) delete fOutputList;
111    if (fCandidates) delete fCandidates;
112    if (fFlowEvent) delete fFlowEvent;
113    if (fBayesianResponse) delete fBayesianResponse;
114    if (fEventMixing) delete fPoolManager;
115    if (fDebug > 0) cout << " === Deleted instance of AliAnalysisTaskPhiFlow === " << endl;
116 }
117 //_____________________________________________________________________________
118 TH1F* AliAnalysisTaskPhiFlow::BookHistogram(const char* name)
119 {
120    // Return a pointer to a TH1 with predefined binning
121    if(fDebug > 0) cout << " *** BookHistogram() *** " << name << endl;
122    TH1F *hist = new TH1F(name, Form("M_{INV} (%s)", name), 60, .99, 1.092);
123    hist->GetXaxis()->SetTitle("M_{INV} (GeV / c^{2})");
124    hist->GetYaxis()->SetTitle("No. of pairs");
125    hist->SetMarkerStyle(kFullCircle);
126    hist->Sumw2();
127    fOutputList->Add(hist);
128    return hist;
129 }
130 //_____________________________________________________________________________
131 TH2F* AliAnalysisTaskPhiFlow::BookPIDHistogram(const char* name, Bool_t TPC)
132 {
133    // Return a pointer to a TH2 with predefined binning
134    if(fDebug > 0) cout << " *** BookPIDHisotgram() *** " << endl;
135    TH2F *hist = 0x0;
136    if(TPC) {
137        hist = new TH2F(name, Form("PID (%s)", name), 100, 0, 5, 100, 0, 1000);
138        hist->GetYaxis()->SetTitle("dE/dx (a.u.)");
139    }
140    if(!TPC) {
141        hist = new TH2F(name, Form("PID (%s)", name), 100, 0, 5, 100, 0, 1.5);
142        hist->GetYaxis()->SetTitle("#beta");
143    }
144    hist->GetXaxis()->SetTitle("P (GeV / c)");
145    fOutputList->Add(hist);
146    return hist;
147 }
148 //_____________________________________________________________________________
149 TH1F* AliAnalysisTaskPhiFlow::InitPtSpectraHistograms(Float_t nmin, Float_t nmax)
150 {
151    // intialize p_t histograms for each p_t bin
152    if(fDebug > 0) cout << " *** InitPtSpectraHistograms() *** " << endl;
153    TH1F* hist = new TH1F(Form("%4.2f p_{t} %4.2f", nmin, nmax), Form("%f p_{t} %f", nmin, nmax), 60, nmin, nmax);
154    hist->GetXaxis()->SetTitle("p_{T} GeV / c");
155    fOutputList->Add(hist);
156    return hist;
157 }
158 //_____________________________________________________________________________
159 TH1F* AliAnalysisTaskPhiFlow::BookPtHistogram(const char* name)
160 {
161    // Return a pointer to a p_T spectrum histogram
162    if(fDebug > 0) cout << " *** BookPtHistogram() *** " << endl;
163    TH1F* ratio = new TH1F(name, name, 100, 0, 7);
164    ratio->GetXaxis()->SetTitle("p_{T} ( GeV / c^{2} )");
165    ratio->GetYaxis()->SetTitle("No. of events");
166    ratio->Sumw2();
167    fOutputList->Add(ratio);
168    return ratio;
169 }
170 //_____________________________________________________________________________
171 void AliAnalysisTaskPhiFlow::AddPhiIdentificationOutputObjects()
172 {
173    // Add Phi Identification Output Objects
174    if(fDebug > 0) cout << " ** AddPhiIdentificationOutputObjects() *** " << endl;
175    if(fQA) {
176        fEventStats = new TH1F("fHistStats", "Event Statistics", 18, -.25, 4.25);
177        fEventStats->GetXaxis()->SetTitle("No. of events");
178        fEventStats->GetYaxis()->SetTitle("Statistics");
179        fOutputList->Add(fEventStats);
180    }
181    fCentralityPass = new TH1F("fCentralityPass", "Centrality Pass", 101, -1, 100);
182    fOutputList->Add(fCentralityPass);
183    if(fQA) {
184        fCentralityNoPass = new TH1F("fCentralityNoPass", "Centrality No Pass", 101, -1, 100);
185        fOutputList->Add(fCentralityNoPass);
186        fNOPID = BookPIDHistogram("TPC signal, all particles", kTRUE);
187        fPIDk = BookPIDHistogram("TPC signal, kaons", kTRUE);
188        fNOPIDTOF = BookPIDHistogram("TOF signal, all particles", kFALSE);
189        fPIDTOF = BookPIDHistogram("TOF signal, kaons", kFALSE);
190    }
191    for(Int_t ptbin(0); ptbin < fNPtBins; ptbin++) {
192        fInvMNP[ptbin] = BookHistogram(Form("NP, %4.2f < p_{T} < %4.2f GeV", fPtBins[ptbin], fPtBins[ptbin+1]));
193        fInvMNN[ptbin] = BookHistogram(Form("NN, %4.2f < p_{T} < %4.2f GeV", fPtBins[ptbin], fPtBins[ptbin+1]));
194        fInvMPP[ptbin] = BookHistogram(Form("PP, %4.2f < p_{T} < %4.2f GeV", fPtBins[ptbin], fPtBins[ptbin+1]));
195    }
196    if(fQA) {
197        for(Int_t ptbin(0); ptbin < fNPtBins; ptbin++) fPtSpectra[ptbin] = InitPtSpectraHistograms(fPtBins[ptbin], fPtBins[ptbin+1]);
198        fPtP = BookPtHistogram("i^{+}");
199        fPtN = BookPtHistogram("i^{-}");
200        fPtKP = BookPtHistogram("K^{+}");
201        fPtKN = BookPtHistogram("K^{-}");
202        fPhi = new TH1F("fPhi", "#phi distribution", 100, -.5, 7);
203        fOutputList->Add(fPhi);
204        fPt = new TH1F("fPt", "p_{T}", 100, 0, 5.5);
205        fOutputList->Add(fPt);
206        fEta = new TH1F("fEta", "#eta distribution", 100, -1.1, 1.1);
207        fOutputList->Add(fEta);
208        fVZEROA = new TH1F("fVZEROA", "VZERO A Multiplicity", 1000, 0, 10000);
209        fOutputList->Add(fVZEROA);
210        fVZEROC = new TH1F("fVZEROC", "VZERO C Multiplicity", 1000, 0, 10000);
211        fOutputList->Add(fVZEROC);
212        fTPCM = new TH1F("fTPCM", "TPC multiplicity", 1000, 0, 10000);
213        fOutputList->Add(fTPCM);
214        fDCAXYQA = new TH1F("fDCAXYQA", "fDCAXYQA", 1000, -5, 5);
215        fOutputList->Add(fDCAXYQA);
216        fDCAZQA = new TH1F("fDCAZQA", "fDCAZQA", 1000, -5, 5);
217        fOutputList->Add(fDCAZQA);
218        if(fCentralityCut2010 || fCentralityCut2011) {
219            fMultCorAfterCuts = new TH2F("fMultCorAfterCuts", "TPC vs Global multiplicity (After cuts); Global multiplicity; TPC multiplicity", 100, 0, 3000, 100, 0, 3000);
220            fOutputList->Add(fMultCorAfterCuts);
221            fMultvsCentr = new TH2F("fMultvsCentr", "Multiplicity vs centrality; centrality; Multiplicity", 9, -0.5, 100.5, 101, 0, 3000);
222            fOutputList->Add(fMultvsCentr);
223        }
224    }
225    if(fIsMC || fQA) {
226        fDCAAll = new TH2F("fDCAAll", "fDCAAll", 1000, 0, 10, 1000, -5, 5);
227        fOutputList->Add(fDCAAll);
228        fDCAPrim = new TH2F("fDCAprim","fDCAprim", 1000, 0, 10, 1000, -5, 5);
229        fOutputList->Add(fDCAPrim);
230        fDCASecondaryWeak = new TH2F("fDCASecondaryWeak","fDCASecondaryWeak", 1000, 0, 10, 1000, -5, 5);
231        fOutputList->Add(fDCASecondaryWeak);
232        fDCAMaterial = new TH2F("fDCAMaterial","fDCAMaterial", 1000, 0, 10, 1000, -5, 5);
233        fOutputList->Add(fDCAMaterial);
234    }
235    if(fV0) {
236        fSubEventDPhiv2 = new TProfile("fSubEventDPhiv2", "fSubEventDPhiv2", 5, 0, 5);
237        fSubEventDPhiv2->GetXaxis()->SetBinLabel(1, "<#Psi_{a} - #Psi_{b}>");
238        fSubEventDPhiv2->GetXaxis()->SetBinLabel(2, "<#Psi_{a} - #Psi_{c}>");
239        fSubEventDPhiv2->GetXaxis()->SetBinLabel(3, "<#Psi_{b} - #Psi_{c}>");
240        fSubEventDPhiv2->GetXaxis()->SetBinLabel(4, "#sqrt{#frac{<#Psi_{a} - #Psi_{b}><#Psi_{a} - #Psi_{c}>}{<#Psi_{b} - #Psi_{c}>}}");
241        fSubEventDPhiv2->GetXaxis()->SetBinLabel(5, "#sqrt{#frac{<#Psi_{a} - #Psi_{c}><#Psi_{b} - #Psi_{c}>}{<#Psi_{a} - #Psi_{b}>}}");
242        fOutputList->Add(fSubEventDPhiv2);
243        const char* V0[] = {"V0A", "V0C"};
244        for(Int_t ptbin(0); ptbin < fNPtBins; ptbin++)
245            for(Int_t iV0(0); iV0 < 2; iV0++) {
246                    fV0Data[ptbin][iV0] = new TProfile(Form("%s v2 %4.2f < p_{T} < %4.2f GeV", V0[iV0], fPtBins[ptbin], fPtBins[ptbin+1]), Form("%s v2 %4.2f < p_{T} < %4.2f GeV", V0[iV0], fPtBins[ptbin], fPtBins[ptbin+1]), fMassBins, fMinMass, fMaxMass);
247                    fOutputList->Add(fV0Data[ptbin][iV0]);
248            }
249    }
250 }
251 //_____________________________________________________________________________
252 void AliAnalysisTaskPhiFlow::UserCreateOutputObjects()
253 {
254    // Create user defined output objects
255    if(fDebug > 0) cout << " *** UserCreateOutputObjects() *** " << endl;
256    fNullCuts = new AliFlowTrackCuts("null_cuts");
257    fBayesianResponse = new AliFlowBayesianPID();
258    Double_t t(0);
259    for(Int_t i = 0; i < 7; i++) t+=TMath::Abs(fPIDConfig[i]);
260    if(t < 0.1) AliFatal("No valid PID procedure recognized -- terminating analysis !!!");
261    if(fNPtBins > 18) AliFatal("Invalid number of pt bins initialied ( > 18 ) -- terminating analysis !!!");
262    AliFlowCommonConstants* cc = AliFlowCommonConstants::GetMaster();
263    cc->SetNbinsQ(500);           cc->SetNbinsPhi(180);           cc->SetNbinsMult(10000);
264    cc->SetQMin(0.0);             cc->SetPhiMin(0.0);             cc->SetMultMin(0);
265    cc->SetQMax(3.0);             cc->SetPhiMax(TMath::TwoPi());  cc->SetMultMax(10000);
266    cc->SetNbinsMass(fMassBins);  cc->SetNbinsEta(200);           (fMassBins == 1) ? cc->SetNbinsPt(15) : cc->SetNbinsPt(100); // high pt
267    cc->SetMassMin(fMinMass);     cc->SetEtaMin(-5.0);            cc->SetPtMin(0);
268    cc->SetMassMax(fMaxMass);     cc->SetEtaMax(+5.0);            (fMassBins == 1) ? cc->SetPtMax(15) : cc->SetPtMax(10); // high pt
269    if (!fOldTrackParam) fBayesianResponse->SetNewTrackParam();
270    AliAnalysisManager *man = AliAnalysisManager::GetAnalysisManager();
271    if (man) {
272       AliInputEventHandler* inputHandler = (AliInputEventHandler*)(man->GetInputEventHandler());
273       if (inputHandler)   fPIDResponse = inputHandler->GetPIDResponse();
274    }
275    // Create all output objects and store them to a list
276    fOutputList = new TList();
277    fOutputList->SetOwner(kTRUE);
278    // Create and post the Phi identification output objects
279    AddPhiIdentificationOutputObjects();
280    PostData(1, fOutputList);
281    // create candidate array
282    fCandidates = new TObjArray(1000);
283    fCandidates->SetOwner(kTRUE);
284    // create and post flowevent
285    fFlowEvent = new AliFlowEvent(10000);
286    PostData(2, fFlowEvent);
287    if(fEventMixing) fPoolManager = InitializeEventMixing();
288 }
289 //_____________________________________________________________________________
290 AliEventPoolManager* AliAnalysisTaskPhiFlow::InitializeEventMixing()
291 {
292    // initialize event mixing
293   if(fDebug > 0) cout << " *** InitializeEventMixing() *** " << endl;
294   Int_t _c(0), _v(0);
295   for(Int_t i(0); i < 19; i++) {
296       if (fCentralityMixingBins[i+1] < fCentralityMixingBins[i]) { _c = i; break; }
297       else _c = 19;
298   }
299   for(Int_t i(0); i < 19; i++) {
300       if (fVertexMixingBins[i+1] < fVertexMixingBins[i]) { _v = i; break; }
301       else _v = 19;
302   }
303   if(fDebug > 0 ) cout << Form("   --> found %d centrality bins for mixing, %d vertex bins for mixing", _c, _v) << endl;
304   Double_t centralityBins[_c];
305   Double_t vertexBins[_v];
306   for(Int_t i(0); i < _c + 1; i++) centralityBins[i] = fCentralityMixingBins[i];
307   for(Int_t i(0); i < _v + 1; i++) vertexBins[i] = fVertexMixingBins[i];
308   return new AliEventPoolManager(fMixingParameters[0], fMixingParameters[1], _c, (Double_t*)centralityBins, _v, (Double_t*)vertexBins);
309 }
310 //_____________________________________________________________________________
311 template <typename T> Double_t AliAnalysisTaskPhiFlow::InvariantMass(const T* track1, const T* track2) const
312 {
313    // Return the invariant mass of two tracks, assuming both tracks are kaons
314    if(fDebug > 1) cout << " *** InvariantMass() *** " << endl;
315    if ((!track2) || (!track1)) return 0.;
316    Double_t masss = TMath::Power(4.93676999999999977e-01, 2);
317    Double_t pxs = TMath::Power((track1->Px() + track2->Px()), 2);
318    Double_t pys = TMath::Power((track1->Py() + track2->Py()), 2);
319    Double_t pzs = TMath::Power((track1->Pz() + track2->Pz()), 2);
320    Double_t e1 = TMath::Sqrt(track1->P() * track1->P() + masss);
321    Double_t e2 = TMath::Sqrt(track2->P() * track2->P() + masss);
322    Double_t es = TMath::Power((e1 + e2), 2);
323    if ((es - (pxs + pys + pzs)) < 0) return 0.;
324    return TMath::Sqrt((es - (pxs + pys + pzs)));
325 }
326 //_____________________________________________________________________________
327 template <typename T> Double_t AliAnalysisTaskPhiFlow::DeltaDipAngle(const T* track1, const T* track2) const
328 {
329    // Calculate the delta dip angle between two particles
330    if(fDebug > 1) cout << " *** DeltaDipAngle() *** " << endl;
331    if (track1->P()*track2->P() == 0) return 999;
332    return TMath::ACos(((track1->Pt() * track2->Pt()) + (track1->Pz() * track2->Pz())) / (track1->P() * track2->P()));
333 }
334 //_____________________________________________________________________________
335 template <typename T> Bool_t AliAnalysisTaskPhiFlow::CheckDeltaDipAngle(const T* track1, const T* track2) const
336 {
337    // Check if pair passes delta dip angle cut within 0 < p_t < fDeltaDipPt
338    if(fDebug > 1) cout << " *** CheckDeltaDipAngle() *** " << endl;
339    if ((TMath::Abs(DeltaDipAngle(track1, track2)) < fDeltaDipAngle) && (PhiPt(track1, track2) < fDeltaDipPt)) return kFALSE;
340    return kTRUE;
341 }
342 //_____________________________________________________________________________
343 template <typename T> Bool_t AliAnalysisTaskPhiFlow::CheckCandidateEtaPtCut(const T* track1, const T* track2) const
344 {
345    // Check if pair passes eta and pt cut
346    if(fDebug > 1) cout << " *** CheckCandidateEtaPtCut() *** " << endl;
347    if (fCandidateMinPt > PhiPt(track1, track2) || fCandidateMaxPt < PhiPt(track1, track2)) return kFALSE;
348    TVector3 a(track1->Px(), track1->Py(), track1->Pz());
349    TVector3 b(track2->Px(), track2->Py(), track2->Pz());
350    TVector3 c = a + b;
351    if (fCandidateMinEta > c.Eta() || fCandidateMaxEta < c.Eta()) return kFALSE;
352    return kTRUE;
353 }
354 //_____________________________________________________________________________
355 template <typename T> Bool_t AliAnalysisTaskPhiFlow::EventCut(T* event)
356 {
357    // Impose event cuts
358    if(fDebug > 0) cout << " *** EventCut() *** " << endl;
359    if (!event) return kFALSE;
360    if (!CheckVertex(event)) return kFALSE;
361    if (!CheckCentrality(event)) return kFALSE;
362    if(fQA) PlotMultiplcities(event);
363    return kTRUE;
364 }
365 //_____________________________________________________________________________
366 template <typename T> void AliAnalysisTaskPhiFlow::PlotMultiplcities(const T* event) const
367 {
368    // QA multiplicity plots
369    if(fDebug > 1) cout << " *** PlotMultiplcities() *** " << endl;
370    fVZEROA->Fill(event->GetVZEROData()->GetMTotV0A());
371    fVZEROC->Fill(event->GetVZEROData()->GetMTotV0C());
372    fTPCM->Fill(event->GetNumberOfTracks());
373 }
374 //_____________________________________________________________________________
375 template <typename T> Bool_t AliAnalysisTaskPhiFlow::CheckVertex(const T* event)
376 {
377    // Check if event vertex is within given range
378    if(fDebug > 0) cout << " *** CheckVertex() *** " << endl;
379    if (!event->GetPrimaryVertex()) return 0x0;
380    fVertex = event->GetPrimaryVertex()->GetZ();
381    if (TMath::Abs(fVertex) > fVertexRange) return 0x0;
382    return kTRUE;
383 }
384 //_____________________________________________________________________________
385 template <typename T> Bool_t AliAnalysisTaskPhiFlow::CheckCentrality(T* event)
386 {
387    // Check if event is within the set centrality range. Falls back to V0 centrality determination if no method is set
388    if(fDebug > 0) cout << " *** CheckCentrality() *** " << endl;
389    if (!fkCentralityMethodA) AliFatal("No centrality method set! FATAL ERROR!");
390    fCentrality = event->GetCentrality()->GetCentralityPercentile(fkCentralityMethodA);
391    Double_t cenB(-999);
392    // if a second centrality estimator is requited, set it
393    (fkCentralityMethodB) ? cenB = event->GetCentrality()->GetCentralityPercentile(fkCentralityMethodB) : cenB = fCentrality;
394    if (TMath::Abs(fCentrality-cenB) > 5 || cenB >= 80 || cenB < 0 || fCentrality <= fCentralityMin || fCentrality > fCentralityMax) {
395       if(fQA) fCentralityNoPass->Fill(fCentrality) ;
396       return kFALSE;
397    }
398    const Int_t nGoodTracks = event->GetNumberOfTracks();
399    if(fCentralityCut2010) { // cut on outliers
400       Float_t multTPC(0.); // tpc mult estimate
401       Float_t multGlob(0.); // global multiplicity
402       for(Int_t iTracks = 0; iTracks < nGoodTracks; iTracks++) { // fill tpc mult
403           AliAODTrack* trackAOD = event->GetTrack(iTracks);
404           if (!trackAOD) continue;
405           if (!(trackAOD->TestFilterBit(1))) continue;
406           if ((trackAOD->Pt() < .2) || (trackAOD->Pt() > 5.0) || (TMath::Abs(trackAOD->Eta()) > .8) || (trackAOD->GetTPCNcls() < 70)  || (trackAOD->GetDetPid()->GetTPCsignal() < 10.0) || (trackAOD->Chi2perNDF() < 0.2)) continue;
407           multTPC++;
408       }
409       for(Int_t iTracks = 0; iTracks < nGoodTracks; iTracks++) { // fill global mult
410           AliAODTrack* trackAOD = event->GetTrack(iTracks);
411           if (!trackAOD) continue;
412           if (!(trackAOD->TestFilterBit(16))) continue;
413           if ((trackAOD->Pt() < .2) || (trackAOD->Pt() > 5.0) || (TMath::Abs(trackAOD->Eta()) > .8) || (trackAOD->GetTPCNcls() < 70) || (trackAOD->GetDetPid()->GetTPCsignal() < 10.0) || (trackAOD->Chi2perNDF() < 0.1)) continue;
414           Double_t b[2] = {-99., -99.};
415           Double_t bCov[3] = {-99., -99., -99.};
416           if (!(trackAOD->PropagateToDCA(event->GetPrimaryVertex(), event->GetMagneticField(), 100., b, bCov))) continue;
417           if ((TMath::Abs(b[0]) > 0.3) || (TMath::Abs(b[1]) > 0.3)) continue;
418           multGlob++;
419       } //track loop
420  //     printf(" mult TPC %.2f, mult Glob %.2f \n", multTPC, multGlob);
421       if(! (multTPC > (-40.3+1.22*multGlob) && multTPC < (32.1+1.59*multGlob))) return kFALSE;
422       if(fQA) {  
423           fMultCorAfterCuts->Fill(multGlob, multTPC);  
424           fMultvsCentr->Fill(fCentrality, multTPC);
425       }
426    }
427
428    if(fCentralityCut2011) { // cut on outliers
429       Float_t multTPC(0.); // tpc mult estimate
430       Float_t multGlob(0.); // global multiplicity
431       for(Int_t iTracks = 0; iTracks < nGoodTracks; iTracks++) { // fill tpc mult
432           AliAODTrack* trackAOD = event->GetTrack(iTracks);
433           if (!trackAOD) continue;
434           if (!(trackAOD->TestFilterBit(1))) continue;
435           if ((trackAOD->Pt() < .2) || (trackAOD->Pt() > 5.0) || (TMath::Abs(trackAOD->Eta()) > .8) || (trackAOD->GetTPCNcls() < 70)  || (trackAOD->GetDetPid()->GetTPCsignal() < 10.0) || (trackAOD->Chi2perNDF() < 0.2)) continue;
436           multTPC++;
437       }
438       for(Int_t iTracks = 0; iTracks < nGoodTracks; iTracks++) { // fill global mult
439           AliAODTrack* trackAOD = event->GetTrack(iTracks);
440           if (!trackAOD) continue;
441           if (!(trackAOD->TestFilterBit(16))) continue;
442           if ((trackAOD->Pt() < .2) || (trackAOD->Pt() > 5.0) || (TMath::Abs(trackAOD->Eta()) > .8) || (trackAOD->GetTPCNcls() < 70) || (trackAOD->GetDetPid()->GetTPCsignal() < 10.0) || (trackAOD->Chi2perNDF() < 0.1)) continue;
443           Double_t b[2] = {-99., -99.};
444           Double_t bCov[3] = {-99., -99., -99.};
445           if (!(trackAOD->PropagateToDCA(event->GetPrimaryVertex(), event->GetMagneticField(), 100., b, bCov))) continue;
446           if ((TMath::Abs(b[0]) > 0.3) || (TMath::Abs(b[1]) > 0.3)) continue;
447           multGlob++;
448       } //track loop
449       //printf(" mult TPC %.2f, mult Glob %.2f \n", multTPC, multGlob);
450       if(! (multTPC > (-36.73 + 1.48*multGlob) && multTPC < (62.87 + 1.78*multGlob))) return kFALSE;
451       if(fQA) {  
452           fMultCorAfterCuts->Fill(multGlob, multTPC);  
453           fMultvsCentr->Fill(fCentrality, multTPC);
454       }
455    }
456
457    fCentralityPass->Fill(fCentrality);
458    return kTRUE;
459 }
460 //_____________________________________________________________________________
461 void AliAnalysisTaskPhiFlow::InitializeBayesianPID(AliAODEvent* event)
462 {
463    // Initialize the Bayesian PID object for AOD
464    if(fDebug > 0) cout << " *** InitializeBayesianPID() *** " << endl;
465    fBayesianResponse->SetDetResponse(event, fCentrality);
466 }
467 //_____________________________________________________________________________
468 template <typename T> Bool_t AliAnalysisTaskPhiFlow::PassesTPCbayesianCut(T* track) const
469 {
470    // Check if the particle passes the TPC TOF bayesian cut.
471    if(fDebug > 1) cout << " *** PassesTPCbayesianCut() *** " << endl;
472    fBayesianResponse->ComputeProb(track);
473    if (!fBayesianResponse->GetCurrentMask(0)) return kFALSE; // return false if TPC has no response
474    Float_t *probabilities = fBayesianResponse->GetProb();
475    if (probabilities[3] > fPIDConfig[6]) {
476       if(fQA) {fPhi->Fill(track->Phi()); fPt->Fill(track->Pt()); fEta->Fill(track->Eta());}
477       return kTRUE;
478    }
479    return kFALSE;
480 }
481 //_____________________________________________________________________________
482 Bool_t AliAnalysisTaskPhiFlow::PassesDCACut(AliAODTrack* track) const
483 {
484     // check if track passes dca cut according to dca routine
485     // setup the routine as follows:
486     // fDCAConfig[0] < -1 no pt dependence
487     // fDCAConfig[0] =  0 do nothing
488     // fDCAConfig[0] >  1 pt dependent dca cut
489     if(fDebug > 1) cout << " *** PassesDCACut() *** " << endl;
490     if(fIsMC) return kTRUE;
491     if( (fDCAConfig[0] < 0.1) && (fDCAConfig[0] > -0.1) ) {
492         if(fQA) {
493             Double_t b[2] = { -99., -99.};
494             Double_t bCov[3] = { -99., -99., -99.};
495             track->PropagateToDCA(fAOD->GetPrimaryVertex(), fAOD->GetMagneticField(), 100., b, bCov);
496             fDCAXYQA->Fill(b[0]);
497             fDCAZQA->Fill(b[1]);
498             fDCAPrim->Fill(track->Pt(), b[0]);
499         }
500         return kTRUE;
501     }
502     Double_t b[2] = { -99., -99.};
503     Double_t bCov[3] = { -99., -99., -99.};
504     track->PropagateToDCA(fAOD->GetPrimaryVertex(), fAOD->GetMagneticField(), 100., b, bCov);
505     if((!fIsMC)&&fQA) fDCAMaterial->Fill(track->Pt(), b[0]);
506     if( (fDCAConfig[0] < -.9) && ( (TMath::Abs(b[0]) > fDCAConfig[1]) || (TMath::Abs(b[1]) > fDCAConfig[2])) ) return kFALSE;
507     if(fDCAConfig[0] > .9) {
508         if(fDCAConfig[4] < TMath::Abs(b[1])) return kFALSE;
509         Double_t denom = TMath::Power(track->Pt(), fDCAConfig[3]);
510         if( denom < 0.0000001 ) return kFALSE; // avoid division by zero
511         if( (fDCAConfig[1] + fDCAConfig[2] / denom) < TMath::Abs(b[0]) ) return kFALSE;
512     }
513     if(fQA) {
514         fDCAXYQA->Fill(b[0]);
515         fDCAZQA->Fill(b[1]);
516         fDCAPrim->Fill(track->Pt(), b[0]);
517     }
518     return kTRUE;
519 }
520 //_____________________________________________________________________________
521 Bool_t AliAnalysisTaskPhiFlow::IsKaon(AliAODTrack* track) const
522 {
523    // Kaon identification routine, based on multiple detectors and approaches
524    if(fDebug > 1) cout << " *** IsKaon() *** " << endl;
525    if(fRequireTPCStandAlone && (!track->TestFilterBit(1))) return kFALSE;
526    if(!PassesDCACut(track)) return kFALSE;
527    if(fQA) {fNOPID->Fill(track->P(), track->GetTPCsignal());fNOPIDTOF->Fill(track->P(), track->GetTOFsignal());}
528    if(track->Pt() < fPIDConfig[1]) {
529        if(fDebug>1) cout << " ITS received track with p_t " << track->Pt() << endl;
530        // if tpc control is disabled, pure its pid
531        if(fPIDConfig[2] < 0.) {
532            if (TMath::Abs(fPIDResponse->NumberOfSigmasITS(track, AliPID::kKaon)) < fPIDConfig[0]) return kTRUE;
533            return kFALSE;
534        }
535        // else, switch to ITS pid with TPC rejection of protons and pions
536        if (TMath::Abs(fPIDResponse->NumberOfSigmasTPC(track, AliPID::kProton)) < fPIDConfig[3]) return kFALSE;
537        else if (TMath::Abs(fPIDResponse->NumberOfSigmasTPC(track, AliPID::kPion)) < fPIDConfig[3]) return kFALSE;
538        else if (TMath::Abs(fPIDResponse->NumberOfSigmasITS(track, AliPID::kKaon)) < fPIDConfig[0]) {
539            if(fQA) {fPIDk->Fill(track->P(), track->GetTPCsignal()); fPIDTOF->Fill(track->P(), track->GetTOFsignal());}
540            return kTRUE;
541        }
542        return kFALSE;
543    }
544    if((track->Pt() > fPIDConfig[1]) && (track->Pt() < fPIDConfig[4])) {
545        if(fDebug>1) cout << " TPC received track with p_t " << track->Pt() << endl;
546        // if its control is disabled, pure tpc pid
547        if(fPIDConfig[5] < 0.) {
548            if (TMath::Abs(fPIDResponse->NumberOfSigmasTPC(track, AliPID::kKaon)) < fPIDConfig[3]) return kTRUE;
549            return kFALSE;
550        }
551        // else, switch to TPC pid with ITS rejection of protons and pions
552        if (TMath::Abs(fPIDResponse->NumberOfSigmasITS(track, AliPID::kProton)) < fPIDConfig[0]) return kFALSE;
553        else if (TMath::Abs(fPIDResponse->NumberOfSigmasITS(track, AliPID::kPion)) < fPIDConfig[0]) return kFALSE;
554        else if (TMath::Abs(fPIDResponse->NumberOfSigmasTPC(track, AliPID::kKaon)) < fPIDConfig[3]) {
555            if(fQA) {fPIDk->Fill(track->P(), track->GetTPCsignal()); fPIDTOF->Fill(track->P(), track->GetTOFsignal());}
556            return kTRUE;
557        }
558        return kFALSE;
559    }
560    if(fDebug>1) cout << " Bayesian method received track with p_t " << track->Pt() << endl;
561    // switch to bayesian PID
562    if (PassesTPCbayesianCut(track)) {
563        if(fQA) {fPIDk->Fill(track->P(), track->GetTPCsignal());fPIDTOF->Fill(track->P(), track->GetTOFsignal());}
564        return kTRUE;
565    }
566    return kFALSE;
567 }
568 //_____________________________________________________________________________
569 template <typename T> Double_t AliAnalysisTaskPhiFlow::PhiPt(const T* track1, const T* track2) const
570 {
571    // return p_t of track pair
572    TVector3 a(track1->Px(), track1->Py(), track1->Pz());
573    TVector3 b(track2->Px(), track2->Py(), track2->Pz());
574    TVector3 c = a + b;
575    return c.Pt();
576 }
577 //_____________________________________________________________________________
578 template <typename T> void AliAnalysisTaskPhiFlow::PtSelector(Int_t tracktype, const T* track1, const T* track2) const
579 {
580    // plot m_inv spectra of like- and unlike-sign kaon pairs for each pt bin
581    Double_t pt = PhiPt(track1, track2);
582    if (tracktype == 0) {
583        for(Int_t ptbin(0); ptbin < fNPtBins; ptbin++) {
584            if ((fPtBins[ptbin] <= pt) && (fPtBins[ptbin+1] > pt )) {
585                fInvMNP[ptbin]->Fill(InvariantMass(track1, track2));
586                if(fQA) fPtSpectra[ptbin]->Fill(pt);
587            }
588        }
589    }
590    if (tracktype == 1) {
591        for(Int_t ptbin(0); ptbin < fNPtBins; ptbin++) {
592            if ((fPtBins[ptbin] <= pt) && (fPtBins[ptbin+1] > pt )) {
593                fInvMPP[ptbin]->Fill(InvariantMass(track1, track2));
594            }
595        }
596    }
597    if (tracktype == 2) {
598        for(Int_t ptbin(0); ptbin < fNPtBins; ptbin++) {
599            if ((fPtBins[ptbin] <= pt) && (fPtBins[ptbin+1] > pt )) {
600                fInvMNN[ptbin]->Fill(InvariantMass(track1, track2));
601            }
602        }
603    }
604 }
605 //_____________________________________________________________________________
606 template <typename T> Bool_t AliAnalysisTaskPhiFlow::PhiTrack(T* track) const
607 {
608    // check if track meets quality cuts
609    if(!track) return kFALSE;
610    return fPOICuts->IsSelected(track);
611 }
612 //_____________________________________________________________________________
613 template <typename T> void AliAnalysisTaskPhiFlow::SetNullCuts(T* event)
614 {
615    // Set null cuts
616    if (fDebug > 0) cout << " *** SetNullCuts() *** " << fCutsRP << endl;
617    fCutsRP->SetEvent(event, MCEvent());
618    fNullCuts->SetParamType(AliFlowTrackCuts::kGlobal);
619    fNullCuts->SetPtRange(+1, -1); // select nothing QUICK
620    fNullCuts->SetEtaRange(+1, -1); // select nothing VZERO
621    fNullCuts->SetEvent(event, MCEvent());
622 }
623 //_____________________________________________________________________________
624 void AliAnalysisTaskPhiFlow::PrepareFlowEvent(Int_t iMulti)
625 {
626    // Prepare flow events
627    if (fDebug > 0 ) cout << " *** PrepareFlowEvent() *** " << endl;
628    fFlowEvent->ClearFast();
629    fFlowEvent->Fill(fCutsRP, fNullCuts);
630    fFlowEvent->SetReferenceMultiplicity(iMulti);
631    fFlowEvent->DefineDeadZone(0, 0, 0, 0);
632 }
633 //_____________________________________________________________________________
634 void AliAnalysisTaskPhiFlow::VZEROSubEventAnalysis()
635 {
636     // vzero event plane analysis using three subevents
637     if(fDebug > 0) cout << " ** VZEROSubEventAnalysis() *** " << endl;
638     if(!AliAnalysisTaskVnV0::IsPsiComputed()) { // AliAnalysisTaskVnV0 must be added to analysis que before this task !!!
639         if(fDebug > 0 ) cout << "Fatal error: unable to retrieve VZERO task output !!! Exiting ..." << endl;
640         return;
641     }
642     // retrieve data
643     Float_t abcPsi2[] = {AliAnalysisTaskVnV0::GetPsi2V0A(), AliAnalysisTaskVnV0::GetPsi2TPC(), AliAnalysisTaskVnV0::GetPsi2V0C()};
644     for(Int_t i(0); i < 3; i++) if(abcPsi2[i] == 0)  {
645         if(fDebug > 0) cout << " Warning: I've encountered 0 in one of the symmetry panes (TPC, VZERO) - skipping event !" << endl;
646             return;
647     }
648     // save info for resolution calculations
649     Float_t qaqb = TMath::Cos(2.*(abcPsi2[0]-abcPsi2[1])); // vzeroa - tpc
650     Float_t qaqc = TMath::Cos(2.*(abcPsi2[0]-abcPsi2[2])); // vzeroa - vzeroc
651     Float_t qbqc = TMath::Cos(2.*(abcPsi2[1]-abcPsi2[2])); // tpc - vzeroc
652     fSubEventDPhiv2->Fill(0.5, qaqb);
653     fSubEventDPhiv2->Fill(1.5, qaqc);
654     fSubEventDPhiv2->Fill(2.5, qbqc);
655     Float_t minv[31];
656     Float_t _dphi[30][fNPtBins][2]; //minv, pt, v0a-c
657     Int_t occurence[30][fNPtBins]; //minv, pt
658     for(Int_t mb(0); mb < 31; mb++) minv[mb] = 0.99 + mb * 0.0034;
659     for(Int_t i(0); i < 30; i++) for (Int_t j(0); j < fNPtBins; j++) for(Int_t k(0); k < 2; k ++) {
660         _dphi[i][j][k] = 0;
661         if(k==0) occurence[i][j] = 0;
662     }
663     for(Int_t iCand(0); iCand < fCandidates->GetEntriesFast(); iCand++) { // loop over unlike sign kaon pairs
664         AliFlowTrackSimple *track = dynamic_cast<AliFlowTrackSimple*>(fCandidates->At(iCand));
665         if(!track) {
666             if(fDebug > 1) cout << " --> dynamic_cast returned null-pointer, skipping candidate <-- " << endl;
667             continue;
668         }
669         for(Int_t mb(0); mb < 30; mb++) { // loop over massbands
670             if(track->Mass() <= minv[mb] || track->Mass() > minv[mb+1]) continue;
671             for(Int_t ptbin(0); ptbin < fNPtBins; ptbin++) { // loop over pt bins
672                 if(track->Pt() <= fPtBins[ptbin] || track->Pt() > fPtBins[ptbin+1]) continue;
673                 _dphi[mb][ptbin][0]+=TMath::Cos(2.*(track->Phi() - abcPsi2[0]));
674                 _dphi[mb][ptbin][1]+=TMath::Cos(2.*(track->Phi() - abcPsi2[2]));
675                 occurence[mb][ptbin]+=1;
676             }
677         }
678         for(Int_t mb(0); mb < 30; mb++) // write vn values to tprofiles
679             for(Int_t ptbin(0); ptbin < fNPtBins; ptbin++) {
680                 if(occurence[mb][ptbin]==0) continue;
681                 fV0Data[ptbin][0]->Fill(mb*0.0034+0.99+0.0017, _dphi[mb][ptbin][0]/(Float_t)occurence[mb][ptbin]);
682                 fV0Data[ptbin][1]->Fill(mb*0.0034+0.99+0.0017, _dphi[mb][ptbin][1]/(Float_t)occurence[mb][ptbin]);
683             }
684     }
685 }
686 //_____________________________________________________________________________
687 void AliAnalysisTaskPhiFlow::UserExec(Option_t *)
688 {
689    // UserExec: called for each event. Commented where necessary
690    if(fDebug > 0 ) cout << " *** UserExec() *** " << endl;
691    TObjArray* MixingCandidates = 0x0; // init null pointer for event mixing
692    if(fEventMixing) {
693        MixingCandidates = new TObjArray();
694        MixingCandidates->SetOwner(kTRUE); // mixing candidates has ownership of objects in array
695    }
696    if (!fPIDResponse) {
697       if(fDebug > 0 ) cout << " Could not get PID response " << endl;
698       return;
699    }
700    fAOD = dynamic_cast<AliAODEvent*>(InputEvent()); // check for aod data type
701    if (fAOD) {
702       fAODAnalysis = kTRUE;
703       if (!EventCut(fAOD)) return; // check for event cuts
704       InitializeBayesianPID(fAOD); // init event objects
705       SetNullCuts(fAOD);
706       PrepareFlowEvent(fAOD->GetNumberOfTracks());
707       fCandidates->SetLast(-1);
708       if(fIsMC) IsMC(); // launch mc stuff
709       if(fQA) fEventStats->Fill(0);
710       Int_t unTracks = fAOD->GetNumberOfTracks();
711       AliAODTrack* un[unTracks];
712       AliAODTrack* up[unTracks];
713       Int_t unp(0);
714       Int_t unn(0);
715       for (Int_t iTracks = 0; iTracks < unTracks; iTracks++) { // select analysis candidates
716          AliAODTrack* track = fAOD->GetTrack(iTracks);
717          if (!PhiTrack(track)) continue;
718          if (fQA) {
719             if(track->Charge() > 0) {fEventStats->Fill(1); fPtP->Fill(track->Pt());}
720             if(track->Charge() < 0) {fEventStats->Fill(2); fPtN->Fill(track->Pt());}
721          }
722          if (IsKaon(track)) {
723             if(fEventMixing) MixingCandidates->Add(new AliPhiMesonHelperTrack(track->Eta(), track->Phi(), track->P(),
724                                                                               track->Px(), track->Py(), track->Pz(),
725                                                                               track->Pt(), track->Charge()));
726             if (track->Charge() > 0) {
727                up[unp] = track;
728                unp++;
729                if(fQA) {fEventStats->Fill(3);fPtKP->Fill(track->Pt());}
730             }
731             else if (track->Charge() < 0) {
732                un[unn] = track;
733                unn++;
734                if(fQA) {fEventStats->Fill(4); fPtKN->Fill(track->Pt());}
735             }
736          }
737       }
738       for (Int_t pTracks = 0; pTracks < unp ; pTracks++) { // perform analysis
739          for (Int_t nTracks = 0; nTracks < unn ; nTracks++) {
740             if (fApplyDeltaDipCut && (!CheckDeltaDipAngle(up[pTracks], un[nTracks]))) continue;
741             if (fCandidateEtaPtCut && (!CheckCandidateEtaPtCut(up[pTracks], un[nTracks]))) continue;
742             PtSelector(0, up[pTracks], un[nTracks]);
743             Double_t pt = PhiPt(up[pTracks], un[nTracks]);
744             Double_t mass = InvariantMass(up[pTracks], un[nTracks]);
745             TVector3 a(up[pTracks]->Px(), up[pTracks]->Py(), up[pTracks]->Pz());
746             TVector3 b(un[nTracks]->Px(), un[nTracks]->Py(), up[pTracks]->Pz());
747             TVector3 c = a + b;
748             Double_t phi = c.Phi();
749             Double_t eta = c.Eta();
750             Int_t nIDs[2];
751             nIDs[0] = up[pTracks]->GetID();
752             nIDs[1] = un[nTracks]->GetID();
753             MakeTrack(mass, pt, phi, eta, 2, nIDs);
754          }
755       }
756       if(fV0) VZEROSubEventAnalysis();
757       if (fDebug > 0)  printf("I received %d candidates\n", fCandidates->GetEntriesFast()); // insert candidates into flow events
758       for (int iCand = 0; iCand != fCandidates->GetEntriesFast(); ++iCand) {
759          AliFlowCandidateTrack *cand = dynamic_cast<AliFlowCandidateTrack*>(fCandidates->At(iCand));
760          if (!cand) continue;
761          if (fDebug > 1) printf(" --> Checking at candidate %d with %d daughters: mass %f\n", iCand, cand->GetNDaughters(), cand->Mass());
762          for (int iDau = 0; iDau != cand->GetNDaughters(); ++iDau) {
763             if (fDebug>1) printf("     *** Daughter %d with fID %d ***", iDau, cand->GetIDDaughter(iDau));
764             for (int iRPs = 0; iRPs != fFlowEvent->NumberOfTracks(); ++iRPs) {
765                AliFlowTrack *iRP = dynamic_cast<AliFlowTrack*>(fFlowEvent->GetTrack(iRPs));
766                if (!iRP) continue;
767                if (!iRP->InRPSelection()) continue;
768                if (cand->GetIDDaughter(iDau) == iRP->GetID()) {
769                   if (fDebug > 1) printf("      was in RP set");
770                   iRP->SetForRPSelection(kFALSE);
771                   fFlowEvent->SetNumberOfRPs(fFlowEvent->GetNumberOfRPs() - 1);
772                }
773             }
774             if (fDebug > 1) printf("\n");
775          }
776          cand->SetForPOISelection(kTRUE);
777          fFlowEvent->InsertTrack(((AliFlowTrack*) cand));
778       }
779       if (fDebug > 0) printf("TPCevent %d\n", fFlowEvent->NumberOfTracks());
780       if(!fEventMixing) { // combinatorial background
781           for (Int_t pTracks = 0; pTracks < unp ; pTracks++) {
782              for (Int_t nTracks = pTracks + 1; nTracks < unp ; nTracks++) {
783                 if (fApplyDeltaDipCut && (!CheckDeltaDipAngle(up[pTracks], up[nTracks]))) continue;
784                 if (fCandidateEtaPtCut && (!CheckCandidateEtaPtCut(up[pTracks], up[nTracks]))) continue;
785                 PtSelector(1, up[pTracks], up[nTracks]);
786              }
787           }
788           for (Int_t nTracks = 0; nTracks < unn ; nTracks++) {
789              for (Int_t pTracks = nTracks + 1; pTracks < unn ; pTracks++) {
790                 if (fApplyDeltaDipCut && (!CheckDeltaDipAngle(un[nTracks], un[pTracks]))) continue;
791                 if (fCandidateEtaPtCut && (!CheckCandidateEtaPtCut(un[nTracks], un[pTracks]))) continue;
792                 PtSelector(2, un[nTracks], un[pTracks]);
793              }
794           }
795       }
796       if(fEventMixing) ReconstructionWithEventMixing(MixingCandidates);
797       PostData(1, fOutputList);
798       PostData(2, fFlowEvent);
799    }
800    if (fESD) AliFatal(" ESD analysis no longer supported ! " );
801 }
802 //_____________________________________________________________________________
803 void AliAnalysisTaskPhiFlow::ReconstructionWithEventMixing(TObjArray* MixingCandidates) const
804 {
805     // perform phi reconstruction with event mixing
806     if(fDebug > 0) cout << " *** ReconstructionWithEventMixing() *** " << endl;
807     AliEventPool* pool = fPoolManager->GetEventPool(fCentrality, fVertex);
808     if(!pool) AliFatal(Form("No pool found for centrality = %f, zVtx = %f", fCentrality, fVertex));
809     if(pool->IsReady() || pool->NTracksInPool() > fMixingParameters[1] / 10 || pool->GetCurrentNEvents() >= fMixingParameters[2]) {
810         Int_t nEvents = pool->GetCurrentNEvents();
811         if(fDebug > 0) cout << " --> " << nEvents << " events in mixing buffer ... " << endl;
812         for (Int_t iEvent(0); iEvent < nEvents; iEvent++) {
813             TObjArray* mixed_candidates = pool->GetEvent(iEvent);
814             if(!mixed_candidates) continue; // this should NEVER happen
815             Int_t bufferTracks = mixed_candidates->GetEntriesFast(); // buffered candidates
816             Int_t candidates = MixingCandidates->GetEntriesFast(); // mixing candidates
817             if(fDebug > 0) cout << Form("   - mixing %d tracks with %d buffered tracks from event %d ... ", candidates, bufferTracks, iEvent) << endl;
818             AliPhiMesonHelperTrack* buffer_un[bufferTracks]; // set values for buffered candidates
819             AliPhiMesonHelperTrack* buffer_up[bufferTracks];
820             Int_t buffer_unp(0);
821             Int_t buffer_unn(0);
822             AliPhiMesonHelperTrack* mix_un[candidates];// set values for mixing candidates
823             AliPhiMesonHelperTrack* mix_up[candidates];
824             Int_t mix_unp(0);
825             Int_t mix_unn(0);
826             for (Int_t iTracks = 0; iTracks < candidates; iTracks++) { // distinguish between k+ and k- for mixing candidates
827                 AliPhiMesonHelperTrack* track = (AliPhiMesonHelperTrack*)MixingCandidates->At(iTracks);
828                 if(!track) continue;
829                 if (track->Charge() > 0) {
830                     mix_up[mix_unp] = track;
831                     mix_unp++;
832                 }
833                 else if (track->Charge() < 0 ) {
834                    mix_un[mix_unn] = track;
835                    mix_unn++;
836                 }
837             }
838             for (Int_t iTracks = 0; iTracks < bufferTracks; iTracks++) { // distinguish between k+ and k- for buffered candidates
839                 AliPhiMesonHelperTrack* track = (AliPhiMesonHelperTrack*)mixed_candidates->At(iTracks);
840                 if(!track) continue;
841                 if (track->Charge() > 0) {
842                     buffer_up[buffer_unp] = track;
843                     buffer_unp++;
844                 }
845                 else if (track->Charge() < 0 ) {
846                    buffer_un[buffer_unn] = track;
847                    buffer_unn++;
848                 }
849             }
850             for (Int_t pMix = 0; pMix < mix_unp; pMix++) { // mix k+ (current event) k+ (buffer)
851                 if(fDebug > 1 ) cout << Form("mixing current k+ track %d with", pMix);
852                 if(!fTypeMixing) { // mix with like-sign kaons
853                     for(Int_t pBuf = 0; pBuf < buffer_unp; pBuf++) {
854                         if(fDebug > 1 ) cout << Form(" buffered k+ track %d", pBuf) << endl;
855                         if (fApplyDeltaDipCut && (!CheckDeltaDipAngle(mix_up[pMix], buffer_up[pBuf]))) continue;
856                         if (fCandidateMinEta && (!CheckCandidateEtaPtCut(mix_up[pMix], buffer_up[pBuf]))) continue;
857                         PtSelector(1, mix_up[pMix], buffer_up[pBuf]);
858                     }
859                 }
860                 else if(fTypeMixing) { // mix with unlike-sign kaons
861                     for(Int_t nBuf = 0; nBuf < buffer_unn; nBuf++) {
862                         if(fDebug > 1 ) cout << Form(" buffered k- track %d", nBuf) << endl;
863                         if (fApplyDeltaDipCut && (!CheckDeltaDipAngle(mix_up[pMix], buffer_un[nBuf]))) continue;
864                         if (fCandidateMinEta && (!CheckCandidateEtaPtCut(mix_up[pMix], buffer_un[nBuf]))) continue;
865                         PtSelector(2, mix_up[pMix], buffer_un[nBuf]);
866                     }
867                 }
868             }
869             for (Int_t nMix = 0; nMix < mix_unn; nMix++) { // mix k- (current event) k- (buffer)
870                 if(fDebug > 1 ) cout << Form("mixing current k- track %d with", nMix);
871                 if(!fTypeMixing) { // mix with like-sign kaons
872                     for(Int_t nBuf = 0; nBuf < buffer_unn; nBuf++) {
873                         if(fDebug > 1 ) cout << Form(" buffered k- track %d", nBuf) << endl;
874                         if (fApplyDeltaDipCut && (!CheckDeltaDipAngle(mix_un[nMix], buffer_un[nBuf]))) continue;
875                         if (fCandidateMinEta && (!CheckCandidateEtaPtCut(mix_un[nMix], buffer_un[nBuf]))) continue;
876                         PtSelector(2, mix_un[nMix], buffer_un[nBuf]);
877                     }
878                 }
879                 else if(fTypeMixing) { // mix with unlike-sign kaons
880                     for(Int_t pBuf = 0; pBuf < buffer_unp; pBuf++) {
881                         if(fDebug > 1 ) cout << Form(" buffered k+ track %d", pBuf) << endl;
882                         if (fApplyDeltaDipCut && (!CheckDeltaDipAngle(mix_un[nMix], buffer_up[pBuf]))) continue;
883                         if (fCandidateMinEta && (!CheckCandidateEtaPtCut(mix_un[nMix], buffer_up[pBuf]))) continue;
884                         PtSelector(1, mix_un[nMix], buffer_up[pBuf]);
885                     }
886                 }
887             }
888         } // end of mixed events loop
889     } // end of checking to see whether pool is filled correctly
890     pool->UpdatePool(MixingCandidates); // update pool with current mixing candidates (for next event)
891     if(fDebug > 0 ) pool->PrintInfo();
892 }
893 //_____________________________________________________________________________
894 void AliAnalysisTaskPhiFlow::Terminate(Option_t *)
895 {
896     // terminate
897     if(fDebug > 0) cout << " *** Terminate() *** " << endl;
898 }
899 //______________________________________________________________________________
900 void  AliAnalysisTaskPhiFlow::MakeTrack(Double_t mass, Double_t pt, Double_t phi, Double_t eta, Int_t nDau, Int_t iID[]) const
901 {
902    // Construct Flow Candidate Track from two selected candidates
903    if(fDebug > 1 ) cout << " *** MakeTrack() *** " << endl;
904    Bool_t overwrite = kTRUE;
905    AliFlowCandidateTrack *sTrack = static_cast<AliFlowCandidateTrack*>(fCandidates->At(fCandidates->GetLast() + 1));
906    if (!sTrack) {
907       sTrack = new AliFlowCandidateTrack(); //deleted by fCandidates
908       overwrite = kFALSE;
909    }
910    else sTrack->ClearMe();
911    sTrack->SetMass(mass);
912    sTrack->SetPt(pt);
913    sTrack->SetPhi(phi);
914    sTrack->SetEta(eta);
915    for (Int_t iDau = 0; iDau != nDau; ++iDau) sTrack->AddDaughter(iID[iDau]);
916    sTrack->SetForPOISelection(kTRUE);
917    sTrack->SetForRPSelection(kFALSE);
918    if (overwrite) fCandidates->SetLast(fCandidates->GetLast() + 1);
919    else fCandidates->AddLast(sTrack);
920    return;
921 }
922 //_____________________________________________________________________________
923 void AliAnalysisTaskPhiFlow::IsMC()
924 {
925     // Fill QA histos for MC analysis
926    TClonesArray *arrayMC = 0;
927    if(fDebug > 0) cout << " -> Switching to MC mode <- " << endl;
928    // fill array with mc tracks 
929    arrayMC = (TClonesArray*) fAOD->GetList()->FindObject(AliAODMCParticle::StdBranchName());
930    if (!arrayMC) AliFatal("Error: MC particles branch not found!\n");
931    for (Int_t iTracks = 0; iTracks < fAOD->GetNumberOfTracks(); iTracks++) {
932      AliAODTrack* track = fAOD->GetTrack(iTracks);
933      if(!PhiTrack(track) || !IsKaon(track)) { // check for kaons
934          if(fDebug>1) cout << " Rejected track" << endl;
935          continue;
936      }
937      if (fDebug>1) cout << " Received MC kaon " << endl;
938      Double_t b[2] = { -99., -99.};
939      Double_t bCov[3] = { -99., -99., -99.};
940      track->PropagateToDCA(fAOD->GetPrimaryVertex(), fAOD->GetMagneticField(), 100., b, bCov);
941      // find corresponding mc particle
942      AliAODMCParticle *partMC = (AliAODMCParticle*) arrayMC->At(TMath::Abs(track->GetLabel()));
943      if (!partMC) {
944          if(fDebug > 1) cout << "Cannot get MC particle" << endl;
945          continue;
946      }
947      // Check if it is primary, secondary from material or secondary from weak decay
948      Bool_t isPrimary           = partMC->IsPhysicalPrimary();
949      Bool_t isSecondaryMaterial = kFALSE;
950      Bool_t isSecondaryWeak     = kFALSE;
951      if (!isPrimary) {
952          Int_t mfl = -999, codemoth = -999;
953          Int_t indexMoth = partMC->GetMother();
954          if (indexMoth >= 0) { // is not fake
955             AliAODMCParticle* moth = (AliAODMCParticle*) arrayMC->At(indexMoth);
956             codemoth = TMath::Abs(moth->GetPdgCode());
957             mfl = Int_t(codemoth / TMath::Power(10, Int_t(TMath::Log10(codemoth))));
958          }
959          if (mfl == 3) isSecondaryWeak     = kTRUE;
960          else       isSecondaryMaterial = kTRUE;
961       }
962       if (isPrimary) {
963           fDCAPrim->Fill(track->Pt(), b[0]);
964           fDCAXYQA->Fill(b[0]);
965           fDCAZQA->Fill(b[1]);
966       }
967       if (isSecondaryWeak)  fDCASecondaryWeak->Fill(track->Pt(), b[0]);
968       if (isSecondaryMaterial) fDCAMaterial->Fill(track->Pt(), b[0]);
969    }
970 }
971 //_____________________________________________________________________________