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