]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/SPECTRA/AliAnalysisTaskCheckV0.cxx
bda9e8bbf85cf1722fc19eab7fcdb360cd5f419c
[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 "TCanvas.h"
22 #include "TLegend.h"
23
24 #include "AliAnalysisTaskSE.h"
25
26 #include "AliESDEvent.h"
27 #include "AliESDVertex.h"
28 #include "AliAODEvent.h"
29 #include "AliAODVertex.h"
30
31 #include "AliESDv0.h"
32
33 #include "AliAnalysisTaskCheckV0.h"
34
35 ClassImp(AliAnalysisTaskCheckV0)
36
37 //________________________________________________________________________
38 AliAnalysisTaskCheckV0::AliAnalysisTaskCheckV0() 
39   : AliAnalysisTaskSE(), fAnalysisType("ESD"), fCollidingSystems(0), fListHist(), 
40     fHistPrimaryVertexPosX(0), fHistPrimaryVertexPosY(0), fHistPrimaryVertexPosZ(0),
41     fHistTrackMultiplicity(0), fHistV0Multiplicity(0), fHistV0OnFlyStatus(0),
42     fHistV0MultiplicityOff(0), fHistV0Chi2Off(0),
43     fHistDcaV0DaughtersOff(0), fHistV0CosineOfPointingAngleOff(0),
44     fHistV0RadiusOff(0),fHistDcaV0ToPrimVertexOff(0),
45     fHistDcaPosToPrimVertexOff(0),fHistDcaNegToPrimVertexOff(0),
46     fHistMassK0Off(0),fHistMassLambdaOff(0),fHistMassAntiLambdaOff(0),
47     fHistV0MultiplicityOn(0), fHistV0Chi2On(0),
48     fHistDcaV0DaughtersOn(0), fHistV0CosineOfPointingAngleOn(0),
49     fHistV0RadiusOn(0),fHistDcaV0ToPrimVertexOn(0),
50     fHistDcaPosToPrimVertexOn(0),fHistDcaNegToPrimVertexOn(0),
51     fHistMassK0On(0),fHistMassLambdaOn(0),fHistMassAntiLambdaOn(0)
52 {
53   // Dummy constructor
54 }
55 //________________________________________________________________________
56 AliAnalysisTaskCheckV0::AliAnalysisTaskCheckV0(const char *name) 
57   : AliAnalysisTaskSE(name), fAnalysisType("ESD"), fCollidingSystems(0), fListHist(), 
58     fHistPrimaryVertexPosX(0), fHistPrimaryVertexPosY(0), fHistPrimaryVertexPosZ(0),
59     fHistTrackMultiplicity(0), fHistV0Multiplicity(0), fHistV0OnFlyStatus(0),
60     fHistV0MultiplicityOff(0), fHistV0Chi2Off(0),
61     fHistDcaV0DaughtersOff(0), fHistV0CosineOfPointingAngleOff(0),
62     fHistV0RadiusOff(0),fHistDcaV0ToPrimVertexOff(0),
63     fHistDcaPosToPrimVertexOff(0),fHistDcaNegToPrimVertexOff(0),
64     fHistMassK0Off(0),fHistMassLambdaOff(0),fHistMassAntiLambdaOff(0),
65     fHistV0MultiplicityOn(0), fHistV0Chi2On(0),
66     fHistDcaV0DaughtersOn(0), fHistV0CosineOfPointingAngleOn(0),
67     fHistV0RadiusOn(0),fHistDcaV0ToPrimVertexOn(0),
68     fHistDcaPosToPrimVertexOn(0),fHistDcaNegToPrimVertexOn(0),
69     fHistMassK0On(0),fHistMassLambdaOn(0),fHistMassAntiLambdaOn(0)
70 {
71   // Constructor
72   // Define output slots only here
73   // Output slot #1 writes into a TList container
74   DefineOutput(1, TList::Class());
75 }
76
77 //________________________________________________________________________
78 void AliAnalysisTaskCheckV0::UserCreateOutputObjects()
79 {
80   // Create histograms
81   // Called once
82
83   // Distinguish Track and V0 Multiplicity !
84
85   fListHist = new TList();
86  
87   if (!fHistPrimaryVertexPosX) {
88     fHistPrimaryVertexPosX = new TH1F("fHistPrimaryVertexPosX", "Primary vertex position in x;Position in x (cm);Events;",100,-1,1);
89     fListHist->Add(fHistPrimaryVertexPosX);
90   }
91   if (!fHistPrimaryVertexPosY) {
92     fHistPrimaryVertexPosY = new TH1F("fHistPrimaryVertexPosY", "Primary vertex position in y;Position in y (cm);Events;",100,-1,1);
93     fListHist->Add(fHistPrimaryVertexPosY);
94   }
95   if (!fHistPrimaryVertexPosZ) {
96     fHistPrimaryVertexPosZ = new TH1F("fHistPrimaryVertexPosZ", "Primary vertex position in z;Position in z (cm);Events;",100,-1,1);
97     fListHist->Add(fHistPrimaryVertexPosZ);
98   }
99
100   if (!fHistTrackMultiplicity) {
101     if (fCollidingSystems)
102       fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 200, 0, 40000);
103     else
104       fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 250, 0, 250);
105     fListHist->Add(fHistTrackMultiplicity);
106   }
107   if (!fHistV0Multiplicity) {
108     if (fCollidingSystems)
109       fHistV0Multiplicity = new TH1F("fHistV0Multiplicity", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
110     else
111       fHistV0Multiplicity = new TH1F("fHistV0Multiplicity", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
112     fListHist->Add(fHistV0Multiplicity);
113   }
114   if (!fHistV0OnFlyStatus) {
115     fHistV0OnFlyStatus = new TH1F("fHistV0OnFlyStatus", "V0 On fly status;status;Number of V0s", 3, 0, 3);
116     fListHist->Add(fHistV0OnFlyStatus);
117   }
118
119   // V0 offline distributions
120   if (!fHistV0MultiplicityOff) {
121     if (fCollidingSystems)
122       fHistV0MultiplicityOff = new TH1F("fHistV0MultiplicityOff", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
123     else
124       fHistV0MultiplicityOff = new TH1F("fHistV0MultiplicityOff", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50); 
125     fListHist->Add(fHistV0MultiplicityOff);
126   }
127   if (!fHistV0Chi2Off) {
128     fHistV0Chi2Off = new TH1F("fHistV0Chi2Off", "V0 chi2;chi2;Number of V0s", 33, 0, 33);
129     fListHist->Add(fHistV0Chi2Off);
130   }
131   if (!fHistDcaV0DaughtersOff) {
132     fHistDcaV0DaughtersOff = new TH1F("fHistDcaV0DaughtersOff", "DCA between V0 daughters;DCA (cm);Number of V0s", 300, 0, 3);
133     fListHist->Add(fHistDcaV0DaughtersOff);
134   }
135   if (!fHistV0CosineOfPointingAngleOff) {
136     fHistV0CosineOfPointingAngleOff = new TH1F("fHistV0CosineOfPointingAngleOff", "V0 Cosine of Pointing Angle;Number of V0s", 200, 0, 1);
137     fListHist->Add(fHistV0CosineOfPointingAngleOff);
138   }
139   if (!fHistV0RadiusOff) {
140     fHistV0RadiusOff = new TH1F("fHistV0RadiusOff", "V0 decay radius;Radius (cm);Number of V0s", 33, 0, 33);
141     fListHist->Add(fHistV0RadiusOff);
142   }
143   if (!fHistDcaV0ToPrimVertexOff) {
144     fHistDcaV0ToPrimVertexOff = new TH1F("fHistDcaV0ToPrimVertexOff", "DCA of V0 to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
145     fListHist->Add(fHistDcaV0ToPrimVertexOff);
146   }
147   if (!fHistDcaPosToPrimVertexOff) {
148     fHistDcaPosToPrimVertexOff = new TH1F("fHistDcaPosToPrimVertexOff", "DCA of V0 neg daughter to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
149     fListHist->Add(fHistDcaPosToPrimVertexOff);
150   }
151   if (!fHistDcaNegToPrimVertexOff) {
152     fHistDcaNegToPrimVertexOff = new TH1F("fHistDcaNegToPrimVertexOff", "DCA of V0 pos daughter to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
153     fListHist->Add(fHistDcaNegToPrimVertexOff);
154   }
155
156   if (!fHistMassK0Off) {
157     fHistMassK0Off = new TH1F("fHistMassK0Off","K^{0} candidates;M(#pi^{+}#pi^{-}) (GeV/c^{2});Counts",100,0.4,0.6);
158     fListHist->Add(fHistMassK0Off);
159   }
160   if (!fHistMassLambdaOff) {
161     fHistMassLambdaOff = new TH1F("fHistMassLambdaOff","#Lambda^{0} candidates;M(p#pi^{-}) (GeV/c^{2});Counts",75,1.05,1.2);
162     fListHist->Add(fHistMassLambdaOff);
163   }
164   if (!fHistMassAntiLambdaOff) {
165     fHistMassAntiLambdaOff = new TH1F("fHistMassAntiLambdaOff","#bar{#Lambda}^{0} candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});Counts",75,1.05,1.2);
166     fListHist->Add(fHistMassAntiLambdaOff);
167   }
168
169   // V0 on-the-fly distributions
170   if (!fHistV0MultiplicityOn) {
171     if (fCollidingSystems)
172       fHistV0MultiplicityOn = new TH1F("fHistV0MultiplicityOn", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
173     else
174       fHistV0MultiplicityOn = new TH1F("fHistV0MultiplicityOn", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
175     fListHist->Add(fHistV0MultiplicityOn);
176   }
177   if (!fHistV0Chi2On) {
178     fHistV0Chi2On = new TH1F("fHistV0Chi2On", "V0 chi2;chi2;Number of V0s", 33, 0, 33);
179     fListHist->Add(fHistV0Chi2On);
180   }
181   if (!fHistDcaV0DaughtersOn) {
182     fHistDcaV0DaughtersOn = new TH1F("fHistDcaV0DaughtersOn", "DCA between V0 daughters;DCA (cm);Number of V0s", 300, 0, 3);
183     fListHist->Add(fHistDcaV0DaughtersOn);
184   }
185   if (!fHistV0CosineOfPointingAngleOn) {
186     fHistV0CosineOfPointingAngleOn = new TH1F("fHistV0CosineOfPointingAngleOn", "V0 Cosine of Pointing Angle;Number of V0s", 200, 0, 1);
187     fListHist->Add(fHistV0CosineOfPointingAngleOn);
188   }
189   if (!fHistV0RadiusOn) {
190     fHistV0RadiusOn = new TH1F("fHistV0RadiusOn", "V0 decay radius;Radius (cm);Number of V0s", 33, 0, 33);
191     fListHist->Add(fHistV0RadiusOn);
192   }
193   if (!fHistDcaV0ToPrimVertexOn) {
194     fHistDcaV0ToPrimVertexOn = new TH1F("fHistDcaV0ToPrimVertexOn", "DCA of V0 to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
195     fListHist->Add(fHistDcaV0ToPrimVertexOn);
196   }
197   if (!fHistDcaPosToPrimVertexOn) {
198     fHistDcaPosToPrimVertexOn = new TH1F("fHistDcaPosToPrimVertexOn", "DCA of V0 neg daughter to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
199     fListHist->Add(fHistDcaPosToPrimVertexOn);
200   }
201   if (!fHistDcaNegToPrimVertexOn) {
202     fHistDcaNegToPrimVertexOn = new TH1F("fHistDcaNegToPrimVertexOn", "DCA of V0 pos daughter to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
203     fListHist->Add(fHistDcaNegToPrimVertexOn);
204   }
205
206   if (!fHistMassK0On) {
207     fHistMassK0On = new TH1F("fHistMassK0On","K^{0} candidates;M(#pi^{+}#pi^{-}) (GeV/c^{2});Counts",100,0.4,0.6);
208     fListHist->Add(fHistMassK0On);
209   }
210   if (!fHistMassLambdaOn) {
211     fHistMassLambdaOn = new TH1F("fHistMassLambdaOn","#Lambda^{0} candidates;M(p#pi^{-}) (GeV/c^{2});Counts",75,1.05,1.2);
212     fListHist->Add(fHistMassLambdaOn);
213   }
214   if (!fHistMassAntiLambdaOn) {
215     fHistMassAntiLambdaOn = new TH1F("fHistMassAntiLambdaOn","#bar{#Lambda}^{0} candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});Counts",75,1.05,1.2);
216     fListHist->Add(fHistMassAntiLambdaOn);
217   }
218
219 }
220
221 //________________________________________________________________________
222 void AliAnalysisTaskCheckV0::UserExec(Option_t *) 
223 {
224   // Main loop
225   // Called for each event
226   AliVEvent* lEvent = InputEvent();
227   if (!lEvent) {
228     Printf("ERROR: Event not available");
229     return;
230   }
231   fHistTrackMultiplicity->Fill(lEvent->GetNumberOfTracks());
232   
233   Double_t tPrimaryVtxPosition[3];
234   Int_t nv0s = 0;
235   nv0s = lEvent->GetNumberOfV0s();
236   Printf("CheckV0 analysis task: There are %d v0s in this event",nv0s);
237
238   Int_t    lOnFlyStatus = 0, nv0sOn = 0, nv0sOff = 0;
239   Double_t lChi2V0 = 0;
240   Double_t lDcaV0Daughters = 0, lDcaV0ToPrimVertex = 0;
241   Double_t lDcaPosToPrimVertex = 0, lDcaNegToPrimVertex = 0;
242   Double_t lV0CosineOfPointingAngle = 0;
243   Double_t lV0Radius = 0;
244   Double_t lInvMassK0 = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
245
246   if(fAnalysisType == "ESD") {
247
248     const AliESDVertex *primaryVtx = ((AliESDEvent*)lEvent)->GetPrimaryVertex();
249     tPrimaryVtxPosition[0] = primaryVtx->GetXv();
250     tPrimaryVtxPosition[1] = primaryVtx->GetYv();
251     tPrimaryVtxPosition[2] = primaryVtx->GetZv();
252
253     fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
254     fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
255     fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
256
257     for (Int_t iV0 = 0; iV0 < nv0s; iV0++) 
258       {// This is the begining of the V0 loop
259         AliESDv0 *v0 = ((AliESDEvent*)lEvent)->GetV0(iV0);
260         if (!v0) continue;
261
262         Double_t tDecayVertexV0[3]; v0->GetXYZ(tDecayVertexV0[0],tDecayVertexV0[1],tDecayVertexV0[2]); 
263
264         lV0Radius = TMath::Sqrt(tDecayVertexV0[0]*tDecayVertexV0[0]+tDecayVertexV0[1]*tDecayVertexV0[1]);
265
266         UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
267         UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
268
269         Double_t lMomPos[3]; v0->GetPPxPyPz(lMomPos[0],lMomPos[1],lMomPos[2]);
270         Double_t lMomNeg[3]; v0->GetNPxPyPz(lMomNeg[0],lMomNeg[1],lMomNeg[2]);
271
272         AliESDtrack *pTrack=((AliESDEvent*)lEvent)->GetTrack(lKeyPos);
273         AliESDtrack *nTrack=((AliESDEvent*)lEvent)->GetTrack(lKeyNeg);
274         if (!pTrack || !nTrack) {
275           Printf("ERROR: Could not retreive one of the daughter track");
276           continue;
277         }
278
279         Float_t tDcaPosToPrimVertex[2];
280         if(pTrack) pTrack->GetImpactParameters(tDcaPosToPrimVertex[0],tDcaPosToPrimVertex[1]);
281         else { tDcaPosToPrimVertex[0]=999.;  tDcaPosToPrimVertex[1]=999.;}
282         lDcaPosToPrimVertex = TMath::Sqrt(tDcaPosToPrimVertex[0]*tDcaPosToPrimVertex[0]+tDcaPosToPrimVertex[1]*tDcaPosToPrimVertex[1]);
283
284         Float_t tDcaNegToPrimVertex[2];
285         if(nTrack) nTrack->GetImpactParameters(tDcaNegToPrimVertex[0],tDcaNegToPrimVertex[1]);
286         else { tDcaNegToPrimVertex[0]=999.;  tDcaNegToPrimVertex[1]=999.;}
287         lDcaNegToPrimVertex = TMath::Sqrt(tDcaNegToPrimVertex[0]*tDcaNegToPrimVertex[0]+tDcaNegToPrimVertex[1]*tDcaNegToPrimVertex[1]);
288
289
290         lOnFlyStatus = v0->GetOnFlyStatus();
291         lChi2V0 = v0->GetChi2V0();
292         lDcaV0Daughters = v0->GetDcaV0Daughters();
293         lDcaV0ToPrimVertex = v0->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
294         lV0CosineOfPointingAngle = v0->GetV0CosineOfPointingAngle(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
295
296         // Getting invariant mass infos directly from ESD
297         v0->ChangeMassHypothesis(310);
298         lInvMassK0 = v0->GetEffMass();
299         v0->ChangeMassHypothesis(3122);
300         lInvMassLambda = v0->GetEffMass();
301         v0->ChangeMassHypothesis(-3122);
302         lInvMassAntiLambda = v0->GetEffMass();
303
304         fHistV0OnFlyStatus->Fill(lOnFlyStatus);
305         if(!lOnFlyStatus){
306           nv0sOff++;
307           fHistV0Chi2Off->Fill(lChi2V0);
308           fHistDcaV0ToPrimVertexOff->Fill(lDcaV0ToPrimVertex);
309           fHistDcaV0DaughtersOff->Fill(lDcaV0Daughters);
310           fHistV0CosineOfPointingAngleOff->Fill(lV0CosineOfPointingAngle);
311
312           fHistV0RadiusOff->Fill(lV0Radius);
313           fHistDcaPosToPrimVertexOff->Fill(lDcaPosToPrimVertex);
314           fHistDcaNegToPrimVertexOff->Fill(lDcaNegToPrimVertex);
315
316           // Filling invariant mass histos for all candidates
317           fHistMassK0Off->Fill(lInvMassK0);
318           fHistMassLambdaOff->Fill(lInvMassLambda);
319           fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
320         }
321         else {
322           nv0sOn++;
323           fHistV0Chi2On->Fill(lChi2V0);
324           fHistDcaV0ToPrimVertexOn->Fill(lDcaV0ToPrimVertex);
325           fHistDcaV0DaughtersOn->Fill(lDcaV0Daughters);
326           fHistV0CosineOfPointingAngleOn->Fill(lV0CosineOfPointingAngle);
327
328           fHistV0RadiusOn->Fill(lV0Radius);
329           fHistDcaPosToPrimVertexOn->Fill(lDcaPosToPrimVertex);
330           fHistDcaNegToPrimVertexOn->Fill(lDcaNegToPrimVertex);
331
332           // Filling invariant mass histos for all candidates
333           fHistMassK0On->Fill(lInvMassK0);
334           fHistMassLambdaOn->Fill(lInvMassLambda);
335           fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
336         }
337       }// This is the end of the V0 loop
338   } // end of "ESD" analysis
339
340   else if(fAnalysisType == "AOD") {
341
342     const AliAODVertex *primaryVtx = ((AliAODEvent*)lEvent)->GetPrimaryVertex();
343     tPrimaryVtxPosition[0] = primaryVtx->GetX();
344     tPrimaryVtxPosition[1] = primaryVtx->GetY();
345     tPrimaryVtxPosition[2] = primaryVtx->GetZ();
346
347     fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
348     fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
349     fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
350
351     for (Int_t iV0 = 0; iV0 < nv0s; iV0++) 
352       {// This is the begining of the V0 loop
353         AliAODv0 *v0 = ((AliAODEvent*)lEvent)->GetV0(iV0);
354         if (!v0) continue;
355
356         lV0Radius = v0->RadiusV0();
357         lDcaPosToPrimVertex = v0->DcaPosToPrimVertex();
358         lDcaNegToPrimVertex = v0->DcaNegToPrimVertex();
359
360
361         lOnFlyStatus = v0->GetOnFlyStatus();
362         lChi2V0 = v0->Chi2V0();
363         lDcaV0Daughters = v0->DcaV0Daughters();
364         lDcaV0ToPrimVertex = v0->DcaV0ToPrimVertex();
365         lV0CosineOfPointingAngle = v0->CosPointingAngle(tPrimaryVtxPosition);
366
367         lInvMassK0 = v0->MassK0Short();
368         lInvMassLambda = v0->MassLambda();
369         lInvMassAntiLambda = v0->MassAntiLambda();
370
371         fHistV0OnFlyStatus->Fill(lOnFlyStatus);
372         if(!lOnFlyStatus){
373           nv0sOff++;
374           fHistV0Chi2Off->Fill(lChi2V0);
375           fHistDcaV0ToPrimVertexOff->Fill(lDcaV0ToPrimVertex);
376           fHistDcaV0DaughtersOff->Fill(lDcaV0Daughters);
377           fHistV0CosineOfPointingAngleOff->Fill(lV0CosineOfPointingAngle);
378
379           fHistV0RadiusOff->Fill(lV0Radius);
380           fHistDcaPosToPrimVertexOff->Fill(lDcaPosToPrimVertex);
381           fHistDcaNegToPrimVertexOff->Fill(lDcaNegToPrimVertex);
382
383           // Filling invariant mass histos for all candidates
384           fHistMassK0Off->Fill(lInvMassK0);
385           fHistMassLambdaOff->Fill(lInvMassLambda);
386           fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
387         }
388         else {
389           nv0sOn++;
390           fHistV0Chi2On->Fill(lChi2V0);
391           fHistDcaV0ToPrimVertexOn->Fill(lDcaV0ToPrimVertex);
392           fHistDcaV0DaughtersOn->Fill(lDcaV0Daughters);
393           fHistV0CosineOfPointingAngleOn->Fill(lV0CosineOfPointingAngle);
394
395           fHistV0RadiusOn->Fill(lV0Radius);
396           fHistDcaPosToPrimVertexOn->Fill(lDcaPosToPrimVertex);
397           fHistDcaNegToPrimVertexOn->Fill(lDcaNegToPrimVertex);
398
399           // Filling invariant mass histos for all candidates
400           fHistMassK0On->Fill(lInvMassK0);
401           fHistMassLambdaOn->Fill(lInvMassLambda);
402           fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
403         }
404       }// This is the end of the V0 loop
405   } // end of "AOD" analysis
406
407   fHistV0Multiplicity->Fill(nv0s);
408   fHistV0MultiplicityOff->Fill(nv0sOff);
409   fHistV0MultiplicityOn->Fill(nv0sOn);
410
411   // Post output data.
412   PostData(1, fListHist);
413 }      
414
415 //________________________________________________________________________
416 void AliAnalysisTaskCheckV0::Terminate(Option_t *) 
417 {
418   // Draw result to the screen
419   // Called once at the end of the query
420
421   fHistTrackMultiplicity = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistTrackMultiplicity"));
422   if (!fHistTrackMultiplicity) {
423     Printf("ERROR: fHistTrackMultiplicity not available");
424     return;
425   }
426   fHistV0Multiplicity = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistV0Multiplicity"));
427   if (!fHistV0Multiplicity) {
428     Printf("ERROR: fHistV0Multiplicity not available");
429     return;
430   }
431   fHistV0MultiplicityOff = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistV0MultiplicityOff"));
432   if (!fHistV0MultiplicityOff) {
433     Printf("ERROR: fHistV0MultiplicityOff not available");
434     return;
435   }
436   fHistV0MultiplicityOn = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistV0MultiplicityOn"));
437   if (!fHistV0MultiplicityOn) {
438     Printf("ERROR: fHistV0MultiplicityOn not available");
439     return;
440   }
441    
442   TCanvas *canCheckV0 = new TCanvas("AliAnalysisTaskCheckV0","Multiplicity",10,10,510,510);
443   canCheckV0->Divide(2,2);
444   if (fHistTrackMultiplicity->GetMaximum() > 0.) canCheckV0->cd(1)->SetLogy();
445   fHistTrackMultiplicity->SetMarkerStyle(26);
446   fHistTrackMultiplicity->DrawCopy("E");
447   fHistV0Multiplicity->SetMarkerStyle(25);
448   fHistV0Multiplicity->DrawCopy("ESAME");
449   fHistV0MultiplicityOff->SetMarkerStyle(24);
450   fHistV0MultiplicityOff->DrawCopy("ESAME");
451   fHistV0MultiplicityOn->SetMarkerStyle(20);
452   fHistV0MultiplicityOn->DrawCopy("ESAME");
453
454   TLegend *legendMultiplicity = new TLegend(0.4,0.4,0.65,0.6);
455   legendMultiplicity->AddEntry(fHistTrackMultiplicity,"tracks");
456   legendMultiplicity->AddEntry(fHistV0Multiplicity,"all V^{0}");
457   legendMultiplicity->AddEntry(fHistV0MultiplicityOff,"offline V^{0}");
458   legendMultiplicity->AddEntry(fHistV0MultiplicityOn,"onthefly V^{0}");
459   legendMultiplicity->Draw();
460
461   fHistMassK0Off = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassK0Off"));
462   if (!fHistMassK0Off) {
463     Printf("ERROR: fHistMassK0Off not available");
464     return;
465   }
466   fHistMassK0On = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassK0On"));
467   if (!fHistMassK0On) {
468     Printf("ERROR: fHistMassK0On not available");
469     return;
470   }
471   fHistMassLambdaOff = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassLambdaOff"));
472   if (!fHistMassLambdaOff) {
473     Printf("ERROR: fHistMassLambdaOff not available");
474     return;
475   }
476   fHistMassLambdaOn = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassLambdaOn"));
477   if (!fHistMassLambdaOn) {
478     Printf("ERROR: fHistMassLambdaOn not available");
479     return;
480   }
481   fHistMassAntiLambdaOff = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassAntiLambdaOff"));
482   if (!fHistMassAntiLambdaOff) {
483     Printf("ERROR: fHistMassAntiLambdaOff not available");
484     return;
485   }
486   fHistMassAntiLambdaOn = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassAntiLambdaOn"));
487   if (!fHistMassAntiLambdaOn) {
488     Printf("ERROR: fHistMassAntiLambdaOn not available");
489     return;
490   }
491
492   canCheckV0->cd(2);
493   fHistMassK0Off->SetMarkerStyle(24);
494   fHistMassK0Off->DrawCopy("E");
495   fHistMassK0On->SetMarkerStyle(20);
496   fHistMassK0On->DrawCopy("ESAME");
497
498   canCheckV0->cd(3);
499   fHistMassLambdaOff->SetMarkerStyle(24);
500   fHistMassLambdaOff->DrawCopy("E");
501   fHistMassLambdaOn->SetMarkerStyle(20);
502   fHistMassLambdaOn->DrawCopy("ESAME");
503
504   canCheckV0->cd(4);
505   fHistMassAntiLambdaOff->SetMarkerStyle(24);
506   fHistMassAntiLambdaOff->DrawCopy("E");
507   fHistMassAntiLambdaOn->SetMarkerStyle(20);
508   fHistMassAntiLambdaOn->DrawCopy("ESAME");
509
510 }