]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/SPECTRA/AliAnalysisTaskCheckV0.cxx
Compilation with CMake: corrected EINCLUDE
[u/mrichter/AliRoot.git] / PWG2 / SPECTRA / AliAnalysisTaskCheckV0.cxx
CommitLineData
96cad837 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//-----------------------------------------------------------------
c8eada4f 15// AliAnalysisTaskCheckV0 class
16// This task is for QAing the V0s from ESD/AOD
96cad837 17// Origin: B.H. Nov2007, hippolyt@in2p3.fr
18//-----------------------------------------------------------------
96cad837 19#include "TList.h"
20#include "TH1F.h"
21#include "TCanvas.h"
c8eada4f 22#include "TLegend.h"
96cad837 23
c8eada4f 24#include "AliAnalysisTaskSE.h"
96cad837 25
26#include "AliESDEvent.h"
b638d7f5 27#include "AliESDVertex.h"
c8eada4f 28#include "AliAODEvent.h"
96cad837 29
30#include "AliESDv0.h"
31
c8eada4f 32#include "AliAnalysisTaskCheckV0.h"
96cad837 33
c8eada4f 34ClassImp(AliAnalysisTaskCheckV0)
96cad837 35
b638d7f5 36//________________________________________________________________________
37AliAnalysisTaskCheckV0::AliAnalysisTaskCheckV0()
976dc218 38 : AliAnalysisTaskSE(), fAnalysisType("ESD"), fCollidingSystems(0), fListHist(),
b638d7f5 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),
9ea746fc 45 fHistMassK0sOff(0),fHistMassLambdaOff(0),fHistMassAntiLambdaOff(0),
b638d7f5 46 fHistV0MultiplicityOn(0), fHistV0Chi2On(0),
47 fHistDcaV0DaughtersOn(0), fHistV0CosineOfPointingAngleOn(0),
48 fHistV0RadiusOn(0),fHistDcaV0ToPrimVertexOn(0),
49 fHistDcaPosToPrimVertexOn(0),fHistDcaNegToPrimVertexOn(0),
9ea746fc 50 fHistMassK0sOn(0),fHistMassLambdaOn(0),fHistMassAntiLambdaOn(0)
b638d7f5 51{
52 // Dummy constructor
53}
96cad837 54//________________________________________________________________________
c8eada4f 55AliAnalysisTaskCheckV0::AliAnalysisTaskCheckV0(const char *name)
976dc218 56 : AliAnalysisTaskSE(name), fAnalysisType("ESD"), fCollidingSystems(0), fListHist(),
b638d7f5 57 fHistPrimaryVertexPosX(0), fHistPrimaryVertexPosY(0), fHistPrimaryVertexPosZ(0),
96cad837 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),
9ea746fc 63 fHistMassK0sOff(0),fHistMassLambdaOff(0),fHistMassAntiLambdaOff(0),
96cad837 64 fHistV0MultiplicityOn(0), fHistV0Chi2On(0),
65 fHistDcaV0DaughtersOn(0), fHistV0CosineOfPointingAngleOn(0),
66 fHistV0RadiusOn(0),fHistDcaV0ToPrimVertexOn(0),
67 fHistDcaPosToPrimVertexOn(0),fHistDcaNegToPrimVertexOn(0),
9ea746fc 68 fHistMassK0sOn(0),fHistMassLambdaOn(0),fHistMassAntiLambdaOn(0)
96cad837 69{
70 // Constructor
976dc218 71 // Define output slots only here
72 // Output slot #1 writes into a TList container
73 DefineOutput(1, TList::Class());
96cad837 74}
75
76//________________________________________________________________________
976dc218 77void AliAnalysisTaskCheckV0::UserCreateOutputObjects()
96cad837 78{
79 // Create histograms
80 // Called once
81
82 // Distinguish Track and V0 Multiplicity !
83
84 fListHist = new TList();
b638d7f5 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 }
96cad837 98
99 if (!fHistTrackMultiplicity) {
c8eada4f 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);
96cad837 104 fListHist->Add(fHistTrackMultiplicity);
105 }
106 if (!fHistV0Multiplicity) {
c8eada4f 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);
96cad837 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) {
c8eada4f 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);
96cad837 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
9ea746fc 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);
96cad837 158 }
159 if (!fHistMassLambdaOff) {
c8eada4f 160 fHistMassLambdaOff = new TH1F("fHistMassLambdaOff","#Lambda^{0} candidates;M(p#pi^{-}) (GeV/c^{2});Counts",75,1.05,1.2);
96cad837 161 fListHist->Add(fHistMassLambdaOff);
162 }
163 if (!fHistMassAntiLambdaOff) {
c8eada4f 164 fHistMassAntiLambdaOff = new TH1F("fHistMassAntiLambdaOff","#bar{#Lambda}^{0} candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});Counts",75,1.05,1.2);
96cad837 165 fListHist->Add(fHistMassAntiLambdaOff);
166 }
167
168 // V0 on-the-fly distributions
169 if (!fHistV0MultiplicityOn) {
c8eada4f 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);
96cad837 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
9ea746fc 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);
96cad837 208 }
209 if (!fHistMassLambdaOn) {
c8eada4f 210 fHistMassLambdaOn = new TH1F("fHistMassLambdaOn","#Lambda^{0} candidates;M(p#pi^{-}) (GeV/c^{2});Counts",75,1.05,1.2);
96cad837 211 fListHist->Add(fHistMassLambdaOn);
212 }
213 if (!fHistMassAntiLambdaOn) {
c8eada4f 214 fHistMassAntiLambdaOn = new TH1F("fHistMassAntiLambdaOn","#bar{#Lambda}^{0} candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});Counts",75,1.05,1.2);
96cad837 215 fListHist->Add(fHistMassAntiLambdaOn);
216 }
217
218}
219
220//________________________________________________________________________
976dc218 221void AliAnalysisTaskCheckV0::UserExec(Option_t *)
96cad837 222{
223 // Main loop
224 // Called for each event
976dc218 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();
83db4a5d 235 // Printf("CheckV0 analysis task: There are %d v0s in this event",nv0s);
976dc218 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;
9ea746fc 243 Double_t lInvMassK0s = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
96cad837 244
c8eada4f 245 if(fAnalysisType == "ESD") {
246
976dc218 247 const AliESDVertex *primaryVtx = ((AliESDEvent*)lEvent)->GetPrimaryVertex();
b638d7f5 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
bd7ae1e4 256 Double_t lMagneticField = ((AliESDEvent*)lEvent)->GetMagneticField();
257
c8eada4f 258 for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
259 {// This is the begining of the V0 loop
976dc218 260 AliESDv0 *v0 = ((AliESDEvent*)lEvent)->GetV0(iV0);
c8eada4f 261 if (!v0) continue;
262
263 Double_t tDecayVertexV0[3]; v0->GetXYZ(tDecayVertexV0[0],tDecayVertexV0[1],tDecayVertexV0[2]);
264
265 lV0Radius = TMath::Sqrt(tDecayVertexV0[0]*tDecayVertexV0[0]+tDecayVertexV0[1]*tDecayVertexV0[1]);
266
267 UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
268 UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
269
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]);
272
976dc218 273 AliESDtrack *pTrack=((AliESDEvent*)lEvent)->GetTrack(lKeyPos);
274 AliESDtrack *nTrack=((AliESDEvent*)lEvent)->GetTrack(lKeyNeg);
c8eada4f 275 if (!pTrack || !nTrack) {
276 Printf("ERROR: Could not retreive one of the daughter track");
277 continue;
278 }
279
bd7ae1e4 280 lDcaPosToPrimVertex = TMath::Abs(pTrack->GetD(tPrimaryVtxPosition[0],
281 tPrimaryVtxPosition[1],
282 lMagneticField) );
c8eada4f 283
bd7ae1e4 284 lDcaNegToPrimVertex = TMath::Abs(nTrack->GetD(tPrimaryVtxPosition[0],
285 tPrimaryVtxPosition[1],
286 lMagneticField) );
c8eada4f 287
288 lOnFlyStatus = v0->GetOnFlyStatus();
289 lChi2V0 = v0->GetChi2V0();
290 lDcaV0Daughters = v0->GetDcaV0Daughters();
60c8c322 291 lDcaV0ToPrimVertex = v0->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
292 lV0CosineOfPointingAngle = v0->GetV0CosineOfPointingAngle(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
c8eada4f 293
294 // Getting invariant mass infos directly from ESD
295 v0->ChangeMassHypothesis(310);
9ea746fc 296 lInvMassK0s = v0->GetEffMass();
c8eada4f 297 v0->ChangeMassHypothesis(3122);
298 lInvMassLambda = v0->GetEffMass();
299 v0->ChangeMassHypothesis(-3122);
300 lInvMassAntiLambda = v0->GetEffMass();
301
302 fHistV0OnFlyStatus->Fill(lOnFlyStatus);
303 if(!lOnFlyStatus){
304 nv0sOff++;
305 fHistV0Chi2Off->Fill(lChi2V0);
306 fHistDcaV0ToPrimVertexOff->Fill(lDcaV0ToPrimVertex);
307 fHistDcaV0DaughtersOff->Fill(lDcaV0Daughters);
308 fHistV0CosineOfPointingAngleOff->Fill(lV0CosineOfPointingAngle);
309
310 fHistV0RadiusOff->Fill(lV0Radius);
311 fHistDcaPosToPrimVertexOff->Fill(lDcaPosToPrimVertex);
312 fHistDcaNegToPrimVertexOff->Fill(lDcaNegToPrimVertex);
313
314 // Filling invariant mass histos for all candidates
9ea746fc 315 fHistMassK0sOff->Fill(lInvMassK0s);
c8eada4f 316 fHistMassLambdaOff->Fill(lInvMassLambda);
317 fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
318 }
319 else {
320 nv0sOn++;
321 fHistV0Chi2On->Fill(lChi2V0);
322 fHistDcaV0ToPrimVertexOn->Fill(lDcaV0ToPrimVertex);
323 fHistDcaV0DaughtersOn->Fill(lDcaV0Daughters);
324 fHistV0CosineOfPointingAngleOn->Fill(lV0CosineOfPointingAngle);
325
326 fHistV0RadiusOn->Fill(lV0Radius);
327 fHistDcaPosToPrimVertexOn->Fill(lDcaPosToPrimVertex);
328 fHistDcaNegToPrimVertexOn->Fill(lDcaNegToPrimVertex);
329
330 // Filling invariant mass histos for all candidates
9ea746fc 331 fHistMassK0sOn->Fill(lInvMassK0s);
c8eada4f 332 fHistMassLambdaOn->Fill(lInvMassLambda);
333 fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
334 }
335 }// This is the end of the V0 loop
c8eada4f 336 } // end of "ESD" analysis
337
338 else if(fAnalysisType == "AOD") {
339
976dc218 340 const AliAODVertex *primaryVtx = ((AliAODEvent*)lEvent)->GetPrimaryVertex();
b638d7f5 341 tPrimaryVtxPosition[0] = primaryVtx->GetX();
342 tPrimaryVtxPosition[1] = primaryVtx->GetY();
343 tPrimaryVtxPosition[2] = primaryVtx->GetZ();
344
345 fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
346 fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
347 fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
348
c8eada4f 349 for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
350 {// This is the begining of the V0 loop
976dc218 351 AliAODv0 *v0 = ((AliAODEvent*)lEvent)->GetV0(iV0);
c8eada4f 352 if (!v0) continue;
353
354 lV0Radius = v0->RadiusV0();
355 lDcaPosToPrimVertex = v0->DcaPosToPrimVertex();
356 lDcaNegToPrimVertex = v0->DcaNegToPrimVertex();
357
358
b638d7f5 359 lOnFlyStatus = v0->GetOnFlyStatus();
c8eada4f 360 lChi2V0 = v0->Chi2V0();
361 lDcaV0Daughters = v0->DcaV0Daughters();
362 lDcaV0ToPrimVertex = v0->DcaV0ToPrimVertex();
b638d7f5 363 lV0CosineOfPointingAngle = v0->CosPointingAngle(tPrimaryVtxPosition);
c8eada4f 364
9ea746fc 365 lInvMassK0s = v0->MassK0Short();
c8eada4f 366 lInvMassLambda = v0->MassLambda();
367 lInvMassAntiLambda = v0->MassAntiLambda();
368
369 fHistV0OnFlyStatus->Fill(lOnFlyStatus);
370 if(!lOnFlyStatus){
371 nv0sOff++;
372 fHistV0Chi2Off->Fill(lChi2V0);
373 fHistDcaV0ToPrimVertexOff->Fill(lDcaV0ToPrimVertex);
374 fHistDcaV0DaughtersOff->Fill(lDcaV0Daughters);
375 fHistV0CosineOfPointingAngleOff->Fill(lV0CosineOfPointingAngle);
376
377 fHistV0RadiusOff->Fill(lV0Radius);
378 fHistDcaPosToPrimVertexOff->Fill(lDcaPosToPrimVertex);
379 fHistDcaNegToPrimVertexOff->Fill(lDcaNegToPrimVertex);
380
381 // Filling invariant mass histos for all candidates
9ea746fc 382 fHistMassK0sOff->Fill(lInvMassK0s);
c8eada4f 383 fHistMassLambdaOff->Fill(lInvMassLambda);
384 fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
385 }
386 else {
387 nv0sOn++;
388 fHistV0Chi2On->Fill(lChi2V0);
389 fHistDcaV0ToPrimVertexOn->Fill(lDcaV0ToPrimVertex);
390 fHistDcaV0DaughtersOn->Fill(lDcaV0Daughters);
391 fHistV0CosineOfPointingAngleOn->Fill(lV0CosineOfPointingAngle);
392
393 fHistV0RadiusOn->Fill(lV0Radius);
394 fHistDcaPosToPrimVertexOn->Fill(lDcaPosToPrimVertex);
395 fHistDcaNegToPrimVertexOn->Fill(lDcaNegToPrimVertex);
396
397 // Filling invariant mass histos for all candidates
9ea746fc 398 fHistMassK0sOn->Fill(lInvMassK0s);
c8eada4f 399 fHistMassLambdaOn->Fill(lInvMassLambda);
400 fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
401 }
402 }// This is the end of the V0 loop
c8eada4f 403 } // end of "AOD" analysis
96cad837 404
976dc218 405 fHistV0Multiplicity->Fill(nv0s);
406 fHistV0MultiplicityOff->Fill(nv0sOff);
407 fHistV0MultiplicityOn->Fill(nv0sOn);
408
96cad837 409 // Post output data.
976dc218 410 PostData(1, fListHist);
96cad837 411}
412
413//________________________________________________________________________
c8eada4f 414void AliAnalysisTaskCheckV0::Terminate(Option_t *)
96cad837 415{
416 // Draw result to the screen
417 // Called once at the end of the query
418
976dc218 419 fHistTrackMultiplicity = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistTrackMultiplicity"));
96cad837 420 if (!fHistTrackMultiplicity) {
421 Printf("ERROR: fHistTrackMultiplicity not available");
422 return;
423 }
976dc218 424 fHistV0Multiplicity = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistV0Multiplicity"));
96cad837 425 if (!fHistV0Multiplicity) {
426 Printf("ERROR: fHistV0Multiplicity not available");
427 return;
428 }
976dc218 429 fHistV0MultiplicityOff = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistV0MultiplicityOff"));
96cad837 430 if (!fHistV0MultiplicityOff) {
431 Printf("ERROR: fHistV0MultiplicityOff not available");
432 return;
433 }
976dc218 434 fHistV0MultiplicityOn = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistV0MultiplicityOn"));
96cad837 435 if (!fHistV0MultiplicityOn) {
436 Printf("ERROR: fHistV0MultiplicityOn not available");
437 return;
438 }
439
c8eada4f 440 TCanvas *canCheckV0 = new TCanvas("AliAnalysisTaskCheckV0","Multiplicity",10,10,510,510);
441 canCheckV0->Divide(2,2);
976dc218 442 if (fHistTrackMultiplicity->GetMaximum() > 0.) canCheckV0->cd(1)->SetLogy();
b638d7f5 443 fHistTrackMultiplicity->SetMarkerStyle(26);
96cad837 444 fHistTrackMultiplicity->DrawCopy("E");
b638d7f5 445 fHistV0Multiplicity->SetMarkerStyle(25);
96cad837 446 fHistV0Multiplicity->DrawCopy("ESAME");
447 fHistV0MultiplicityOff->SetMarkerStyle(24);
448 fHistV0MultiplicityOff->DrawCopy("ESAME");
449 fHistV0MultiplicityOn->SetMarkerStyle(20);
450 fHistV0MultiplicityOn->DrawCopy("ESAME");
c8eada4f 451
a95976c8 452 TLegend *legendMultiplicity = new TLegend(0.5,0.5,0.75,0.75);
c8eada4f 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();
458
9ea746fc 459 fHistMassK0sOff = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassK0sOff"));
460 if (!fHistMassK0sOff) {
461 Printf("ERROR: fHistMassK0sOff not available");
c8eada4f 462 return;
463 }
9ea746fc 464 fHistMassK0sOn = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassK0sOn"));
465 if (!fHistMassK0sOn) {
466 Printf("ERROR: fHistMassK0sOn not available");
c8eada4f 467 return;
468 }
976dc218 469 fHistMassLambdaOff = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassLambdaOff"));
c8eada4f 470 if (!fHistMassLambdaOff) {
471 Printf("ERROR: fHistMassLambdaOff not available");
472 return;
473 }
976dc218 474 fHistMassLambdaOn = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassLambdaOn"));
c8eada4f 475 if (!fHistMassLambdaOn) {
476 Printf("ERROR: fHistMassLambdaOn not available");
477 return;
478 }
976dc218 479 fHistMassAntiLambdaOff = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassAntiLambdaOff"));
c8eada4f 480 if (!fHistMassAntiLambdaOff) {
481 Printf("ERROR: fHistMassAntiLambdaOff not available");
482 return;
483 }
976dc218 484 fHistMassAntiLambdaOn = dynamic_cast<TH1F*> (((TList*)GetOutputData(1))->FindObject("fHistMassAntiLambdaOn"));
c8eada4f 485 if (!fHistMassAntiLambdaOn) {
486 Printf("ERROR: fHistMassAntiLambdaOn not available");
487 return;
488 }
489
490 canCheckV0->cd(2);
9ea746fc 491 fHistMassK0sOn->SetMarkerStyle(20);
492 fHistMassK0sOn->DrawCopy("E");
493 fHistMassK0sOff->SetMarkerStyle(24);
494 fHistMassK0sOff->DrawCopy("ESAME");
c8eada4f 495
496 canCheckV0->cd(3);
c8eada4f 497 fHistMassLambdaOn->SetMarkerStyle(20);
a95976c8 498 fHistMassLambdaOn->DrawCopy("E");
499 fHistMassLambdaOff->SetMarkerStyle(24);
500 fHistMassLambdaOff->DrawCopy("ESAME");
c8eada4f 501
502 canCheckV0->cd(4);
c8eada4f 503 fHistMassAntiLambdaOn->SetMarkerStyle(20);
a95976c8 504 fHistMassAntiLambdaOn->DrawCopy("E");
505 fHistMassAntiLambdaOff->SetMarkerStyle(24);
506 fHistMassAntiLambdaOff->DrawCopy("ESAME");
976dc218 507
96cad837 508}