]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/SPECTRA/AliAnalysisTaskCheckV0.cxx
Request by Martin: added flag for big output
[u/mrichter/AliRoot.git] / PWG2 / SPECTRA / AliAnalysisTaskCheckV0.cxx
1 /**************************************************************************
2  * Author: Boris Hippolyte.                                               *
3  * Contributors are mentioned in the code where appropriate.              *
4  *                                                                        *
5  * Permission to use, copy, modify and distribute this software and its   *
6  * documentation strictly for non-commercial purposes is hereby granted   *
7  * without fee, provided that the above copyright notice appears in all   *
8  * copies and that both the copyright notice and this permission notice   *
9  * appear in the supporting documentation. The authors make no claims     *
10  * about the suitability of this software for any purpose. It is          *
11  * provided "as is" without express or implied warranty.                  *
12  **************************************************************************/
13
14 //-----------------------------------------------------------------
15 //                 AliAnalysisTaskCheckV0 class
16 //            This task is for QAing the V0s from ESD/AOD
17 //              Origin: B.H. Nov2007, hippolyt@in2p3.fr
18 //-----------------------------------------------------------------
19 #include "TList.h"
20 #include "TH1F.h"
21 #include "TH2F.h"
22
23 #include "TStyle.h"
24 #include "TCanvas.h"
25 #include "TLegend.h"
26
27 #include "AliAnalysisTaskSE.h"
28 #include "AliAnalysisManager.h"
29 #include "AliInputEventHandler.h"
30
31 #include "AliESDEvent.h"
32 #include "AliESDVertex.h"
33 #include "AliAODEvent.h"
34
35 #include "AliESDv0.h"
36
37 #include "AliAnalysisTaskCheckV0.h"
38
39 ClassImp(AliAnalysisTaskCheckV0)
40
41 //________________________________________________________________________
42 AliAnalysisTaskCheckV0::AliAnalysisTaskCheckV0() 
43   : AliAnalysisTaskSE(), fAnalysisType("ESD"), fCollidingSystems(0), fUsePhysicsSelection(0),
44     fMaxPrimaryVtxPosZ(100.), fMinV0Pt(0.), fMaxV0Pt(100.), fMaxV0Rapidity(1.), fMinDaughterTpcClusters(80),
45     fListHist(), 
46     fHistPrimaryVertexPosX(0), fHistPrimaryVertexPosY(0), fHistPrimaryVertexPosZ(0),
47     fHistKeptPrimaryVertexPosX(0), fHistKeptPrimaryVertexPosY(0), fHistKeptPrimaryVertexPosZ(0),
48     fHistTrackMultiplicity(0), fHistV0Multiplicity(0), fHistV0OnFlyStatus(0),
49     fHistV0MultiplicityOff(0),
50     fHistV0Chi2Off(0),
51     fHistDcaV0DaughtersOff(0), fHistV0CosineOfPointingAngleOff(0),
52     fHistV0RadiusOff(0),fHistDcaV0ToPrimVertexOff(0),
53     fHistDcaPosToPrimVertexOff(0),fHistDcaNegToPrimVertexOff(0),
54     fHistMassK0sOff(0),fHistMassLambdaOff(0),fHistMassAntiLambdaOff(0),
55     fHistMassK0sOffVsPt(0),fHistMassLambdaOffVsPt(0),fHistMassAntiLambdaOffVsPt(0),
56     fHistArmenterosPodolanskiOff(0),
57     fHistV0MultiplicityOn(0),
58     fHistV0Chi2On(0),
59     fHistDcaV0DaughtersOn(0), fHistV0CosineOfPointingAngleOn(0),
60     fHistV0RadiusOn(0),fHistDcaV0ToPrimVertexOn(0),
61     fHistDcaPosToPrimVertexOn(0),fHistDcaNegToPrimVertexOn(0),
62     fHistMassK0sOn(0),fHistMassLambdaOn(0),fHistMassAntiLambdaOn(0),
63     fHistMassK0sOnVsPt(0),fHistMassLambdaOnVsPt(0),fHistMassAntiLambdaOnVsPt(0),
64     fHistArmenterosPodolanskiOn(0)
65 {
66   // Dummy constructor
67 }
68 //________________________________________________________________________
69 AliAnalysisTaskCheckV0::AliAnalysisTaskCheckV0(const char *name) 
70   : AliAnalysisTaskSE(name), fAnalysisType("ESD"), fCollidingSystems(0), fUsePhysicsSelection(0),
71     fMaxPrimaryVtxPosZ(100.), fMinV0Pt(0.), fMaxV0Pt(100.), fMaxV0Rapidity(1.), fMinDaughterTpcClusters(80),
72     fListHist(), 
73     fHistPrimaryVertexPosX(0), fHistPrimaryVertexPosY(0), fHistPrimaryVertexPosZ(0),
74     fHistKeptPrimaryVertexPosX(0), fHistKeptPrimaryVertexPosY(0), fHistKeptPrimaryVertexPosZ(0),
75     fHistTrackMultiplicity(0), fHistV0Multiplicity(0), fHistV0OnFlyStatus(0),
76     fHistV0MultiplicityOff(0),
77     fHistV0Chi2Off(0),
78     fHistDcaV0DaughtersOff(0), fHistV0CosineOfPointingAngleOff(0),
79     fHistV0RadiusOff(0),fHistDcaV0ToPrimVertexOff(0),
80     fHistDcaPosToPrimVertexOff(0),fHistDcaNegToPrimVertexOff(0),
81     fHistMassK0sOff(0),fHistMassLambdaOff(0),fHistMassAntiLambdaOff(0),
82     fHistMassK0sOffVsPt(0),fHistMassLambdaOffVsPt(0),fHistMassAntiLambdaOffVsPt(0),
83     fHistArmenterosPodolanskiOff(0),
84     fHistV0MultiplicityOn(0),
85     fHistV0Chi2On(0),
86     fHistDcaV0DaughtersOn(0), fHistV0CosineOfPointingAngleOn(0),
87     fHistV0RadiusOn(0),fHistDcaV0ToPrimVertexOn(0),
88     fHistDcaPosToPrimVertexOn(0),fHistDcaNegToPrimVertexOn(0),
89     fHistMassK0sOn(0),fHistMassLambdaOn(0),fHistMassAntiLambdaOn(0),
90     fHistMassK0sOnVsPt(0),fHistMassLambdaOnVsPt(0),fHistMassAntiLambdaOnVsPt(0),
91     fHistArmenterosPodolanskiOn(0)
92 {
93   // Constructor
94   // Define output slots only here
95   // Output slot #1 writes into a TList container
96   DefineOutput(1, TList::Class());
97 }
98 //________________________________________________________________________
99 AliAnalysisTaskCheckV0::~AliAnalysisTaskCheckV0(){
100   // Destructor
101   if (fListHist) { delete fListHist; fListHist = 0x0; }
102 }
103 //________________________________________________________________________
104 void AliAnalysisTaskCheckV0::UserCreateOutputObjects()
105 {
106   // Create histograms
107   // Called once
108
109   // Distinguish Track and V0 Multiplicity !
110
111   fListHist = new TList();
112   fListHist->SetOwner();
113  
114   if (!fHistPrimaryVertexPosX) {
115     fHistPrimaryVertexPosX = new TH1F("fHistPrimaryVertexPosX", "Primary vertex position in x;Position in x (cm);Events;",100,-1,1);
116     fListHist->Add(fHistPrimaryVertexPosX);
117   }
118   if (!fHistPrimaryVertexPosY) {
119     fHistPrimaryVertexPosY = new TH1F("fHistPrimaryVertexPosY", "Primary vertex position in y;Position in y (cm);Events;",100,-1,1);
120     fListHist->Add(fHistPrimaryVertexPosY);
121   }
122   if (!fHistPrimaryVertexPosZ) {
123     fHistPrimaryVertexPosZ = new TH1F("fHistPrimaryVertexPosZ", "Primary vertex position in z;Position in z (cm);Events;",100,-1,1);
124     fListHist->Add(fHistPrimaryVertexPosZ);
125   }
126   if (!fHistKeptPrimaryVertexPosX) {
127     fHistKeptPrimaryVertexPosX = new TH1F("fHistKeptPrimaryVertexPosX", "Kept Primary vertex position in x;Position in x (cm);Events;",100,-1,1);
128     fListHist->Add(fHistKeptPrimaryVertexPosX);
129   }
130   if (!fHistKeptPrimaryVertexPosY) {
131     fHistKeptPrimaryVertexPosY = new TH1F("fHistKeptPrimaryVertexPosY", "Kept Primary vertex position in y;Position in y (cm);Events;",100,-1,1);
132     fListHist->Add(fHistKeptPrimaryVertexPosY);
133   }
134   if (!fHistKeptPrimaryVertexPosZ) {
135     fHistKeptPrimaryVertexPosZ = new TH1F("fHistKeptPrimaryVertexPosZ", "Kept Primary vertex position in z;Position in z (cm);Events;",100,-1,1);
136     fListHist->Add(fHistKeptPrimaryVertexPosZ);
137   }
138
139   if (!fHistTrackMultiplicity) {
140     if (fCollidingSystems)
141       fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 200, 0, 40000);
142     else
143       fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 250, 0, 250);
144     fListHist->Add(fHistTrackMultiplicity);
145   }
146   if (!fHistV0Multiplicity) {
147     if (fCollidingSystems)
148       fHistV0Multiplicity = new TH1F("fHistV0Multiplicity", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
149     else
150       fHistV0Multiplicity = new TH1F("fHistV0Multiplicity", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
151     fListHist->Add(fHistV0Multiplicity);
152   }
153   if (!fHistV0OnFlyStatus) {
154     fHistV0OnFlyStatus = new TH1F("fHistV0OnFlyStatus", "V0 On fly status;status;Number of V0s", 3, 0, 3);
155     fListHist->Add(fHistV0OnFlyStatus);
156   }
157
158   // V0 offline distributions
159   if (!fHistV0MultiplicityOff) {
160     if (fCollidingSystems)
161       fHistV0MultiplicityOff = new TH1F("fHistV0MultiplicityOff", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
162     else
163       fHistV0MultiplicityOff = new TH1F("fHistV0MultiplicityOff", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50); 
164     fListHist->Add(fHistV0MultiplicityOff);
165   }
166   if (!fHistV0Chi2Off) {
167     fHistV0Chi2Off = new TH1F("fHistV0Chi2Off", "V0 chi2;chi2;Number of V0s", 33, 0, 33);
168     fListHist->Add(fHistV0Chi2Off);
169   }
170   if (!fHistDcaV0DaughtersOff) {
171     fHistDcaV0DaughtersOff = new TH1F("fHistDcaV0DaughtersOff", "DCA between V0 daughters;DCA (cm);Number of V0s", 300, 0, 3);
172     fListHist->Add(fHistDcaV0DaughtersOff);
173   }
174   if (!fHistV0CosineOfPointingAngleOff) {
175     fHistV0CosineOfPointingAngleOff = new TH1F("fHistV0CosineOfPointingAngleOff", "V0 Cosine of Pointing Angle;Number of V0s", 200, 0, 1);
176     fListHist->Add(fHistV0CosineOfPointingAngleOff);
177   }
178   if (!fHistV0RadiusOff) {
179     fHistV0RadiusOff = new TH1F("fHistV0RadiusOff", "V0 decay radius;Radius (cm);Number of V0s", 33, 0, 33);
180     fListHist->Add(fHistV0RadiusOff);
181   }
182   if (!fHistDcaV0ToPrimVertexOff) {
183     fHistDcaV0ToPrimVertexOff = new TH1F("fHistDcaV0ToPrimVertexOff", "DCA of V0 to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
184     fListHist->Add(fHistDcaV0ToPrimVertexOff);
185   }
186   if (!fHistDcaPosToPrimVertexOff) {
187     fHistDcaPosToPrimVertexOff = new TH1F("fHistDcaPosToPrimVertexOff", "DCA of V0 neg daughter to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
188     fListHist->Add(fHistDcaPosToPrimVertexOff);
189   }
190   if (!fHistDcaNegToPrimVertexOff) {
191     fHistDcaNegToPrimVertexOff = new TH1F("fHistDcaNegToPrimVertexOff", "DCA of V0 pos daughter to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
192     fListHist->Add(fHistDcaNegToPrimVertexOff);
193   }
194
195   if (!fHistMassK0sOff) {
196     fHistMassK0sOff = new TH1F("fHistMassK0sOff","K^{0} candidates;M(#pi^{+}#pi^{-}) (GeV/c^{2});Counts",100,0.4,0.6);
197     fListHist->Add(fHistMassK0sOff);
198   }
199   if (!fHistMassLambdaOff) {
200     fHistMassLambdaOff = new TH1F("fHistMassLambdaOff","#Lambda^{0} candidates;M(p#pi^{-}) (GeV/c^{2});Counts",75,1.05,1.2);
201     fListHist->Add(fHistMassLambdaOff);
202   }
203   if (!fHistMassAntiLambdaOff) {
204     fHistMassAntiLambdaOff = new TH1F("fHistMassAntiLambdaOff","#bar{#Lambda}^{0} candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});Counts",75,1.05,1.2);
205     fListHist->Add(fHistMassAntiLambdaOff);
206   }
207   if (!fHistMassK0sOffVsPt) {
208     fHistMassK0sOffVsPt = new TH2F("fHistMassK0sOffVsPt","K^{0} candidates;p_{T} (GeV/c);M(#pi^{+}#pi^{-}) (GeV/c^{2})",200,0,10,100,0.4,0.6);
209     fListHist->Add(fHistMassK0sOffVsPt);
210   }
211   if (!fHistMassLambdaOffVsPt) {
212     fHistMassLambdaOffVsPt = new TH2F("fHistMassLambdaOffVsPt","#Lambda^{0} candidates;p_{T} (GeV/c);M(p#pi^{-}) (GeV/c^{2})",200,0,10,75,1.05,1.2);
213     fListHist->Add(fHistMassLambdaOffVsPt);
214   }
215   if (!fHistMassAntiLambdaOffVsPt) {
216     fHistMassAntiLambdaOffVsPt = new TH2F("fHistMassAntiLambdaOffVsPt","#bar{#Lambda}^{0} candidates;p_{T} (GeV/c);M(#bar{p}#pi^{+}) (GeV/c^{2})",200,0,10,75,1.05,1.2);
217     fListHist->Add(fHistMassAntiLambdaOffVsPt);
218   }
219   if (!fHistArmenterosPodolanskiOff) {
220     fHistArmenterosPodolanskiOff   = new TH2F("fHistArmenterosPodolanskiOff","Armenteros-Podolanski Offline phase space;#alpha;p_{t} arm",200,-1.0,1.0,150,0,0.3);
221     fListHist->Add(fHistArmenterosPodolanskiOff);
222   }
223
224   // V0 on-the-fly distributions
225   if (!fHistV0MultiplicityOn) {
226     if (fCollidingSystems)
227       fHistV0MultiplicityOn = new TH1F("fHistV0MultiplicityOn", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
228     else
229       fHistV0MultiplicityOn = new TH1F("fHistV0MultiplicityOn", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
230     fListHist->Add(fHistV0MultiplicityOn);
231   }
232   if (!fHistV0Chi2On) {
233     fHistV0Chi2On = new TH1F("fHistV0Chi2On", "V0 chi2;chi2;Number of V0s", 33, 0, 33);
234     fListHist->Add(fHistV0Chi2On);
235   }
236   if (!fHistDcaV0DaughtersOn) {
237     fHistDcaV0DaughtersOn = new TH1F("fHistDcaV0DaughtersOn", "DCA between V0 daughters;DCA (cm);Number of V0s", 300, 0, 3);
238     fListHist->Add(fHistDcaV0DaughtersOn);
239   }
240   if (!fHistV0CosineOfPointingAngleOn) {
241     fHistV0CosineOfPointingAngleOn = new TH1F("fHistV0CosineOfPointingAngleOn", "V0 Cosine of Pointing Angle;Number of V0s", 200, 0, 1);
242     fListHist->Add(fHistV0CosineOfPointingAngleOn);
243   }
244   if (!fHistV0RadiusOn) {
245     fHistV0RadiusOn = new TH1F("fHistV0RadiusOn", "V0 decay radius;Radius (cm);Number of V0s", 33, 0, 33);
246     fListHist->Add(fHistV0RadiusOn);
247   }
248   if (!fHistDcaV0ToPrimVertexOn) {
249     fHistDcaV0ToPrimVertexOn = new TH1F("fHistDcaV0ToPrimVertexOn", "DCA of V0 to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
250     fListHist->Add(fHistDcaV0ToPrimVertexOn);
251   }
252   if (!fHistDcaPosToPrimVertexOn) {
253     fHistDcaPosToPrimVertexOn = new TH1F("fHistDcaPosToPrimVertexOn", "DCA of V0 neg daughter to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
254     fListHist->Add(fHistDcaPosToPrimVertexOn);
255   }
256   if (!fHistDcaNegToPrimVertexOn) {
257     fHistDcaNegToPrimVertexOn = new TH1F("fHistDcaNegToPrimVertexOn", "DCA of V0 pos daughter to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
258     fListHist->Add(fHistDcaNegToPrimVertexOn);
259   }
260
261   if (!fHistMassK0sOn) {
262     fHistMassK0sOn = new TH1F("fHistMassK0sOn","K^{0} candidates;M(#pi^{+}#pi^{-}) (GeV/c^{2});Counts",100,0.4,0.6);
263     fListHist->Add(fHistMassK0sOn);
264   }
265   if (!fHistMassLambdaOn) {
266     fHistMassLambdaOn = new TH1F("fHistMassLambdaOn","#Lambda^{0} candidates;M(p#pi^{-}) (GeV/c^{2});Counts",75,1.05,1.2);
267     fListHist->Add(fHistMassLambdaOn);
268   }
269   if (!fHistMassAntiLambdaOn) {
270     fHistMassAntiLambdaOn = new TH1F("fHistMassAntiLambdaOn","#bar{#Lambda}^{0} candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});Counts",75,1.05,1.2);
271     fListHist->Add(fHistMassAntiLambdaOn);
272   }
273   if (!fHistMassK0sOnVsPt) {
274     fHistMassK0sOnVsPt = new TH2F("fHistMassK0sOnVsPt","K^{0} candidates;p_{T} (GeV/c);M(#pi^{+}#pi^{-}) (GeV/c^{2})",200,0,10,100,0.4,0.6);
275     fListHist->Add(fHistMassK0sOnVsPt);
276   }
277   if (!fHistMassLambdaOnVsPt) {
278     fHistMassLambdaOnVsPt = new TH2F("fHistMassLambdaOnVsPt","#Lambda^{0} candidates;p_{T} (GeV/c);M(p#pi^{-}) (GeV/c^{2})",200,0,10,75,1.05,1.2);
279     fListHist->Add(fHistMassLambdaOnVsPt);
280   }
281   if (!fHistMassAntiLambdaOnVsPt) {
282     fHistMassAntiLambdaOnVsPt = new TH2F("fHistMassAntiLambdaOnVsPt","#bar{#Lambda}^{0} candidates;p_{T} (GeV/c);M(#bar{p}#pi^{+}) (GeV/c^{2})",200,0,10,75,1.05,1.2);
283     fListHist->Add(fHistMassAntiLambdaOnVsPt);
284   }
285   if (!fHistArmenterosPodolanskiOn) {
286     fHistArmenterosPodolanskiOn    = new TH2F("fHistArmenterosPodolanskiOn","Armenteros-Podolanski Onthefly phase space;#alpha;p_{t} arm",200,-1.0,1.0,150,0,0.3);
287     fListHist->Add(fHistArmenterosPodolanskiOn);
288   }
289   // Post output data.
290   PostData(1, fListHist);
291 }
292
293 //________________________________________________________________________
294 void AliAnalysisTaskCheckV0::UserExec(Option_t *) 
295 {
296   // Main loop
297   // Called for each event
298   AliVEvent* lEvent = InputEvent();
299   if (!lEvent) {
300     Printf("ERROR: Event not available");
301     return;
302   }
303
304   Bool_t isSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
305   if ((fUsePhysicsSelection)&&(!isSelected)) return;
306
307   fHistTrackMultiplicity->Fill(lEvent->GetNumberOfTracks());
308   
309   Double_t tPrimaryVtxPosition[3];
310   Int_t nv0s = 0;
311   nv0s = lEvent->GetNumberOfV0s();
312   //  Printf("CheckV0 analysis task: There are %d v0s in this event",nv0s);
313
314   Int_t    lOnFlyStatus = 0, nv0sOn = 0, nv0sOff = 0;
315   Double_t lChi2V0 = 0;
316   Double_t lDcaV0Daughters = 0, lDcaV0ToPrimVertex = 0;
317   Double_t lDcaPosToPrimVertex = 0, lDcaNegToPrimVertex = 0;
318   Double_t lV0CosineOfPointingAngle = 0;
319   Double_t lV0Radius = 0, lPt = 0;
320   Double_t lRapK0Short = 0, lRapLambda = 0;
321   Double_t lInvMassK0s = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
322   Double_t lAlphaV0 = 0, lPtArmV0 = 0;
323
324   const AliVVertex *primaryVtx = lEvent->GetPrimaryVertex();
325   tPrimaryVtxPosition[0] = primaryVtx->GetX();
326   tPrimaryVtxPosition[1] = primaryVtx->GetY();
327   tPrimaryVtxPosition[2] = primaryVtx->GetZ();
328   fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
329   fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
330   fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
331
332   if (TMath::Abs(tPrimaryVtxPosition[2])<fMaxPrimaryVtxPosZ){// event selections
333
334     fHistKeptPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
335     fHistKeptPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
336     fHistKeptPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
337
338     if(fAnalysisType == "ESD") {
339
340       Double_t lMagneticField = ((AliESDEvent*)lEvent)->GetMagneticField();
341
342       for (Int_t iV0 = 0; iV0 < nv0s; iV0++) 
343         {// This is the begining of the V0 loop
344           AliESDv0 *v0 = ((AliESDEvent*)lEvent)->GetV0(iV0);
345           if (!v0) continue;
346
347           Double_t tDecayVertexV0[3]; v0->GetXYZ(tDecayVertexV0[0],tDecayVertexV0[1],tDecayVertexV0[2]); 
348
349           lV0Radius = TMath::Sqrt(tDecayVertexV0[0]*tDecayVertexV0[0]+tDecayVertexV0[1]*tDecayVertexV0[1]);
350           lPt = v0->Pt();
351           lRapK0Short = v0->RapK0Short();
352           lRapLambda  = v0->RapLambda();
353           if ((lPt<fMinV0Pt)||(fMaxV0Pt<lPt)) continue;
354
355           UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
356           UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
357
358           Double_t lMomPos[3]; v0->GetPPxPyPz(lMomPos[0],lMomPos[1],lMomPos[2]);
359           Double_t lMomNeg[3]; v0->GetNPxPyPz(lMomNeg[0],lMomNeg[1],lMomNeg[2]);
360
361           AliESDtrack *pTrack=((AliESDEvent*)lEvent)->GetTrack(lKeyPos);
362           AliESDtrack *nTrack=((AliESDEvent*)lEvent)->GetTrack(lKeyNeg);
363           if (!pTrack || !nTrack) {
364             Printf("ERROR: Could not retreive one of the daughter track");
365             continue;
366           }
367
368           // Filter like-sign V0 (next: add counter and distribution)
369           if ( pTrack->GetSign() == nTrack->GetSign()){
370             continue;
371           } 
372
373           // WARNING: the following selections cannot be done for AOD yet...
374
375           // Tracks quality cuts 
376           if ( ( (pTrack->GetTPCNcls()) < fMinDaughterTpcClusters ) || ( (nTrack->GetTPCNcls()) < fMinDaughterTpcClusters ) ) continue;
377         
378           // TPC refit condition (done during reconstruction for Offline but not for On-the-fly)
379           if( !(pTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;      
380           if( !(nTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
381
382
383           lDcaPosToPrimVertex = TMath::Abs(pTrack->GetD(tPrimaryVtxPosition[0],
384                                                         tPrimaryVtxPosition[1],
385                                                         lMagneticField) );
386
387           lDcaNegToPrimVertex = TMath::Abs(nTrack->GetD(tPrimaryVtxPosition[0],
388                                                         tPrimaryVtxPosition[1],
389                                                         lMagneticField) );
390
391           lOnFlyStatus = v0->GetOnFlyStatus();
392           lChi2V0 = v0->GetChi2V0();
393           lDcaV0Daughters = v0->GetDcaV0Daughters();
394           lDcaV0ToPrimVertex = v0->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
395           lV0CosineOfPointingAngle = v0->GetV0CosineOfPointingAngle(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
396
397           // Getting invariant mass infos directly from ESD
398           v0->ChangeMassHypothesis(310);
399           lInvMassK0s = v0->GetEffMass();
400           v0->ChangeMassHypothesis(3122);
401           lInvMassLambda = v0->GetEffMass();
402           v0->ChangeMassHypothesis(-3122);
403           lInvMassAntiLambda = v0->GetEffMass();
404           lAlphaV0 = v0->AlphaV0();
405           lPtArmV0 = v0->PtArmV0();
406
407           fHistV0OnFlyStatus->Fill(lOnFlyStatus);
408           if(!lOnFlyStatus){
409             nv0sOff++;
410             fHistV0Chi2Off->Fill(lChi2V0);
411             fHistDcaV0ToPrimVertexOff->Fill(lDcaV0ToPrimVertex);
412             fHistDcaV0DaughtersOff->Fill(lDcaV0Daughters);
413             fHistV0CosineOfPointingAngleOff->Fill(lV0CosineOfPointingAngle);
414
415             fHistV0RadiusOff->Fill(lV0Radius);
416             fHistDcaPosToPrimVertexOff->Fill(lDcaPosToPrimVertex);
417             fHistDcaNegToPrimVertexOff->Fill(lDcaNegToPrimVertex);
418
419             // Filling invariant mass histos for all candidates
420             fHistMassK0sOff->Fill(lInvMassK0s);
421             fHistMassLambdaOff->Fill(lInvMassLambda);
422             fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
423             if (TMath::Abs(lRapK0Short)<fMaxV0Rapidity)
424               fHistMassK0sOffVsPt->Fill(lPt,lInvMassK0s);
425             if (TMath::Abs(lRapLambda)<fMaxV0Rapidity){
426               fHistMassLambdaOffVsPt->Fill(lPt,lInvMassLambda);
427               fHistMassAntiLambdaOffVsPt->Fill(lPt,lInvMassAntiLambda);
428             }
429             fHistArmenterosPodolanskiOff->Fill(lAlphaV0,lPtArmV0);
430           }
431           else {
432             nv0sOn++;
433             fHistV0Chi2On->Fill(lChi2V0);
434             fHistDcaV0ToPrimVertexOn->Fill(lDcaV0ToPrimVertex);
435             fHistDcaV0DaughtersOn->Fill(lDcaV0Daughters);
436             fHistV0CosineOfPointingAngleOn->Fill(lV0CosineOfPointingAngle);
437
438             fHistV0RadiusOn->Fill(lV0Radius);
439             fHistDcaPosToPrimVertexOn->Fill(lDcaPosToPrimVertex);
440             fHistDcaNegToPrimVertexOn->Fill(lDcaNegToPrimVertex);
441
442             // Filling invariant mass histos for all candidates
443             fHistMassK0sOn->Fill(lInvMassK0s);
444             fHistMassLambdaOn->Fill(lInvMassLambda);
445             fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
446             if (TMath::Abs(lRapK0Short)<fMaxV0Rapidity)
447               fHistMassK0sOnVsPt->Fill(lPt,lInvMassK0s);
448             if (TMath::Abs(lRapLambda)<fMaxV0Rapidity){
449               fHistMassLambdaOnVsPt->Fill(lPt,lInvMassLambda);
450               fHistMassAntiLambdaOnVsPt->Fill(lPt,lInvMassAntiLambda);
451             }
452             fHistArmenterosPodolanskiOn->Fill(lAlphaV0,lPtArmV0);
453           }
454         }// This is the end of the V0 loop
455     } // end of "ESD" analysis
456
457     else if(fAnalysisType == "AOD") {
458
459       for (Int_t iV0 = 0; iV0 < nv0s; iV0++) 
460         {// This is the begining of the V0 loop
461           AliAODv0 *v0 = ((AliAODEvent*)lEvent)->GetV0(iV0);
462           if (!v0) continue;
463
464           lV0Radius = v0->RadiusV0();
465           lPt = v0->Pt();
466           lRapK0Short = v0->RapK0Short();
467           lRapLambda  = v0->RapLambda();
468           if ((lPt<fMinV0Pt)||(fMaxV0Pt<lPt)) continue;
469           lDcaPosToPrimVertex = v0->DcaPosToPrimVertex();
470           lDcaNegToPrimVertex = v0->DcaNegToPrimVertex();
471
472           // propose getOnFlyStatus() = -1 for like-sign
473
474           lOnFlyStatus = v0->GetOnFlyStatus();
475           lChi2V0 = v0->Chi2V0();
476           lDcaV0Daughters = v0->DcaV0Daughters();
477           lDcaV0ToPrimVertex = v0->DcaV0ToPrimVertex();
478           lV0CosineOfPointingAngle = v0->CosPointingAngle(tPrimaryVtxPosition);
479
480           lInvMassK0s = v0->MassK0Short();
481           lInvMassLambda = v0->MassLambda();
482           lInvMassAntiLambda = v0->MassAntiLambda();
483           lAlphaV0 = v0->AlphaV0();
484           lPtArmV0 = v0->PtArmV0();
485
486           fHistV0OnFlyStatus->Fill(lOnFlyStatus);
487           if(!lOnFlyStatus){
488             nv0sOff++;
489             fHistV0Chi2Off->Fill(lChi2V0);
490             fHistDcaV0ToPrimVertexOff->Fill(lDcaV0ToPrimVertex);
491             fHistDcaV0DaughtersOff->Fill(lDcaV0Daughters);
492             fHistV0CosineOfPointingAngleOff->Fill(lV0CosineOfPointingAngle);
493
494             fHistV0RadiusOff->Fill(lV0Radius);
495             fHistDcaPosToPrimVertexOff->Fill(lDcaPosToPrimVertex);
496             fHistDcaNegToPrimVertexOff->Fill(lDcaNegToPrimVertex);
497
498             // Filling invariant mass histos for all candidates
499             fHistMassK0sOff->Fill(lInvMassK0s);
500             fHistMassLambdaOff->Fill(lInvMassLambda);
501             fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
502             if (TMath::Abs(lRapK0Short)<fMaxV0Rapidity)
503               fHistMassK0sOffVsPt->Fill(lPt,lInvMassK0s);
504             if (TMath::Abs(lRapLambda)<fMaxV0Rapidity){
505               fHistMassLambdaOffVsPt->Fill(lPt,lInvMassLambda);
506               fHistMassAntiLambdaOffVsPt->Fill(lPt,lInvMassAntiLambda);
507             }
508             fHistArmenterosPodolanskiOff->Fill(lAlphaV0,lPtArmV0);
509           }
510           else {
511             nv0sOn++;
512             fHistV0Chi2On->Fill(lChi2V0);
513             fHistDcaV0ToPrimVertexOn->Fill(lDcaV0ToPrimVertex);
514             fHistDcaV0DaughtersOn->Fill(lDcaV0Daughters);
515             fHistV0CosineOfPointingAngleOn->Fill(lV0CosineOfPointingAngle);
516
517             fHistV0RadiusOn->Fill(lV0Radius);
518             fHistDcaPosToPrimVertexOn->Fill(lDcaPosToPrimVertex);
519             fHistDcaNegToPrimVertexOn->Fill(lDcaNegToPrimVertex);
520
521             // Filling invariant mass histos for all candidates
522             fHistMassK0sOn->Fill(lInvMassK0s);
523             fHistMassLambdaOn->Fill(lInvMassLambda);
524             fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
525             if (TMath::Abs(lRapK0Short)<fMaxV0Rapidity)
526               fHistMassK0sOnVsPt->Fill(lPt,lInvMassK0s);
527             if (TMath::Abs(lRapLambda)<fMaxV0Rapidity){
528               fHistMassLambdaOnVsPt->Fill(lPt,lInvMassLambda);
529               fHistMassAntiLambdaOnVsPt->Fill(lPt,lInvMassAntiLambda);
530             }
531             fHistArmenterosPodolanskiOn->Fill(lAlphaV0,lPtArmV0);
532           }
533         }// This is the end of the V0 loop
534     } // end of "AOD" analysis
535
536     fHistV0Multiplicity->Fill(nv0s);
537     fHistV0MultiplicityOff->Fill(nv0sOff);
538     fHistV0MultiplicityOn->Fill(nv0sOn);
539
540   }// end of event selection
541
542   // Post output data.
543   PostData(1, fListHist);
544 }
545
546 //________________________________________________________________________
547 void AliAnalysisTaskCheckV0::Terminate(Option_t *) 
548 {
549   // Draw result to the screen
550   // Called once at the end of the query
551
552   TList *cRetrievedList = 0x0;
553   cRetrievedList = (TList*)GetOutputData(1);
554
555   if(!cRetrievedList){
556     AliWarning("ERROR - AliAnalysisTaskCheckV0: output data container list not available\n"); return;
557   }
558
559   // Implement a decent style
560   TStyle *myStyle = new TStyle("myStyle","my style");
561   Int_t font = 42;
562   myStyle->SetCanvasColor(10);
563   myStyle->SetStatColor(10);
564   myStyle->SetPadColor(10);
565   myStyle->SetDrawBorder(0);
566   myStyle->SetCanvasBorderMode(0);
567   myStyle->SetPadBorderMode(0);
568   myStyle->SetTextFont(font);
569   myStyle->SetStatFont(font);
570   myStyle->SetLabelFont(font,"xyz"); 
571   myStyle->SetTitleFont(font);  
572   myStyle->SetPalette(1,0); 
573   myStyle->cd();
574
575   fHistTrackMultiplicity = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistTrackMultiplicity"));
576   if (!fHistTrackMultiplicity) {
577     Printf("ERROR: fHistTrackMultiplicity not available");
578     return;
579   }
580   fHistV0Multiplicity = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistV0Multiplicity"));
581   if (!fHistV0Multiplicity) {
582     Printf("ERROR: fHistV0Multiplicity not available");
583     return;
584   }
585   fHistV0MultiplicityOff = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistV0MultiplicityOff"));
586   if (!fHistV0MultiplicityOff) {
587     Printf("ERROR: fHistV0MultiplicityOff not available");
588     return;
589   }
590   fHistV0MultiplicityOn = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistV0MultiplicityOn"));
591   if (!fHistV0MultiplicityOn) {
592     Printf("ERROR: fHistV0MultiplicityOn not available");
593     return;
594   }
595
596   TCanvas *canCheckV0 = new TCanvas("AliAnalysisTaskCheckV0","Check V0",10,10,510,700);
597   canCheckV0->Divide(2,3);
598   if (fHistTrackMultiplicity->GetMaximum() > 0.) canCheckV0->cd(1)->SetLogy();
599   fHistTrackMultiplicity->SetMarkerStyle(26);
600   fHistTrackMultiplicity->DrawCopy("E");
601   fHistV0Multiplicity->SetMarkerStyle(25);
602   fHistV0Multiplicity->DrawCopy("ESAME");
603   fHistV0MultiplicityOff->SetMarkerStyle(24);
604   fHistV0MultiplicityOff->DrawCopy("ESAME");
605   fHistV0MultiplicityOn->SetMarkerStyle(20);
606   fHistV0MultiplicityOn->DrawCopy("ESAME");
607
608   TLegend *legendMultiplicity = new TLegend(0.5,0.5,0.75,0.75);
609   legendMultiplicity->AddEntry(fHistTrackMultiplicity,"tracks");
610   legendMultiplicity->AddEntry(fHistV0Multiplicity,"all V^{0}");
611   legendMultiplicity->AddEntry(fHistV0MultiplicityOff,"offline V^{0}");
612   legendMultiplicity->AddEntry(fHistV0MultiplicityOn,"onthefly V^{0}");
613   legendMultiplicity->Draw();
614
615   fHistMassK0sOff = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistMassK0sOff"));
616   if (!fHistMassK0sOff) {
617     Printf("ERROR: fHistMassK0sOff not available");
618     return;
619   }
620   fHistMassK0sOn = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistMassK0sOn"));
621   if (!fHistMassK0sOn) {
622     Printf("ERROR: fHistMassK0sOn not available");
623     return;
624   }
625   fHistMassLambdaOff = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistMassLambdaOff"));
626   if (!fHistMassLambdaOff) {
627     Printf("ERROR: fHistMassLambdaOff not available");
628     return;
629   }
630   fHistMassLambdaOn = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistMassLambdaOn"));
631   if (!fHistMassLambdaOn) {
632     Printf("ERROR: fHistMassLambdaOn not available");
633     return;
634   }
635   fHistMassAntiLambdaOff = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistMassAntiLambdaOff"));
636   if (!fHistMassAntiLambdaOff) {
637     Printf("ERROR: fHistMassAntiLambdaOff not available");
638     return;
639   }
640   fHistMassAntiLambdaOn = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistMassAntiLambdaOn"));
641   if (!fHistMassAntiLambdaOn) {
642     Printf("ERROR: fHistMassAntiLambdaOn not available");
643     return;
644   }
645   fHistArmenterosPodolanskiOff = dynamic_cast<TH2F*> (cRetrievedList->FindObject("fHistArmenterosPodolanskiOff"));
646   if (!fHistArmenterosPodolanskiOff) {
647     Printf("ERROR: fHistArmenterosPodolanskiOff not available");
648     return;
649   }
650   fHistArmenterosPodolanskiOn = dynamic_cast<TH2F*> (cRetrievedList->FindObject("fHistArmenterosPodolanskiOn"));
651   if (!fHistArmenterosPodolanskiOn) {
652     Printf("ERROR: fHistArmenterosPodolanskiOn not available");
653     return;
654   }
655
656   canCheckV0->cd(2);
657   fHistMassK0sOn->SetMarkerStyle(20);
658   fHistMassK0sOn->DrawCopy("E");
659   fHistMassK0sOff->SetMarkerStyle(24);
660   fHistMassK0sOff->DrawCopy("ESAME");
661
662   canCheckV0->cd(3);
663   fHistMassLambdaOn->SetMarkerStyle(20);
664   fHistMassLambdaOn->DrawCopy("E");
665   fHistMassLambdaOff->SetMarkerStyle(24);
666   fHistMassLambdaOff->DrawCopy("ESAME");
667
668   canCheckV0->cd(4);
669   fHistMassAntiLambdaOn->SetMarkerStyle(20);
670   fHistMassAntiLambdaOn->DrawCopy("E");
671   fHistMassAntiLambdaOff->SetMarkerStyle(24);
672   fHistMassAntiLambdaOff->DrawCopy("ESAME");
673
674   canCheckV0->cd(5);
675   fHistArmenterosPodolanskiOff->DrawCopy("COL2Z");
676   canCheckV0->cd(6);
677   fHistArmenterosPodolanskiOn->DrawCopy("COL2Z");
678
679 }