1 /**************************************************************************
2 * Author: Boris Hippolyte. *
3 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
14 //-----------------------------------------------------------------
15 // AliAnalysisTaskCheckV0 class
16 // This task is for QAing the V0s from ESD/AOD
17 // Origin: B.H. Nov2007, hippolyt@in2p3.fr
18 //-----------------------------------------------------------------
24 #include "AliAnalysisTaskSE.h"
26 #include "AliESDEvent.h"
27 #include "AliESDVertex.h"
28 #include "AliAODEvent.h"
32 #include "AliAnalysisTaskCheckV0.h"
34 ClassImp(AliAnalysisTaskCheckV0)
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)
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)
71 // Define output slots only here
72 // Output slot #1 writes into a TList container
73 DefineOutput(1, TList::Class());
76 //________________________________________________________________________
77 void AliAnalysisTaskCheckV0::UserCreateOutputObjects()
82 // Distinguish Track and V0 Multiplicity !
84 fListHist = new TList();
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);
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);
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);
99 if (!fHistTrackMultiplicity) {
100 if (fCollidingSystems)
101 fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 200, 0, 40000);
103 fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 250, 0, 250);
104 fListHist->Add(fHistTrackMultiplicity);
106 if (!fHistV0Multiplicity) {
107 if (fCollidingSystems)
108 fHistV0Multiplicity = new TH1F("fHistV0Multiplicity", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
110 fHistV0Multiplicity = new TH1F("fHistV0Multiplicity", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
111 fListHist->Add(fHistV0Multiplicity);
113 if (!fHistV0OnFlyStatus) {
114 fHistV0OnFlyStatus = new TH1F("fHistV0OnFlyStatus", "V0 On fly status;status;Number of V0s", 3, 0, 3);
115 fListHist->Add(fHistV0OnFlyStatus);
118 // V0 offline distributions
119 if (!fHistV0MultiplicityOff) {
120 if (fCollidingSystems)
121 fHistV0MultiplicityOff = new TH1F("fHistV0MultiplicityOff", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
123 fHistV0MultiplicityOff = new TH1F("fHistV0MultiplicityOff", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
124 fListHist->Add(fHistV0MultiplicityOff);
126 if (!fHistV0Chi2Off) {
127 fHistV0Chi2Off = new TH1F("fHistV0Chi2Off", "V0 chi2;chi2;Number of V0s", 33, 0, 33);
128 fListHist->Add(fHistV0Chi2Off);
130 if (!fHistDcaV0DaughtersOff) {
131 fHistDcaV0DaughtersOff = new TH1F("fHistDcaV0DaughtersOff", "DCA between V0 daughters;DCA (cm);Number of V0s", 300, 0, 3);
132 fListHist->Add(fHistDcaV0DaughtersOff);
134 if (!fHistV0CosineOfPointingAngleOff) {
135 fHistV0CosineOfPointingAngleOff = new TH1F("fHistV0CosineOfPointingAngleOff", "V0 Cosine of Pointing Angle;Number of V0s", 200, 0, 1);
136 fListHist->Add(fHistV0CosineOfPointingAngleOff);
138 if (!fHistV0RadiusOff) {
139 fHistV0RadiusOff = new TH1F("fHistV0RadiusOff", "V0 decay radius;Radius (cm);Number of V0s", 33, 0, 33);
140 fListHist->Add(fHistV0RadiusOff);
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);
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);
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);
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);
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);
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);
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);
173 fHistV0MultiplicityOn = new TH1F("fHistV0MultiplicityOn", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
174 fListHist->Add(fHistV0MultiplicityOn);
176 if (!fHistV0Chi2On) {
177 fHistV0Chi2On = new TH1F("fHistV0Chi2On", "V0 chi2;chi2;Number of V0s", 33, 0, 33);
178 fListHist->Add(fHistV0Chi2On);
180 if (!fHistDcaV0DaughtersOn) {
181 fHistDcaV0DaughtersOn = new TH1F("fHistDcaV0DaughtersOn", "DCA between V0 daughters;DCA (cm);Number of V0s", 300, 0, 3);
182 fListHist->Add(fHistDcaV0DaughtersOn);
184 if (!fHistV0CosineOfPointingAngleOn) {
185 fHistV0CosineOfPointingAngleOn = new TH1F("fHistV0CosineOfPointingAngleOn", "V0 Cosine of Pointing Angle;Number of V0s", 200, 0, 1);
186 fListHist->Add(fHistV0CosineOfPointingAngleOn);
188 if (!fHistV0RadiusOn) {
189 fHistV0RadiusOn = new TH1F("fHistV0RadiusOn", "V0 decay radius;Radius (cm);Number of V0s", 33, 0, 33);
190 fListHist->Add(fHistV0RadiusOn);
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);
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);
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);
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);
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);
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);
220 //________________________________________________________________________
221 void AliAnalysisTaskCheckV0::UserExec(Option_t *)
224 // Called for each event
225 AliVEvent* lEvent = InputEvent();
227 Printf("ERROR: Event not available");
230 fHistTrackMultiplicity->Fill(lEvent->GetNumberOfTracks());
232 Double_t tPrimaryVtxPosition[3];
234 nv0s = lEvent->GetNumberOfV0s();
235 // Printf("CheckV0 analysis task: There are %d v0s in this event",nv0s);
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;
245 if(fAnalysisType == "ESD") {
247 const AliESDVertex *primaryVtx = ((AliESDEvent*)lEvent)->GetPrimaryVertex();
248 tPrimaryVtxPosition[0] = primaryVtx->GetXv();
249 tPrimaryVtxPosition[1] = primaryVtx->GetYv();
250 tPrimaryVtxPosition[2] = primaryVtx->GetZv();
252 fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
253 fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
254 fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
256 Double_t lMagneticField = ((AliESDEvent*)lEvent)->GetMagneticField();
258 for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
259 {// This is the begining of the V0 loop
260 AliESDv0 *v0 = ((AliESDEvent*)lEvent)->GetV0(iV0);
263 Double_t tDecayVertexV0[3]; v0->GetXYZ(tDecayVertexV0[0],tDecayVertexV0[1],tDecayVertexV0[2]);
265 lV0Radius = TMath::Sqrt(tDecayVertexV0[0]*tDecayVertexV0[0]+tDecayVertexV0[1]*tDecayVertexV0[1]);
267 UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
268 UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
270 Double_t lMomPos[3]; v0->GetPPxPyPz(lMomPos[0],lMomPos[1],lMomPos[2]);
271 Double_t lMomNeg[3]; v0->GetNPxPyPz(lMomNeg[0],lMomNeg[1],lMomNeg[2]);
273 AliESDtrack *pTrack=((AliESDEvent*)lEvent)->GetTrack(lKeyPos);
274 AliESDtrack *nTrack=((AliESDEvent*)lEvent)->GetTrack(lKeyNeg);
275 if (!pTrack || !nTrack) {
276 Printf("ERROR: Could not retreive one of the daughter track");
280 lDcaPosToPrimVertex = TMath::Abs(pTrack->GetD(tPrimaryVtxPosition[0],
281 tPrimaryVtxPosition[1],
284 lDcaNegToPrimVertex = TMath::Abs(nTrack->GetD(tPrimaryVtxPosition[0],
285 tPrimaryVtxPosition[1],
288 lOnFlyStatus = v0->GetOnFlyStatus();
289 lChi2V0 = v0->GetChi2V0();
290 lDcaV0Daughters = v0->GetDcaV0Daughters();
291 lDcaV0ToPrimVertex = v0->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
292 lV0CosineOfPointingAngle = v0->GetV0CosineOfPointingAngle(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
294 // Getting invariant mass infos directly from ESD
295 v0->ChangeMassHypothesis(310);
296 lInvMassK0s = v0->GetEffMass();
297 v0->ChangeMassHypothesis(3122);
298 lInvMassLambda = v0->GetEffMass();
299 v0->ChangeMassHypothesis(-3122);
300 lInvMassAntiLambda = v0->GetEffMass();
302 fHistV0OnFlyStatus->Fill(lOnFlyStatus);
305 fHistV0Chi2Off->Fill(lChi2V0);
306 fHistDcaV0ToPrimVertexOff->Fill(lDcaV0ToPrimVertex);
307 fHistDcaV0DaughtersOff->Fill(lDcaV0Daughters);
308 fHistV0CosineOfPointingAngleOff->Fill(lV0CosineOfPointingAngle);
310 fHistV0RadiusOff->Fill(lV0Radius);
311 fHistDcaPosToPrimVertexOff->Fill(lDcaPosToPrimVertex);
312 fHistDcaNegToPrimVertexOff->Fill(lDcaNegToPrimVertex);
314 // Filling invariant mass histos for all candidates
315 fHistMassK0sOff->Fill(lInvMassK0s);
316 fHistMassLambdaOff->Fill(lInvMassLambda);
317 fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
321 fHistV0Chi2On->Fill(lChi2V0);
322 fHistDcaV0ToPrimVertexOn->Fill(lDcaV0ToPrimVertex);
323 fHistDcaV0DaughtersOn->Fill(lDcaV0Daughters);
324 fHistV0CosineOfPointingAngleOn->Fill(lV0CosineOfPointingAngle);
326 fHistV0RadiusOn->Fill(lV0Radius);
327 fHistDcaPosToPrimVertexOn->Fill(lDcaPosToPrimVertex);
328 fHistDcaNegToPrimVertexOn->Fill(lDcaNegToPrimVertex);
330 // Filling invariant mass histos for all candidates
331 fHistMassK0sOn->Fill(lInvMassK0s);
332 fHistMassLambdaOn->Fill(lInvMassLambda);
333 fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
335 }// This is the end of the V0 loop
336 } // end of "ESD" analysis
338 else if(fAnalysisType == "AOD") {
340 const AliAODVertex *primaryVtx = ((AliAODEvent*)lEvent)->GetPrimaryVertex();
341 tPrimaryVtxPosition[0] = primaryVtx->GetX();
342 tPrimaryVtxPosition[1] = primaryVtx->GetY();
343 tPrimaryVtxPosition[2] = primaryVtx->GetZ();
345 fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
346 fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
347 fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
349 for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
350 {// This is the begining of the V0 loop
351 AliAODv0 *v0 = ((AliAODEvent*)lEvent)->GetV0(iV0);
354 lV0Radius = v0->RadiusV0();
355 lDcaPosToPrimVertex = v0->DcaPosToPrimVertex();
356 lDcaNegToPrimVertex = v0->DcaNegToPrimVertex();
359 lOnFlyStatus = v0->GetOnFlyStatus();
360 lChi2V0 = v0->Chi2V0();
361 lDcaV0Daughters = v0->DcaV0Daughters();
362 lDcaV0ToPrimVertex = v0->DcaV0ToPrimVertex();
363 lV0CosineOfPointingAngle = v0->CosPointingAngle(tPrimaryVtxPosition);
365 lInvMassK0s = v0->MassK0Short();
366 lInvMassLambda = v0->MassLambda();
367 lInvMassAntiLambda = v0->MassAntiLambda();
369 fHistV0OnFlyStatus->Fill(lOnFlyStatus);
372 fHistV0Chi2Off->Fill(lChi2V0);
373 fHistDcaV0ToPrimVertexOff->Fill(lDcaV0ToPrimVertex);
374 fHistDcaV0DaughtersOff->Fill(lDcaV0Daughters);
375 fHistV0CosineOfPointingAngleOff->Fill(lV0CosineOfPointingAngle);
377 fHistV0RadiusOff->Fill(lV0Radius);
378 fHistDcaPosToPrimVertexOff->Fill(lDcaPosToPrimVertex);
379 fHistDcaNegToPrimVertexOff->Fill(lDcaNegToPrimVertex);
381 // Filling invariant mass histos for all candidates
382 fHistMassK0sOff->Fill(lInvMassK0s);
383 fHistMassLambdaOff->Fill(lInvMassLambda);
384 fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
388 fHistV0Chi2On->Fill(lChi2V0);
389 fHistDcaV0ToPrimVertexOn->Fill(lDcaV0ToPrimVertex);
390 fHistDcaV0DaughtersOn->Fill(lDcaV0Daughters);
391 fHistV0CosineOfPointingAngleOn->Fill(lV0CosineOfPointingAngle);
393 fHistV0RadiusOn->Fill(lV0Radius);
394 fHistDcaPosToPrimVertexOn->Fill(lDcaPosToPrimVertex);
395 fHistDcaNegToPrimVertexOn->Fill(lDcaNegToPrimVertex);
397 // Filling invariant mass histos for all candidates
398 fHistMassK0sOn->Fill(lInvMassK0s);
399 fHistMassLambdaOn->Fill(lInvMassLambda);
400 fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
402 }// This is the end of the V0 loop
403 } // end of "AOD" analysis
405 fHistV0Multiplicity->Fill(nv0s);
406 fHistV0MultiplicityOff->Fill(nv0sOff);
407 fHistV0MultiplicityOn->Fill(nv0sOn);
410 PostData(1, fListHist);
413 //________________________________________________________________________
414 void AliAnalysisTaskCheckV0::Terminate(Option_t *)
416 // Draw result to the screen
417 // Called once at the end of the query
419 fHistTrackMultiplicity = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistTrackMultiplicity"));
420 if (!fHistTrackMultiplicity) {
421 Printf("ERROR: fHistTrackMultiplicity not available");
424 fHistV0Multiplicity = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistV0Multiplicity"));
425 if (!fHistV0Multiplicity) {
426 Printf("ERROR: fHistV0Multiplicity not available");
429 fHistV0MultiplicityOff = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistV0MultiplicityOff"));
430 if (!fHistV0MultiplicityOff) {
431 Printf("ERROR: fHistV0MultiplicityOff not available");
434 fHistV0MultiplicityOn = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistV0MultiplicityOn"));
435 if (!fHistV0MultiplicityOn) {
436 Printf("ERROR: fHistV0MultiplicityOn not available");
440 TCanvas *canCheckV0 = new TCanvas("AliAnalysisTaskCheckV0","Multiplicity",10,10,510,510);
441 canCheckV0->Divide(2,2);
442 if (fHistTrackMultiplicity->GetMaximum() > 0.) canCheckV0->cd(1)->SetLogy();
443 fHistTrackMultiplicity->SetMarkerStyle(26);
444 fHistTrackMultiplicity->DrawCopy("E");
445 fHistV0Multiplicity->SetMarkerStyle(25);
446 fHistV0Multiplicity->DrawCopy("ESAME");
447 fHistV0MultiplicityOff->SetMarkerStyle(24);
448 fHistV0MultiplicityOff->DrawCopy("ESAME");
449 fHistV0MultiplicityOn->SetMarkerStyle(20);
450 fHistV0MultiplicityOn->DrawCopy("ESAME");
452 TLegend *legendMultiplicity = new TLegend(0.5,0.5,0.75,0.75);
453 legendMultiplicity->AddEntry(fHistTrackMultiplicity,"tracks");
454 legendMultiplicity->AddEntry(fHistV0Multiplicity,"all V^{0}");
455 legendMultiplicity->AddEntry(fHistV0MultiplicityOff,"offline V^{0}");
456 legendMultiplicity->AddEntry(fHistV0MultiplicityOn,"onthefly V^{0}");
457 legendMultiplicity->Draw();
459 fHistMassK0sOff = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassK0sOff"));
460 if (!fHistMassK0sOff) {
461 Printf("ERROR: fHistMassK0sOff not available");
464 fHistMassK0sOn = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassK0sOn"));
465 if (!fHistMassK0sOn) {
466 Printf("ERROR: fHistMassK0sOn not available");
469 fHistMassLambdaOff = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassLambdaOff"));
470 if (!fHistMassLambdaOff) {
471 Printf("ERROR: fHistMassLambdaOff not available");
474 fHistMassLambdaOn = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassLambdaOn"));
475 if (!fHistMassLambdaOn) {
476 Printf("ERROR: fHistMassLambdaOn not available");
479 fHistMassAntiLambdaOff = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassAntiLambdaOff"));
480 if (!fHistMassAntiLambdaOff) {
481 Printf("ERROR: fHistMassAntiLambdaOff not available");
484 fHistMassAntiLambdaOn = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassAntiLambdaOn"));
485 if (!fHistMassAntiLambdaOn) {
486 Printf("ERROR: fHistMassAntiLambdaOn not available");
491 fHistMassK0sOn->SetMarkerStyle(20);
492 fHistMassK0sOn->DrawCopy("E");
493 fHistMassK0sOff->SetMarkerStyle(24);
494 fHistMassK0sOff->DrawCopy("ESAME");
497 fHistMassLambdaOn->SetMarkerStyle(20);
498 fHistMassLambdaOn->DrawCopy("E");
499 fHistMassLambdaOff->SetMarkerStyle(24);
500 fHistMassLambdaOff->DrawCopy("ESAME");
503 fHistMassAntiLambdaOn->SetMarkerStyle(20);
504 fHistMassAntiLambdaOn->DrawCopy("E");
505 fHistMassAntiLambdaOff->SetMarkerStyle(24);
506 fHistMassAntiLambdaOff->DrawCopy("ESAME");