1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 // +-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
18 // Modified version of AliAnalysisTaskCheckCascade.cxx.
19 // This is a 'hybrid' output version, in that it uses a classic TTree
20 // ROOT object to store the candidates, plus a couple of histograms filled on
21 // a per-event basis for storing variables too numerous to put in a tree.
24 // --- Algorithm Description
25 // 1. Loop over primaries in stack to acquire generated charged Xi
26 // 2. Loop over stack to find Cascades, fill TH3Fs "PrimRawPt"s for Efficiency
27 // 3. Perform Physics Selection
28 // 4. Perform Primary Vertex |z|<10cm selection
29 // 5. Perform Primary Vertex NoTPCOnly vertexing selection
30 // 6. Perform Pileup Rejection
32 // 7a. Fill TH3Fs "PrimAnalysisPt" for control purposes only
34 // Please Report Any Bugs!
36 // --- David Dobrigkeit Chinellato
37 // (david.chinellato@gmail.com)
39 // +-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
45 //class AliMCEventHandler;
54 #include <Riostream.h>
60 #include "THnSparse.h"
67 #include "AliESDEvent.h"
68 #include "AliAODEvent.h"
69 #include "AliV0vertexer.h"
70 #include "AliCascadeVertexer.h"
71 #include "AliESDpid.h"
72 #include "AliESDtrack.h"
73 #include "AliESDtrackCuts.h"
74 #include "AliInputEventHandler.h"
75 #include "AliAnalysisManager.h"
76 #include "AliMCEventHandler.h"
77 #include "AliMCEvent.h"
80 #include "AliV0vertexer.h"
81 #include "AliCascadeVertexer.h"
83 #include "AliCFContainer.h"
84 #include "AliMultiplicity.h"
85 #include "AliAODMCParticle.h"
86 #include "AliESDcascade.h"
87 #include "AliAODcascade.h"
88 #include "AliESDUtils.h"
89 #include "AliGenEventHeader.h"
90 #include "AliAnalysisUtils.h"
92 #include "AliAnalysisTaskExtractCascade.h"
97 ClassImp(AliAnalysisTaskExtractCascade)
99 AliAnalysisTaskExtractCascade::AliAnalysisTaskExtractCascade()
100 : AliAnalysisTaskSE(), fListHist(0), fTreeCascade(0), fPIDResponse(0), fESDtrackCuts(0), fUtils(0),
101 fkIsNuclear ( kFALSE ),
102 fkSwitchINT7 ( kFALSE ),
103 fCentralityEstimator("V0M"),
104 fkpAVertexSelection( kFALSE ),
106 fkRunVertexers ( kFALSE ),
107 //------------------------------------------------
109 //------------------------------------------------
111 fTreeCascVarCharge(0),
112 fTreeCascVarMassAsXi(0),
113 fTreeCascVarMassAsOmega(0),
116 fTreeCascVarRapMC(0),
117 fTreeCascVarRapXi(0),
118 fTreeCascVarRapOmega(0),
119 fTreeCascVarNegEta(0),
120 fTreeCascVarPosEta(0),
121 fTreeCascVarBachEta(0),
122 fTreeCascVarDCACascDaughters(0),
123 fTreeCascVarDCABachToPrimVtx(0),
124 fTreeCascVarDCAV0Daughters(0),
125 fTreeCascVarDCAV0ToPrimVtx(0),
126 fTreeCascVarDCAPosToPrimVtx(0),
127 fTreeCascVarDCANegToPrimVtx(0),
128 fTreeCascVarCascCosPointingAngle(0),
129 fTreeCascVarCascRadius(0),
130 fTreeCascVarV0Mass(0),
131 fTreeCascVarV0CosPointingAngle(0),
132 fTreeCascVarV0CosPointingAngleSpecial(0),
133 fTreeCascVarV0Radius(0),
134 fTreeCascVarLeastNbrClusters(0),
135 fTreeCascVarMultiplicity(0),
136 fTreeCascVarMultiplicityV0A(0),
137 fTreeCascVarMultiplicityZNA(0),
138 fTreeCascVarMultiplicityTRK(0),
139 fTreeCascVarMultiplicitySPD(0),
140 fTreeCascVarDistOverTotMom(0),
142 fTreeCascVarPIDBachelor(0),
143 fTreeCascVarPIDNegative(0),
144 fTreeCascVarPIDPositive(0),
145 fTreeCascVarBachTransMom(0),
146 fTreeCascVarPosTransMom(0),
147 fTreeCascVarNegTransMom(0),
148 fTreeCascVarPosTransMomMC(0),
149 fTreeCascVarNegTransMomMC(0),
150 fTreeCascVarNegNSigmaPion(0),
151 fTreeCascVarNegNSigmaProton(0),
152 fTreeCascVarPosNSigmaPion(0),
153 fTreeCascVarPosNSigmaProton(0),
154 fTreeCascVarBachNSigmaPion(0),
155 fTreeCascVarBachNSigmaKaon(0),
157 fTreeCascVarkITSRefitBachelor(0),
158 fTreeCascVarkITSRefitNegative(0),
159 fTreeCascVarkITSRefitPositive(0),
162 //Part A: EbyE info, Run number
163 fTreeCascVarRunNumber(0),
164 fTreeCascVarEventNumber(0),
166 //Part B: Shared Clusters
167 fTreeCascVarNegSharedClusters(0),
168 fTreeCascVarPosSharedClusters(0),
169 fTreeCascVarBachSharedClusters(0),
171 //Part C: All momenta
172 fTreeCascVarNegPx(0),
173 fTreeCascVarNegPy(0),
174 fTreeCascVarNegPz(0),
175 fTreeCascVarPosPx(0),
176 fTreeCascVarPosPy(0),
177 fTreeCascVarPosPz(0),
178 fTreeCascVarBachPx(0),
179 fTreeCascVarBachPy(0),
180 fTreeCascVarBachPz(0),
182 //------------------------------------------------
184 // --- Filled on an Event-by-event basis
185 //------------------------------------------------
186 fHistV0MultiplicityBeforeTrigSel(0),
187 fHistV0MultiplicityForTrigEvt(0),
188 fHistV0MultiplicityForSelEvt(0),
189 fHistV0MultiplicityForSelEvtNoTPCOnly(0),
190 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup(0),
191 fHistMultiplicityBeforeTrigSel(0),
192 fHistMultiplicityForTrigEvt(0),
193 fHistMultiplicity(0),
194 fHistMultiplicityNoTPCOnly(0),
195 fHistMultiplicityNoTPCOnlyNoPileup(0),
198 fHistMultiplicityV0ABeforeTrigSel(0),
199 fHistMultiplicityV0AForTrigEvt(0),
200 fHistMultiplicityV0A(0),
201 fHistMultiplicityV0ANoTPCOnly(0),
202 fHistMultiplicityV0ANoTPCOnlyNoPileup(0),
205 fHistMultiplicityZNABeforeTrigSel(0),
206 fHistMultiplicityZNAForTrigEvt(0),
207 fHistMultiplicityZNA(0),
208 fHistMultiplicityZNANoTPCOnly(0),
209 fHistMultiplicityZNANoTPCOnlyNoPileup(0),
212 fHistMultiplicityTRKBeforeTrigSel(0),
213 fHistMultiplicityTRKForTrigEvt(0),
214 fHistMultiplicityTRK(0),
215 fHistMultiplicityTRKNoTPCOnly(0),
216 fHistMultiplicityTRKNoTPCOnlyNoPileup(0),
219 fHistMultiplicitySPDBeforeTrigSel(0),
220 fHistMultiplicitySPDForTrigEvt(0),
221 fHistMultiplicitySPD(0),
222 fHistMultiplicitySPDNoTPCOnly(0),
223 fHistMultiplicitySPDNoTPCOnlyNoPileup(0),
235 AliAnalysisTaskExtractCascade::AliAnalysisTaskExtractCascade(const char *name)
236 : AliAnalysisTaskSE(name), fListHist(0), fTreeCascade(0), fPIDResponse(0), fESDtrackCuts(0), fUtils(0),
237 fkIsNuclear ( kFALSE ),
238 fkSwitchINT7 ( kFALSE ),
239 fCentralityEstimator("V0M"),
240 fkpAVertexSelection( kFALSE ),
242 fkRunVertexers ( kFALSE ),
243 //------------------------------------------------
245 //------------------------------------------------
247 fTreeCascVarCharge(0),
248 fTreeCascVarMassAsXi(0),
249 fTreeCascVarMassAsOmega(0),
252 fTreeCascVarRapMC(0),
253 fTreeCascVarRapXi(0),
254 fTreeCascVarRapOmega(0),
255 fTreeCascVarNegEta(0),
256 fTreeCascVarPosEta(0),
257 fTreeCascVarBachEta(0),
258 fTreeCascVarDCACascDaughters(0),
259 fTreeCascVarDCABachToPrimVtx(0),
260 fTreeCascVarDCAV0Daughters(0),
261 fTreeCascVarDCAV0ToPrimVtx(0),
262 fTreeCascVarDCAPosToPrimVtx(0),
263 fTreeCascVarDCANegToPrimVtx(0),
264 fTreeCascVarCascCosPointingAngle(0),
265 fTreeCascVarCascRadius(0),
266 fTreeCascVarV0Mass(0),
267 fTreeCascVarV0CosPointingAngle(0),
268 fTreeCascVarV0CosPointingAngleSpecial(0),
269 fTreeCascVarV0Radius(0),
270 fTreeCascVarLeastNbrClusters(0),
271 fTreeCascVarMultiplicity(0),
272 fTreeCascVarMultiplicityV0A(0),
273 fTreeCascVarMultiplicityZNA(0),
274 fTreeCascVarMultiplicityTRK(0),
275 fTreeCascVarMultiplicitySPD(0),
276 fTreeCascVarDistOverTotMom(0),
278 fTreeCascVarPIDBachelor(0),
279 fTreeCascVarPIDNegative(0),
280 fTreeCascVarPIDPositive(0),
281 fTreeCascVarBachTransMom(0),
282 fTreeCascVarPosTransMom(0),
283 fTreeCascVarNegTransMom(0),
284 fTreeCascVarPosTransMomMC(0),
285 fTreeCascVarNegTransMomMC(0),
286 fTreeCascVarNegNSigmaPion(0),
287 fTreeCascVarNegNSigmaProton(0),
288 fTreeCascVarPosNSigmaPion(0),
289 fTreeCascVarPosNSigmaProton(0),
290 fTreeCascVarBachNSigmaPion(0),
291 fTreeCascVarBachNSigmaKaon(0),
293 fTreeCascVarkITSRefitBachelor(0),
294 fTreeCascVarkITSRefitNegative(0),
295 fTreeCascVarkITSRefitPositive(0),
298 //Part A: EbyE info, Run number
299 fTreeCascVarRunNumber(0),
300 fTreeCascVarEventNumber(0),
302 //Part B: Shared Clusters
303 fTreeCascVarNegSharedClusters(0),
304 fTreeCascVarPosSharedClusters(0),
305 fTreeCascVarBachSharedClusters(0),
307 //Part C: All momenta
308 fTreeCascVarNegPx(0),
309 fTreeCascVarNegPy(0),
310 fTreeCascVarNegPz(0),
311 fTreeCascVarPosPx(0),
312 fTreeCascVarPosPy(0),
313 fTreeCascVarPosPz(0),
314 fTreeCascVarBachPx(0),
315 fTreeCascVarBachPy(0),
316 fTreeCascVarBachPz(0),
318 //------------------------------------------------
320 // --- Filled on an Event-by-event basis
321 //------------------------------------------------
322 fHistV0MultiplicityBeforeTrigSel(0),
323 fHistV0MultiplicityForTrigEvt(0),
324 fHistV0MultiplicityForSelEvt(0),
325 fHistV0MultiplicityForSelEvtNoTPCOnly(0),
326 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup(0),
327 fHistMultiplicityBeforeTrigSel(0),
328 fHistMultiplicityForTrigEvt(0),
329 fHistMultiplicity(0),
330 fHistMultiplicityNoTPCOnly(0),
331 fHistMultiplicityNoTPCOnlyNoPileup(0),
334 fHistMultiplicityV0ABeforeTrigSel(0),
335 fHistMultiplicityV0AForTrigEvt(0),
336 fHistMultiplicityV0A(0),
337 fHistMultiplicityV0ANoTPCOnly(0),
338 fHistMultiplicityV0ANoTPCOnlyNoPileup(0),
341 fHistMultiplicityZNABeforeTrigSel(0),
342 fHistMultiplicityZNAForTrigEvt(0),
343 fHistMultiplicityZNA(0),
344 fHistMultiplicityZNANoTPCOnly(0),
345 fHistMultiplicityZNANoTPCOnlyNoPileup(0),
348 fHistMultiplicityTRKBeforeTrigSel(0),
349 fHistMultiplicityTRKForTrigEvt(0),
350 fHistMultiplicityTRK(0),
351 fHistMultiplicityTRKNoTPCOnly(0),
352 fHistMultiplicityTRKNoTPCOnlyNoPileup(0),
355 fHistMultiplicitySPDBeforeTrigSel(0),
356 fHistMultiplicitySPDForTrigEvt(0),
357 fHistMultiplicitySPD(0),
358 fHistMultiplicitySPDNoTPCOnly(0),
359 fHistMultiplicitySPDNoTPCOnlyNoPileup(0),
370 //Set Variables for re-running the cascade vertexers (as done for MS paper)
372 // New Loose : 1st step for the 7 TeV pp analysis
374 fV0VertexerSels[0] = 33. ; // max allowed chi2
375 fV0VertexerSels[1] = 0.02; // min allowed impact parameter for the 1st daughter (LHC09a4 : 0.05)
376 fV0VertexerSels[2] = 0.02; // min allowed impact parameter for the 2nd daughter (LHC09a4 : 0.05)
377 fV0VertexerSels[3] = 2.0 ; // max allowed DCA between the daughter tracks (LHC09a4 : 0.5)
378 fV0VertexerSels[4] = 0.95; // min allowed cosine of V0's pointing angle (LHC09a4 : 0.99)
379 fV0VertexerSels[5] = 1.0 ; // min radius of the fiducial volume (LHC09a4 : 0.2)
380 fV0VertexerSels[6] = 200. ; // max radius of the fiducial volume (LHC09a4 : 100.0)
382 fCascadeVertexerSels[0] = 33. ; // max allowed chi2 (same as PDC07)
383 fCascadeVertexerSels[1] = 0.05 ; // min allowed V0 impact parameter (PDC07 : 0.05 / LHC09a4 : 0.025 )
384 fCascadeVertexerSels[2] = 0.010; // "window" around the Lambda mass (PDC07 : 0.008 / LHC09a4 : 0.010 )
385 fCascadeVertexerSels[3] = 0.03 ; // min allowed bachelor's impact parameter (PDC07 : 0.035 / LHC09a4 : 0.025 )
386 fCascadeVertexerSels[4] = 2.0 ; // max allowed DCA between the V0 and the bachelor (PDC07 : 0.1 / LHC09a4 : 0.2 )
387 fCascadeVertexerSels[5] = 0.95 ; // min allowed cosine of the cascade pointing angle (PDC07 : 0.9985 / LHC09a4 : 0.998 )
388 fCascadeVertexerSels[6] = 0.4 ; // min radius of the fiducial volume (PDC07 : 0.9 / LHC09a4 : 0.2 )
389 fCascadeVertexerSels[7] = 100. ; // max radius of the fiducial volume (PDC07 : 100 / LHC09a4 : 100 )
391 // Output slot #0 writes into a TList container (Cascade)
392 DefineOutput(1, TList::Class());
393 DefineOutput(2, TTree::Class());
397 AliAnalysisTaskExtractCascade::~AliAnalysisTaskExtractCascade()
399 //------------------------------------------------
401 //------------------------------------------------
411 //cleanup esd track cuts object too...
413 delete fESDtrackCuts;
423 //________________________________________________________________________
424 void AliAnalysisTaskExtractCascade::UserCreateOutputObjects()
429 //------------------------------------------------
431 fTreeCascade = new TTree("fTreeCascade","CascadeCandidates");
433 //------------------------------------------------
434 // fTreeCascade Branch definitions - Cascade Tree
435 //------------------------------------------------
437 //------------------------------------------------
438 // fTreeCascade Branch definitions
439 //------------------------------------------------
441 //-----------BASIC-INFO---------------------------
442 /* 1*/ fTreeCascade->Branch("fTreeCascVarCharge",&fTreeCascVarCharge,"fTreeCascVarCharge/I");
443 /* 2*/ fTreeCascade->Branch("fTreeCascVarMassAsXi",&fTreeCascVarMassAsXi,"fTreeCascVarMassAsXi/F");
444 /* 3*/ fTreeCascade->Branch("fTreeCascVarMassAsOmega",&fTreeCascVarMassAsOmega,"fTreeCascVarMassAsOmega/F");
445 /* 4*/ fTreeCascade->Branch("fTreeCascVarPt",&fTreeCascVarPt,"fTreeCascVarPt/F");
446 /* 5*/ fTreeCascade->Branch("fTreeCascVarRapXi",&fTreeCascVarRapXi,"fTreeCascVarRapXi/F");
447 /* 6*/ fTreeCascade->Branch("fTreeCascVarRapOmega",&fTreeCascVarRapOmega,"fTreeCascVarRapOmega/F");
448 /* 7*/ fTreeCascade->Branch("fTreeCascVarNegEta",&fTreeCascVarNegEta,"fTreeCascVarNegEta/F");
449 /* 8*/ fTreeCascade->Branch("fTreeCascVarPosEta",&fTreeCascVarPosEta,"fTreeCascVarPosEta/F");
450 /* 9*/ fTreeCascade->Branch("fTreeCascVarBachEta",&fTreeCascVarBachEta,"fTreeCascVarBachEta/F");
451 //-----------INFO-FOR-CUTS------------------------
452 /*10*/ fTreeCascade->Branch("fTreeCascVarDCACascDaughters",&fTreeCascVarDCACascDaughters,"fTreeCascVarDCACascDaughters/F");
453 /*11*/ fTreeCascade->Branch("fTreeCascVarDCABachToPrimVtx",&fTreeCascVarDCABachToPrimVtx,"fTreeCascVarDCABachToPrimVtx/F");
454 /*12*/ fTreeCascade->Branch("fTreeCascVarDCAV0Daughters",&fTreeCascVarDCAV0Daughters,"fTreeCascVarDCAV0Daughters/F");
455 /*13*/ fTreeCascade->Branch("fTreeCascVarDCAV0ToPrimVtx",&fTreeCascVarDCAV0ToPrimVtx,"fTreeCascVarDCAV0ToPrimVtx/F");
456 /*14*/ fTreeCascade->Branch("fTreeCascVarDCAPosToPrimVtx",&fTreeCascVarDCAPosToPrimVtx,"fTreeCascVarDCAPosToPrimVtx/F");
457 /*15*/ fTreeCascade->Branch("fTreeCascVarDCANegToPrimVtx",&fTreeCascVarDCANegToPrimVtx,"fTreeCascVarDCANegToPrimVtx/F");
458 /*16*/ fTreeCascade->Branch("fTreeCascVarCascCosPointingAngle",&fTreeCascVarCascCosPointingAngle,"fTreeCascVarCascCosPointingAngle/F");
459 /*17*/ fTreeCascade->Branch("fTreeCascVarCascRadius",&fTreeCascVarCascRadius,"fTreeCascVarCascRadius/F");
460 /*18*/ fTreeCascade->Branch("fTreeCascVarV0Mass",&fTreeCascVarV0Mass,"fTreeCascVarV0Mass/F");
461 /*19*/ fTreeCascade->Branch("fTreeCascVarV0CosPointingAngle",&fTreeCascVarV0CosPointingAngle,"fTreeCascVarV0CosPointingAngle/F");
462 /*19*/ fTreeCascade->Branch("fTreeCascVarV0CosPointingAngleSpecial",&fTreeCascVarV0CosPointingAngleSpecial,"fTreeCascVarV0CosPointingAngleSpecial/F");
463 /*20*/ fTreeCascade->Branch("fTreeCascVarV0Radius",&fTreeCascVarV0Radius,"fTreeCascVarV0Radius/F");
464 /*21*/ fTreeCascade->Branch("fTreeCascVarLeastNbrClusters",&fTreeCascVarLeastNbrClusters,"fTreeCascVarLeastNbrClusters/I");
465 //-----------MULTIPLICITY-INFO--------------------
466 /*22*/ fTreeCascade->Branch("fTreeCascVarMultiplicity",&fTreeCascVarMultiplicity,"fTreeCascVarMultiplicity/I");
467 /*22*/ fTreeCascade->Branch("fTreeCascVarMultiplicityV0A",&fTreeCascVarMultiplicityV0A,"fTreeCascVarMultiplicityV0A/I");
468 /*22*/ fTreeCascade->Branch("fTreeCascVarMultiplicityZNA",&fTreeCascVarMultiplicityZNA,"fTreeCascVarMultiplicityZNA/I");
469 /*22*/ fTreeCascade->Branch("fTreeCascVarMultiplicityTRK",&fTreeCascVarMultiplicityTRK,"fTreeCascVarMultiplicityTRK/I");
470 /*22*/ fTreeCascade->Branch("fTreeCascVarMultiplicitySPD",&fTreeCascVarMultiplicitySPD,"fTreeCascVarMultiplicitySPD/I");
471 //-----------DECAY-LENGTH-INFO--------------------
472 /*23*/ fTreeCascade->Branch("fTreeCascVarDistOverTotMom",&fTreeCascVarDistOverTotMom,"fTreeCascVarDistOverTotMom/F");
473 //------------------------------------------------
474 /*24*/ fTreeCascade->Branch("fTreeCascVarNegNSigmaPion",&fTreeCascVarNegNSigmaPion,"fTreeCascVarNegNSigmaPion/F");
475 /*25*/ fTreeCascade->Branch("fTreeCascVarNegNSigmaProton",&fTreeCascVarNegNSigmaProton,"fTreeCascVarNegNSigmaProton/F");
476 /*26*/ fTreeCascade->Branch("fTreeCascVarPosNSigmaPion",&fTreeCascVarPosNSigmaPion,"fTreeCascVarPosNSigmaPion/F");
477 /*27*/ fTreeCascade->Branch("fTreeCascVarPosNSigmaProton",&fTreeCascVarPosNSigmaProton,"fTreeCascVarPosNSigmaProton/F");
478 /*28*/ fTreeCascade->Branch("fTreeCascVarBachNSigmaPion",&fTreeCascVarBachNSigmaPion,"fTreeCascVarBachNSigmaPion/F");
479 /*29*/ fTreeCascade->Branch("fTreeCascVarBachNSigmaKaon",&fTreeCascVarBachNSigmaKaon,"fTreeCascVarBachNSigmaKaon/F");
481 //Commented out: not needed since all momenta provided! (less info)
482 /*30*/ //fTreeCascade->Branch("fTreeCascVarBachTransMom",&fTreeCascVarBachTransMom,"fTreeCascVarBachTransMom/F");
483 /*30*/ //fTreeCascade->Branch("fTreeCascVarPosTransMom",&fTreeCascVarPosTransMom,"fTreeCascVarPosTransMom/F");
484 /*31*/ //fTreeCascade->Branch("fTreeCascVarNegTransMom",&fTreeCascVarNegTransMom,"fTreeCascVarNegTransMom/F");
486 /*29*/ fTreeCascade->Branch("fTreeCascVarkITSRefitBachelor",&fTreeCascVarkITSRefitBachelor,"fTreeCascVarkITSRefitBachelor/O");
487 /*29*/ fTreeCascade->Branch("fTreeCascVarkITSRefitNegative",&fTreeCascVarkITSRefitNegative,"fTreeCascVarkITSRefitNegative/O");
488 /*29*/ fTreeCascade->Branch("fTreeCascVarkITSRefitPositive",&fTreeCascVarkITSRefitPositive,"fTreeCascVarkITSRefitPositive/O");
490 //-----------Debugging information----------------
491 //Part A: Event-by-event, run-by-run debugging
492 fTreeCascade->Branch("fTreeCascVarRunNumber",&fTreeCascVarRunNumber,"fTreeCascVarRunNumber/I");
493 fTreeCascade->Branch("fTreeCascVarEventNumber",&fTreeCascVarEventNumber,"fTreeCascVarEventNumber/l");
495 //Part B: Shared Clusters for all daughter tracks
496 fTreeCascade->Branch("fTreeCascVarNegSharedClusters",&fTreeCascVarNegSharedClusters,"fTreeCascVarNegSharedClusters/I");
497 fTreeCascade->Branch("fTreeCascVarPosSharedClusters",&fTreeCascVarPosSharedClusters,"fTreeCascVarPosSharedClusters/I");
498 fTreeCascade->Branch("fTreeCascVarBachSharedClusters",&fTreeCascVarBachSharedClusters,"fTreeCascVarBachSharedClusters/I");
500 //Part C: All Momenta of all daughters
501 fTreeCascade->Branch("fTreeCascVarNegPx",&fTreeCascVarNegPx,"fTreeCascVarNegPx/F");
502 fTreeCascade->Branch("fTreeCascVarNegPy",&fTreeCascVarNegPy,"fTreeCascVarNegPy/F");
503 fTreeCascade->Branch("fTreeCascVarNegPz",&fTreeCascVarNegPz,"fTreeCascVarNegPz/F");
504 fTreeCascade->Branch("fTreeCascVarPosPx",&fTreeCascVarPosPx,"fTreeCascVarPosPx/F");
505 fTreeCascade->Branch("fTreeCascVarPosPy",&fTreeCascVarPosPy,"fTreeCascVarPosPy/F");
506 fTreeCascade->Branch("fTreeCascVarPosPz",&fTreeCascVarPosPz,"fTreeCascVarPosPz/F");
507 fTreeCascade->Branch("fTreeCascVarBachPx",&fTreeCascVarBachPx,"fTreeCascVarBachPx/F");
508 fTreeCascade->Branch("fTreeCascVarBachPy",&fTreeCascVarBachPy,"fTreeCascVarBachPy/F");
509 fTreeCascade->Branch("fTreeCascVarBachPz",&fTreeCascVarBachPz,"fTreeCascVarBachPz/F");
510 //------------------------------------------------
512 //------------------------------------------------
513 // Particle Identification Setup
514 //------------------------------------------------
516 AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
517 AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
518 fPIDResponse = inputHandler->GetPIDResponse();
522 if(! fESDtrackCuts ){
523 fESDtrackCuts = new AliESDtrackCuts();
526 fUtils = new AliAnalysisUtils();
529 //------------------------------------------------
530 // V0 Multiplicity Histograms
531 //------------------------------------------------
535 fListHist = new TList();
536 fListHist->SetOwner(); // See http://root.cern.ch/root/html/TCollection.html#TCollection:SetOwner
539 if(! fHistV0MultiplicityBeforeTrigSel) {
540 fHistV0MultiplicityBeforeTrigSel = new TH1F("fHistV0MultiplicityBeforeTrigSel",
541 "V0s per event (before Trig. Sel.);Nbr of V0s/Evt;Events",
543 fListHist->Add(fHistV0MultiplicityBeforeTrigSel);
546 if(! fHistV0MultiplicityForTrigEvt) {
547 fHistV0MultiplicityForTrigEvt = new TH1F("fHistV0MultiplicityForTrigEvt",
548 "V0s per event (for triggered evt);Nbr of V0s/Evt;Events",
550 fListHist->Add(fHistV0MultiplicityForTrigEvt);
553 if(! fHistV0MultiplicityForSelEvt) {
554 fHistV0MultiplicityForSelEvt = new TH1F("fHistV0MultiplicityForSelEvt",
555 "V0s per event;Nbr of V0s/Evt;Events",
557 fListHist->Add(fHistV0MultiplicityForSelEvt);
560 if(! fHistV0MultiplicityForSelEvtNoTPCOnly) {
561 fHistV0MultiplicityForSelEvtNoTPCOnly = new TH1F("fHistV0MultiplicityForSelEvtNoTPCOnly",
562 "V0s per event;Nbr of V0s/Evt;Events",
564 fListHist->Add(fHistV0MultiplicityForSelEvtNoTPCOnly);
566 if(! fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup) {
567 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup = new TH1F("fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup",
568 "V0s per event;Nbr of V0s/Evt;Events",
570 fListHist->Add(fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup);
573 //------------------------------------------------
574 // Track Multiplicity Histograms
575 //------------------------------------------------
577 if(! fHistMultiplicityBeforeTrigSel) {
578 fHistMultiplicityBeforeTrigSel = new TH1F("fHistMultiplicityBeforeTrigSel",
579 "Tracks per event;Nbr of Tracks;Events",
581 fListHist->Add(fHistMultiplicityBeforeTrigSel);
583 if(! fHistMultiplicityForTrigEvt) {
584 fHistMultiplicityForTrigEvt = new TH1F("fHistMultiplicityForTrigEvt",
585 "Tracks per event;Nbr of Tracks;Events",
587 fListHist->Add(fHistMultiplicityForTrigEvt);
589 if(! fHistMultiplicity) {
590 fHistMultiplicity = new TH1F("fHistMultiplicity",
591 "Tracks per event;Nbr of Tracks;Events",
593 fListHist->Add(fHistMultiplicity);
595 if(! fHistMultiplicityNoTPCOnly) {
596 fHistMultiplicityNoTPCOnly = new TH1F("fHistMultiplicityNoTPCOnly",
597 "Tracks per event;Nbr of Tracks;Events",
599 fListHist->Add(fHistMultiplicityNoTPCOnly);
601 if(! fHistMultiplicityNoTPCOnlyNoPileup) {
602 fHistMultiplicityNoTPCOnlyNoPileup = new TH1F("fHistMultiplicityNoTPCOnlyNoPileup",
603 "Tracks per event;Nbr of Tracks;Events",
605 fListHist->Add(fHistMultiplicityNoTPCOnlyNoPileup);
609 //V0A Centrality (if PbPb / pPb)
610 if(! fHistMultiplicityV0ABeforeTrigSel) {
611 fHistMultiplicityV0ABeforeTrigSel = new TH1F("fHistMultiplicityV0ABeforeTrigSel",
612 "Centrality Distribution: V0A;V0A Centrality;Events",
614 fListHist->Add(fHistMultiplicityV0ABeforeTrigSel);
616 if(! fHistMultiplicityV0AForTrigEvt) {
617 fHistMultiplicityV0AForTrigEvt = new TH1F("fHistMultiplicityV0AForTrigEvt",
618 "Centrality Distribution: V0A;V0A Centrality;Events",
620 fListHist->Add(fHistMultiplicityV0AForTrigEvt);
622 if(! fHistMultiplicityV0A) {
623 fHistMultiplicityV0A = new TH1F("fHistMultiplicityV0A",
624 "Centrality Distribution: V0A;V0A Centrality;Events",
626 fListHist->Add(fHistMultiplicityV0A);
628 if(! fHistMultiplicityV0ANoTPCOnly) {
629 fHistMultiplicityV0ANoTPCOnly = new TH1F("fHistMultiplicityV0ANoTPCOnly",
630 "Centrality Distribution: V0A;V0A Centrality;Events",
632 fListHist->Add(fHistMultiplicityV0ANoTPCOnly);
634 if(! fHistMultiplicityV0ANoTPCOnlyNoPileup) {
635 fHistMultiplicityV0ANoTPCOnlyNoPileup = new TH1F("fHistMultiplicityV0ANoTPCOnlyNoPileup",
636 "Centrality Distribution: V0A;V0A Centrality;Events",
638 fListHist->Add(fHistMultiplicityV0ANoTPCOnlyNoPileup);
641 //ZNA Centrality (if PbPb / pPb)
642 if(! fHistMultiplicityZNABeforeTrigSel) {
643 fHistMultiplicityZNABeforeTrigSel = new TH1F("fHistMultiplicityZNABeforeTrigSel",
644 "Centrality Distribution: ZNA;ZNA Centrality;Events",
646 fListHist->Add(fHistMultiplicityZNABeforeTrigSel);
648 if(! fHistMultiplicityZNAForTrigEvt) {
649 fHistMultiplicityZNAForTrigEvt = new TH1F("fHistMultiplicityZNAForTrigEvt",
650 "Centrality Distribution: ZNA;ZNA Centrality;Events",
652 fListHist->Add(fHistMultiplicityZNAForTrigEvt);
654 if(! fHistMultiplicityZNA) {
655 fHistMultiplicityZNA = new TH1F("fHistMultiplicityZNA",
656 "Centrality Distribution: ZNA;ZNA Centrality;Events",
658 fListHist->Add(fHistMultiplicityZNA);
660 if(! fHistMultiplicityZNANoTPCOnly) {
661 fHistMultiplicityZNANoTPCOnly = new TH1F("fHistMultiplicityZNANoTPCOnly",
662 "Centrality Distribution: ZNA;ZNA Centrality;Events",
664 fListHist->Add(fHistMultiplicityZNANoTPCOnly);
666 if(! fHistMultiplicityZNANoTPCOnlyNoPileup) {
667 fHistMultiplicityZNANoTPCOnlyNoPileup = new TH1F("fHistMultiplicityZNANoTPCOnlyNoPileup",
668 "Centrality Distribution: ZNA;ZNA Centrality;Events",
670 fListHist->Add(fHistMultiplicityZNANoTPCOnlyNoPileup);
673 //TRK Centrality (if PbPb / pPb)
674 if(! fHistMultiplicityTRKBeforeTrigSel) {
675 fHistMultiplicityTRKBeforeTrigSel = new TH1F("fHistMultiplicityTRKBeforeTrigSel",
676 "Centrality Distribution: TRK;TRK Centrality;Events",
678 fListHist->Add(fHistMultiplicityTRKBeforeTrigSel);
680 if(! fHistMultiplicityTRKForTrigEvt) {
681 fHistMultiplicityTRKForTrigEvt = new TH1F("fHistMultiplicityTRKForTrigEvt",
682 "Centrality Distribution: TRK;TRK Centrality;Events",
684 fListHist->Add(fHistMultiplicityTRKForTrigEvt);
686 if(! fHistMultiplicityTRK) {
687 fHistMultiplicityTRK = new TH1F("fHistMultiplicityTRK",
688 "Centrality Distribution: TRK;TRK Centrality;Events",
690 fListHist->Add(fHistMultiplicityTRK);
692 if(! fHistMultiplicityTRKNoTPCOnly) {
693 fHistMultiplicityTRKNoTPCOnly = new TH1F("fHistMultiplicityTRKNoTPCOnly",
694 "Centrality Distribution: TRK;TRK Centrality;Events",
696 fListHist->Add(fHistMultiplicityTRKNoTPCOnly);
698 if(! fHistMultiplicityTRKNoTPCOnlyNoPileup) {
699 fHistMultiplicityTRKNoTPCOnlyNoPileup = new TH1F("fHistMultiplicityTRKNoTPCOnlyNoPileup",
700 "Centrality Distribution: TRK;TRK Centrality;Events",
702 fListHist->Add(fHistMultiplicityTRKNoTPCOnlyNoPileup);
705 //SPD Centrality (if PbPb / pPb)
706 if(! fHistMultiplicitySPDBeforeTrigSel) {
707 fHistMultiplicitySPDBeforeTrigSel = new TH1F("fHistMultiplicitySPDBeforeTrigSel",
708 "Centrality Distribution: SPD;SPD Centrality;Events",
710 fListHist->Add(fHistMultiplicitySPDBeforeTrigSel);
712 if(! fHistMultiplicitySPDForTrigEvt) {
713 fHistMultiplicitySPDForTrigEvt = new TH1F("fHistMultiplicitySPDForTrigEvt",
714 "Centrality Distribution: SPD;SPD Centrality;Events",
716 fListHist->Add(fHistMultiplicitySPDForTrigEvt);
718 if(! fHistMultiplicitySPD) {
719 fHistMultiplicitySPD = new TH1F("fHistMultiplicitySPD",
720 "Centrality Distribution: SPD;SPD Centrality;Events",
722 fListHist->Add(fHistMultiplicitySPD);
724 if(! fHistMultiplicitySPDNoTPCOnly) {
725 fHistMultiplicitySPDNoTPCOnly = new TH1F("fHistMultiplicitySPDNoTPCOnly",
726 "Centrality Distribution: SPD;SPD Centrality;Events",
728 fListHist->Add(fHistMultiplicitySPDNoTPCOnly);
730 if(! fHistMultiplicitySPDNoTPCOnlyNoPileup) {
731 fHistMultiplicitySPDNoTPCOnlyNoPileup = new TH1F("fHistMultiplicitySPDNoTPCOnlyNoPileup",
732 "Centrality Distribution: SPD;SPD Centrality;Events",
734 fListHist->Add(fHistMultiplicitySPDNoTPCOnlyNoPileup);
738 //----------------------------------
739 // Primary Vertex Position Histos
740 //----------------------------------
743 fHistPVx = new TH1F("fHistPVx",
744 "PV x position;Nbr of Evts;x",
746 fListHist->Add(fHistPVx);
749 fHistPVy = new TH1F("fHistPVy",
750 "PV y position;Nbr of Evts;y",
752 fListHist->Add(fHistPVy);
755 fHistPVz = new TH1F("fHistPVz",
756 "PV z position;Nbr of Evts;z",
758 fListHist->Add(fHistPVz);
761 if(! fHistPVxAnalysis) {
762 fHistPVxAnalysis = new TH1F("fHistPVxAnalysis",
763 "PV x position;Nbr of Evts;x",
765 fListHist->Add(fHistPVxAnalysis);
767 if(! fHistPVyAnalysis) {
768 fHistPVyAnalysis = new TH1F("fHistPVyAnalysis",
769 "PV y position;Nbr of Evts;y",
771 fListHist->Add(fHistPVyAnalysis);
773 if(! fHistPVzAnalysis) {
774 fHistPVzAnalysis = new TH1F("fHistPVzAnalysis",
775 "PV z position;Nbr of Evts;z",
777 fListHist->Add(fHistPVzAnalysis);
780 //List of Histograms: Normal
781 PostData(1, fListHist);
783 //TTree Object: Saved to base directory. Should cache to disk while saving.
784 //(Important to avoid excessive memory usage, particularly when merging)
785 PostData(2, fTreeCascade);
787 }// end UserCreateOutputObjects
790 //________________________________________________________________________
791 void AliAnalysisTaskExtractCascade::UserExec(Option_t *)
794 // Called for each event
796 AliESDEvent *lESDevent = 0x0;
798 Int_t lNumberOfV0s = -1;
799 Double_t lTrkgPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
800 Double_t lBestPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
801 Double_t lMagneticField = -10.;
803 // Connect to the InputEvent
804 // After these lines, we should have an ESD/AOD event + the number of V0s in it.
806 // Appropriate for ESD analysis!
808 lESDevent = dynamic_cast<AliESDEvent*>( InputEvent() );
810 AliWarning("ERROR: lESDevent not available \n");
814 //--- Acquisition of exact event ID
815 fTreeCascVarRunNumber = lESDevent->GetRunNumber();
816 fTreeCascVarEventNumber =
817 ( ( ((ULong64_t)lESDevent->GetPeriodNumber() ) << 36 ) |
818 ( ((ULong64_t)lESDevent->GetOrbitNumber () ) << 12 ) |
819 ((ULong64_t)lESDevent->GetBunchCrossNumber() ) );
822 //------------------------------------------------
823 // Multiplicity Information Acquistion
824 //------------------------------------------------
826 //REVISED multiplicity estimator after 'multiplicity day' (2011)
827 Int_t lMultiplicity = -100;
828 Int_t lMultiplicityV0A = -100;
829 Int_t lMultiplicityZNA = -100;
830 Int_t lMultiplicityTRK = -100;
831 Int_t lMultiplicitySPD = -100;
834 if(fkIsNuclear == kFALSE) lMultiplicity = fESDtrackCuts->GetReferenceMultiplicity(lESDevent, AliESDtrackCuts::kTrackletsITSTPC, fEtaRefMult );
836 //---> If this is a nuclear collision, then go nuclear on "multiplicity" variable...
837 //---> Warning: Experimental
838 if(fkIsNuclear == kTRUE){
839 AliCentrality* centrality;
840 centrality = lESDevent->GetCentrality();
841 lMultiplicity = ( ( Int_t ) ( centrality->GetCentralityPercentile( fCentralityEstimator.Data() ) ) );
842 lMultiplicityV0A = ( ( Int_t ) ( centrality->GetCentralityPercentile( "V0A" ) ) );
843 lMultiplicityZNA = ( ( Int_t ) ( centrality->GetCentralityPercentile( "ZNA" ) ) );
844 lMultiplicityTRK = ( ( Int_t ) ( centrality->GetCentralityPercentile( "TRK" ) ) );
845 lMultiplicitySPD = ( ( Int_t ) ( centrality->GetCentralityPercentile( "CL1" ) ) );
846 if (centrality->GetQuality()>1) {
847 PostData(1, fListHist);
848 PostData(2, fTreeCascade);
853 //Set variable for filling tree afterwards!
854 //---> pp case......: GetReferenceMultiplicity
855 //---> Pb-Pb case...: Centrality by V0M
857 fTreeCascVarMultiplicity = lMultiplicity;
858 fTreeCascVarMultiplicityV0A = lMultiplicityV0A;
859 fTreeCascVarMultiplicityZNA = lMultiplicityZNA;
860 fTreeCascVarMultiplicityTRK = lMultiplicityTRK;
861 fTreeCascVarMultiplicitySPD = lMultiplicitySPD;
863 fHistV0MultiplicityBeforeTrigSel->Fill ( lESDevent->GetNumberOfV0s() );
864 fHistMultiplicityBeforeTrigSel->Fill ( lMultiplicity );
865 fHistMultiplicityV0ABeforeTrigSel->Fill ( lMultiplicityV0A );
866 fHistMultiplicityZNABeforeTrigSel->Fill ( lMultiplicityZNA );
867 fHistMultiplicityTRKBeforeTrigSel->Fill ( lMultiplicityTRK );
868 fHistMultiplicitySPDBeforeTrigSel->Fill ( lMultiplicitySPD );
870 //------------------------------------------------
872 //------------------------------------------------
874 UInt_t maskIsSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
875 Bool_t isSelected = 0;
876 isSelected = (maskIsSelected & AliVEvent::kMB) == AliVEvent::kMB;
878 //pA triggering: CINT7
879 if( fkSwitchINT7 ) isSelected = (maskIsSelected & AliVEvent::kINT7) == AliVEvent::kINT7;
881 //Standard Min-Bias Selection
882 if ( ! isSelected ) {
883 PostData(1, fListHist);
884 PostData(2, fTreeCascade);
888 //------------------------------------------------
889 // Rerun cascade vertexer!
890 //------------------------------------------------
892 if( fkRunVertexers ){
893 lESDevent->ResetCascades();
894 lESDevent->ResetV0s();
896 AliV0vertexer lV0vtxer;
897 AliCascadeVertexer lCascVtxer;
899 lV0vtxer.SetDefaultCuts(fV0VertexerSels);
900 lCascVtxer.SetDefaultCuts(fCascadeVertexerSels);
902 lV0vtxer.Tracks2V0vertices(lESDevent);
903 lCascVtxer.V0sTracks2CascadeVertices(lESDevent);
905 //------------------------------------------------
906 // After Trigger Selection
907 //------------------------------------------------
909 lNumberOfV0s = lESDevent->GetNumberOfV0s();
911 //Set variable for filling tree afterwards!
912 fHistV0MultiplicityForTrigEvt->Fill(lNumberOfV0s);
913 fHistMultiplicityForTrigEvt->Fill ( lMultiplicity );
914 fHistMultiplicityV0AForTrigEvt ->Fill( lMultiplicityV0A );
915 fHistMultiplicityZNAForTrigEvt ->Fill( lMultiplicityZNA );
916 fHistMultiplicityTRKForTrigEvt ->Fill( lMultiplicityTRK );
917 fHistMultiplicitySPDForTrigEvt ->Fill( lMultiplicitySPD );
919 //------------------------------------------------
920 // Getting: Primary Vertex + MagField Info
921 //------------------------------------------------
923 const AliESDVertex *lPrimaryTrackingESDVtx = lESDevent->GetPrimaryVertexTracks();
924 // get the vtx stored in ESD found with tracks
925 lPrimaryTrackingESDVtx->GetXYZ( lTrkgPrimaryVtxPos );
927 const AliESDVertex *lPrimaryBestESDVtx = lESDevent->GetPrimaryVertex();
928 // get the best primary vertex available for the event
929 // As done in AliCascadeVertexer, we keep the one which is the best one available.
930 // between : Tracking vertex > SPD vertex > TPC vertex > default SPD vertex
931 // This one will be used for next calculations (DCA essentially)
932 lPrimaryBestESDVtx->GetXYZ( lBestPrimaryVtxPos );
934 Double_t lPrimaryVtxPosition[3];
935 const AliVVertex *primaryVtx = lESDevent->GetPrimaryVertex();
936 lPrimaryVtxPosition[0] = primaryVtx->GetX();
937 lPrimaryVtxPosition[1] = primaryVtx->GetY();
938 lPrimaryVtxPosition[2] = primaryVtx->GetZ();
939 fHistPVx->Fill( lPrimaryVtxPosition[0] );
940 fHistPVy->Fill( lPrimaryVtxPosition[1] );
941 fHistPVz->Fill( lPrimaryVtxPosition[2] );
943 //------------------------------------------------
944 // Primary Vertex Requirements Section:
945 // ---> pp and PbPb: Only requires |z|<10cm
946 // ---> pPb: all requirements checked at this stage
947 //------------------------------------------------
949 //Roberto's PV selection criteria, implemented 17th April 2013
951 /* vertex selection */
952 Bool_t fHasVertex = kFALSE;
954 const AliESDVertex *vertex = lESDevent->GetPrimaryVertexTracks();
955 if (vertex->GetNContributors() < 1) {
956 vertex = lESDevent->GetPrimaryVertexSPD();
957 if (vertex->GetNContributors() < 1) fHasVertex = kFALSE;
958 else fHasVertex = kTRUE;
959 TString vtxTyp = vertex->GetTitle();
961 vertex->GetCovarianceMatrix(cov);
962 Double_t zRes = TMath::Sqrt(cov[5]);
963 if (vtxTyp.Contains("vertexer:Z") && (zRes>0.25)) fHasVertex = kFALSE;
965 else fHasVertex = kTRUE;
967 //Is First event in chunk rejection: Still present!
968 if(fkpAVertexSelection==kTRUE && fHasVertex == kFALSE) {
969 AliWarning("Pb / | PV does not satisfy selection criteria!");
970 PostData(1, fListHist);
971 PostData(2, fTreeCascade);
975 //Is First event in chunk rejection: Still present!
976 if(fkpAVertexSelection==kTRUE && fUtils->IsFirstEventInChunk(lESDevent)) {
977 AliWarning("Pb / | This is the first event in the chunk!");
978 PostData(1, fListHist);
979 PostData(2, fTreeCascade);
983 //17 April Fix: Always do primary vertex Z selection, after pA vertex selection from Roberto
984 if(TMath::Abs(lBestPrimaryVtxPos[2]) > 10.0) {
985 AliWarning("Pb / | Z position of Best Prim Vtx | > 10.0 cm ... return !");
986 PostData(1, fListHist);
987 PostData(2, fTreeCascade);
992 lMagneticField = lESDevent->GetMagneticField( );
993 fHistV0MultiplicityForSelEvt ->Fill( lNumberOfV0s );
994 fHistMultiplicity->Fill(lMultiplicity);
995 fHistMultiplicityV0A->Fill(lMultiplicityV0A);
996 fHistMultiplicityZNA->Fill(lMultiplicityZNA);
997 fHistMultiplicityTRK->Fill(lMultiplicityTRK);
998 fHistMultiplicitySPD->Fill(lMultiplicitySPD);
1000 //------------------------------------------------
1001 // SKIP: Events with well-established PVtx
1002 //------------------------------------------------
1004 const AliESDVertex *lPrimaryTrackingESDVtxCheck = lESDevent->GetPrimaryVertexTracks();
1005 const AliESDVertex *lPrimarySPDVtx = lESDevent->GetPrimaryVertexSPD();
1006 if (!lPrimarySPDVtx->GetStatus() && !lPrimaryTrackingESDVtxCheck->GetStatus() && fkpAVertexSelection==kFALSE ){
1007 AliWarning("Pb / No SPD prim. vertex nor prim. Tracking vertex ... return !");
1008 PostData(1, fListHist);
1009 PostData(2, fTreeCascade);
1012 fHistV0MultiplicityForSelEvtNoTPCOnly ->Fill( lNumberOfV0s );
1013 fHistMultiplicityNoTPCOnly->Fill(lMultiplicity);
1014 fHistMultiplicityV0ANoTPCOnly->Fill(lMultiplicityV0A);
1015 fHistMultiplicityZNANoTPCOnly->Fill(lMultiplicityZNA);
1016 fHistMultiplicityTRKNoTPCOnly->Fill(lMultiplicityTRK);
1017 fHistMultiplicitySPDNoTPCOnly->Fill(lMultiplicitySPD);
1019 //------------------------------------------------
1020 // Pileup Rejection Studies
1021 //------------------------------------------------
1023 // FIXME : quality selection regarding pile-up rejection
1024 if(lESDevent->IsPileupFromSPD() && !fkIsNuclear ){// minContributors=3, minZdist=0.8, nSigmaZdist=3., nSigmaDiamXY=2., nSigmaDiamZ=5. -> see http://alisoft.cern.ch/viewvc/trunk/STEER/AliESDEvent.h?root=AliRoot&r1=41914&r2=42199&pathrev=42199
1025 AliWarning("Pb / This is tagged as Pileup from SPD... return !");
1026 PostData(1, fListHist);
1027 PostData(2, fTreeCascade);
1030 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup ->Fill( lNumberOfV0s );
1031 fHistMultiplicityNoTPCOnlyNoPileup->Fill(lMultiplicity);
1032 fHistMultiplicityV0ANoTPCOnlyNoPileup->Fill(lMultiplicityV0A);
1033 fHistMultiplicityZNANoTPCOnlyNoPileup->Fill(lMultiplicityZNA);
1034 fHistMultiplicityTRKNoTPCOnlyNoPileup->Fill(lMultiplicityTRK);
1035 fHistMultiplicitySPDNoTPCOnlyNoPileup->Fill(lMultiplicitySPD);
1037 //Do control histograms without the IsFromVertexerZ events, but consider them in analysis...
1038 if( ! (lESDevent->GetPrimaryVertex()->IsFromVertexerZ() ) ){
1039 fHistPVxAnalysis->Fill( lPrimaryVtxPosition[0] );
1040 fHistPVyAnalysis->Fill( lPrimaryVtxPosition[1] );
1041 fHistPVzAnalysis->Fill( lPrimaryVtxPosition[2] );
1044 //------------------------------------------------
1045 // MAIN CASCADE LOOP STARTS HERE
1046 //------------------------------------------------
1047 // Code Credit: Antonin Maire (thanks^100)
1048 // ---> This is an adaptation
1050 Long_t ncascades = 0;
1051 ncascades = lESDevent->GetNumberOfCascades();
1053 for (Int_t iXi = 0; iXi < ncascades; iXi++){
1054 //------------------------------------------------
1056 //------------------------------------------------
1057 //Double_t lTrkgPrimaryVtxRadius3D = -500.0;
1058 //Double_t lBestPrimaryVtxRadius3D = -500.0;
1060 // - 1st part of initialisation : variables needed to store AliESDCascade data members
1061 Double_t lEffMassXi = 0. ;
1062 //Double_t lChi2Xi = -1. ;
1063 Double_t lDcaXiDaughters = -1. ;
1064 Double_t lXiCosineOfPointingAngle = -1. ;
1065 Double_t lPosXi[3] = { -1000.0, -1000.0, -1000.0 };
1066 Double_t lXiRadius = -1000. ;
1068 // - 2nd part of initialisation : Nbr of clusters within TPC for the 3 daughter cascade tracks
1069 Int_t lPosTPCClusters = -1; // For ESD only ...//FIXME : wait for availability in AOD
1070 Int_t lNegTPCClusters = -1; // For ESD only ...
1071 Int_t lBachTPCClusters = -1; // For ESD only ...
1073 // - 3rd part of initialisation : about V0 part in cascades
1074 Double_t lInvMassLambdaAsCascDghter = 0.;
1075 //Double_t lV0Chi2Xi = -1. ;
1076 Double_t lDcaV0DaughtersXi = -1.;
1078 Double_t lDcaBachToPrimVertexXi = -1., lDcaV0ToPrimVertexXi = -1.;
1079 Double_t lDcaPosToPrimVertexXi = -1.;
1080 Double_t lDcaNegToPrimVertexXi = -1.;
1081 Double_t lV0CosineOfPointingAngleXi = -1. ;
1082 Double_t lV0CosineOfPointingAngleXiSpecial = -1. ;
1083 Double_t lPosV0Xi[3] = { -1000. , -1000., -1000. }; // Position of VO coming from cascade
1084 Double_t lV0RadiusXi = -1000.0;
1085 Double_t lV0quality = 0.;
1087 // - 4th part of initialisation : Effective masses
1088 Double_t lInvMassXiMinus = 0.;
1089 Double_t lInvMassXiPlus = 0.;
1090 Double_t lInvMassOmegaMinus = 0.;
1091 Double_t lInvMassOmegaPlus = 0.;
1093 // - 6th part of initialisation : extra info for QA
1094 Double_t lXiMomX = 0. , lXiMomY = 0., lXiMomZ = 0.;
1095 Double_t lXiTransvMom = 0. ;
1096 Double_t lXiTransvMomMC= 0. ;
1097 Double_t lXiTotMom = 0. ;
1099 Double_t lBachMomX = 0., lBachMomY = 0., lBachMomZ = 0.;
1100 //Double_t lBachTransvMom = 0.;
1101 //Double_t lBachTotMom = 0.;
1103 fTreeCascVarNegNSigmaPion = -100;
1104 fTreeCascVarNegNSigmaProton = -100;
1105 fTreeCascVarPosNSigmaPion = -100;
1106 fTreeCascVarPosNSigmaProton = -100;
1107 fTreeCascVarBachNSigmaPion = -100;
1108 fTreeCascVarBachNSigmaKaon = -100;
1110 Short_t lChargeXi = -2;
1111 //Double_t lV0toXiCosineOfPointingAngle = 0. ;
1113 Double_t lRapXi = -20.0, lRapOmega = -20.0, lRapMC = -20.0; // lEta = -20.0, lTheta = 360., lPhi = 720. ;
1114 //Double_t lAlphaXi = -200., lPtArmXi = -200.0;
1116 // -------------------------------------
1117 // II.ESD - Calculation Part dedicated to Xi vertices (ESD)
1119 AliESDcascade *xi = lESDevent->GetCascade(iXi);
1123 // - II.Step 1 : around primary vertex
1125 //lTrkgPrimaryVtxRadius3D = TMath::Sqrt( lTrkgPrimaryVtxPos[0] * lTrkgPrimaryVtxPos[0] +
1126 // lTrkgPrimaryVtxPos[1] * lTrkgPrimaryVtxPos[1] +
1127 // lTrkgPrimaryVtxPos[2] * lTrkgPrimaryVtxPos[2] );
1129 //lBestPrimaryVtxRadius3D = TMath::Sqrt( lBestPrimaryVtxPos[0] * lBestPrimaryVtxPos[0] +
1130 // lBestPrimaryVtxPos[1] * lBestPrimaryVtxPos[1] +
1131 // lBestPrimaryVtxPos[2] * lBestPrimaryVtxPos[2] );
1133 // - II.Step 2 : Assigning the necessary variables for specific AliESDcascade data members (ESD)
1136 xi->ChangeMassHypothesis(lV0quality , 3312); // default working hypothesis : cascade = Xi- decay
1138 lEffMassXi = xi->GetEffMassXi();
1139 //lChi2Xi = xi->GetChi2Xi();
1140 lDcaXiDaughters = xi->GetDcaXiDaughters();
1141 lXiCosineOfPointingAngle = xi->GetCascadeCosineOfPointingAngle( lBestPrimaryVtxPos[0],
1142 lBestPrimaryVtxPos[1],
1143 lBestPrimaryVtxPos[2] );
1144 // Take care : the best available vertex should be used (like in AliCascadeVertexer)
1146 xi->GetXYZcascade( lPosXi[0], lPosXi[1], lPosXi[2] );
1147 lXiRadius = TMath::Sqrt( lPosXi[0]*lPosXi[0] + lPosXi[1]*lPosXi[1] );
1149 // - II.Step 3 : around the tracks : Bach + V0 (ESD)
1150 // ~ Necessary variables for ESDcascade data members coming from the ESDv0 part (inheritance)
1153 UInt_t lIdxPosXi = (UInt_t) TMath::Abs( xi->GetPindex() );
1154 UInt_t lIdxNegXi = (UInt_t) TMath::Abs( xi->GetNindex() );
1155 UInt_t lBachIdx = (UInt_t) TMath::Abs( xi->GetBindex() );
1156 // Care track label can be negative in MC production (linked with the track quality)
1157 // However = normally, not the case for track index ...
1159 // FIXME : rejection of a double use of a daughter track (nothing but just a crosscheck of what is done in the cascade vertexer)
1160 if(lBachIdx == lIdxNegXi) {
1161 AliWarning("Pb / Idx(Bach. track) = Idx(Neg. track) ... continue!"); continue;
1163 if(lBachIdx == lIdxPosXi) {
1164 AliWarning("Pb / Idx(Bach. track) = Idx(Pos. track) ... continue!"); continue;
1167 AliESDtrack *pTrackXi = lESDevent->GetTrack( lIdxPosXi );
1168 AliESDtrack *nTrackXi = lESDevent->GetTrack( lIdxNegXi );
1169 AliESDtrack *bachTrackXi = lESDevent->GetTrack( lBachIdx );
1171 if (!pTrackXi || !nTrackXi || !bachTrackXi ) {
1172 AliWarning("ERROR: Could not retrieve one of the 3 ESD daughter tracks of the cascade ...");
1176 fTreeCascVarPosEta = pTrackXi->Eta();
1177 fTreeCascVarNegEta = nTrackXi->Eta();
1178 fTreeCascVarBachEta = bachTrackXi->Eta();
1180 //Save shared clusters information
1181 fTreeCascVarNegSharedClusters = nTrackXi->GetTPCnclsS(0,159);
1182 fTreeCascVarPosSharedClusters = pTrackXi->GetTPCnclsS(0,159);
1183 fTreeCascVarBachSharedClusters = bachTrackXi->GetTPCnclsS(0,159);
1185 Double_t lBMom[3], lNMom[3], lPMom[3];
1186 xi->GetBPxPyPz( lBMom[0], lBMom[1], lBMom[2] );
1187 xi->GetPPxPyPz( lPMom[0], lPMom[1], lPMom[2] );
1188 xi->GetNPxPyPz( lNMom[0], lNMom[1], lNMom[2] );
1190 //Save all momentum information
1191 fTreeCascVarNegPx = lNMom[0];
1192 fTreeCascVarNegPy = lNMom[1];
1193 fTreeCascVarNegPz = lNMom[2];
1194 fTreeCascVarPosPx = lPMom[0];
1195 fTreeCascVarPosPy = lPMom[1];
1196 fTreeCascVarPosPz = lPMom[2];
1197 fTreeCascVarBachPx = lBMom[0];
1198 fTreeCascVarBachPy = lBMom[1];
1199 fTreeCascVarBachPz = lBMom[2];
1201 fTreeCascVarBachTransMom = TMath::Sqrt( lBMom[0]*lBMom[0] + lBMom[1]*lBMom[1] );
1202 fTreeCascVarPosTransMom = TMath::Sqrt( lPMom[0]*lPMom[0] + lPMom[1]*lPMom[1] );
1203 fTreeCascVarNegTransMom = TMath::Sqrt( lNMom[0]*lNMom[0] + lNMom[1]*lNMom[1] );
1205 //------------------------------------------------
1206 // TPC dEdx information
1207 //------------------------------------------------
1208 fTreeCascVarNegNSigmaPion = fPIDResponse->NumberOfSigmasTPC( nTrackXi, AliPID::kPion );
1209 fTreeCascVarNegNSigmaProton = fPIDResponse->NumberOfSigmasTPC( nTrackXi, AliPID::kProton );
1210 fTreeCascVarPosNSigmaPion = fPIDResponse->NumberOfSigmasTPC( pTrackXi, AliPID::kPion );
1211 fTreeCascVarPosNSigmaProton = fPIDResponse->NumberOfSigmasTPC( pTrackXi, AliPID::kProton );
1212 fTreeCascVarBachNSigmaPion = fPIDResponse->NumberOfSigmasTPC( bachTrackXi, AliPID::kPion );
1213 fTreeCascVarBachNSigmaKaon = fPIDResponse->NumberOfSigmasTPC( bachTrackXi, AliPID::kKaon );
1215 //------------------------------------------------
1216 // TPC Number of clusters info
1217 // --- modified to save the smallest number
1218 // --- of TPC clusters for the 3 tracks
1219 //------------------------------------------------
1221 lPosTPCClusters = pTrackXi->GetTPCNcls();
1222 lNegTPCClusters = nTrackXi->GetTPCNcls();
1223 lBachTPCClusters = bachTrackXi->GetTPCNcls();
1225 // 1 - Poor quality related to TPCrefit
1226 ULong_t pStatus = pTrackXi->GetStatus();
1227 ULong_t nStatus = nTrackXi->GetStatus();
1228 ULong_t bachStatus = bachTrackXi->GetStatus();
1230 fTreeCascVarkITSRefitBachelor = kTRUE;
1231 fTreeCascVarkITSRefitNegative = kTRUE;
1232 fTreeCascVarkITSRefitPositive = kTRUE;
1234 if ((pStatus&AliESDtrack::kTPCrefit) == 0) { AliWarning("Pb / V0 Pos. track has no TPCrefit ... continue!"); continue; }
1235 if ((nStatus&AliESDtrack::kTPCrefit) == 0) { AliWarning("Pb / V0 Neg. track has no TPCrefit ... continue!"); continue; }
1236 if ((bachStatus&AliESDtrack::kTPCrefit) == 0) { AliWarning("Pb / Bach. track has no TPCrefit ... continue!"); continue; }
1238 //Extra Debug Information: booleans for ITS refit
1239 if ((pStatus&AliESDtrack::kITSrefit) == 0) { fTreeCascVarkITSRefitPositive = kFALSE; }
1240 if ((nStatus&AliESDtrack::kITSrefit) == 0) { fTreeCascVarkITSRefitNegative = kFALSE; }
1241 if ((bachStatus&AliESDtrack::kITSrefit) == 0) { fTreeCascVarkITSRefitBachelor = kFALSE; }
1243 // 2 - Poor quality related to TPC clusters: lowest cut of 70 clusters
1244 if(lPosTPCClusters < 70) { AliWarning("Pb / V0 Pos. track has less than 70 TPC clusters ... continue!"); continue; }
1245 if(lNegTPCClusters < 70) { AliWarning("Pb / V0 Neg. track has less than 70 TPC clusters ... continue!"); continue; }
1246 if(lBachTPCClusters < 70) { AliWarning("Pb / Bach. track has less than 70 TPC clusters ... continue!"); continue; }
1247 Int_t leastnumberofclusters = 1000;
1248 if( lPosTPCClusters < leastnumberofclusters ) leastnumberofclusters = lPosTPCClusters;
1249 if( lNegTPCClusters < leastnumberofclusters ) leastnumberofclusters = lNegTPCClusters;
1250 if( lBachTPCClusters < leastnumberofclusters ) leastnumberofclusters = lBachTPCClusters;
1252 lInvMassLambdaAsCascDghter = xi->GetEffMass();
1253 // This value shouldn't change, whatever the working hyp. is : Xi-, Xi+, Omega-, Omega+
1254 lDcaV0DaughtersXi = xi->GetDcaV0Daughters();
1255 //lV0Chi2Xi = xi->GetChi2V0();
1257 lV0CosineOfPointingAngleXi = xi->GetV0CosineOfPointingAngle( lBestPrimaryVtxPos[0],
1258 lBestPrimaryVtxPos[1],
1259 lBestPrimaryVtxPos[2] );
1260 //Modification: V0 CosPA wrt to Cascade decay vertex
1261 lV0CosineOfPointingAngleXiSpecial = xi->GetV0CosineOfPointingAngle( lPosXi[0],
1265 lDcaV0ToPrimVertexXi = xi->GetD( lBestPrimaryVtxPos[0],
1266 lBestPrimaryVtxPos[1],
1267 lBestPrimaryVtxPos[2] );
1269 lDcaBachToPrimVertexXi = TMath::Abs( bachTrackXi->GetD( lBestPrimaryVtxPos[0],
1270 lBestPrimaryVtxPos[1],
1272 // Note : AliExternalTrackParam::GetD returns an algebraic value ...
1274 xi->GetXYZ( lPosV0Xi[0], lPosV0Xi[1], lPosV0Xi[2] );
1275 lV0RadiusXi = TMath::Sqrt( lPosV0Xi[0]*lPosV0Xi[0] + lPosV0Xi[1]*lPosV0Xi[1] );
1277 lDcaPosToPrimVertexXi = TMath::Abs( pTrackXi ->GetD( lBestPrimaryVtxPos[0],
1278 lBestPrimaryVtxPos[1],
1281 lDcaNegToPrimVertexXi = TMath::Abs( nTrackXi ->GetD( lBestPrimaryVtxPos[0],
1282 lBestPrimaryVtxPos[1],
1285 // - II.Step 4 : around effective masses (ESD)
1286 // ~ change mass hypotheses to cover all the possibilities : Xi-/+, Omega -/+
1288 if( bachTrackXi->Charge() < 0 ) {
1290 xi->ChangeMassHypothesis(lV0quality , 3312);
1291 // Calculate the effective mass of the Xi- candidate.
1292 // pdg code 3312 = Xi-
1293 lInvMassXiMinus = xi->GetEffMassXi();
1296 xi->ChangeMassHypothesis(lV0quality , 3334);
1297 // Calculate the effective mass of the Xi- candidate.
1298 // pdg code 3334 = Omega-
1299 lInvMassOmegaMinus = xi->GetEffMassXi();
1302 xi->ChangeMassHypothesis(lV0quality , 3312); // Back to default hyp.
1303 }// end if negative bachelor
1306 if( bachTrackXi->Charge() > 0 ){
1308 xi->ChangeMassHypothesis(lV0quality , -3312);
1309 // Calculate the effective mass of the Xi+ candidate.
1310 // pdg code -3312 = Xi+
1311 lInvMassXiPlus = xi->GetEffMassXi();
1314 xi->ChangeMassHypothesis(lV0quality , -3334);
1315 // Calculate the effective mass of the Xi+ candidate.
1316 // pdg code -3334 = Omega+
1317 lInvMassOmegaPlus = xi->GetEffMassXi();
1320 xi->ChangeMassHypothesis(lV0quality , -3312); // Back to "default" hyp.
1321 }// end if positive bachelor
1322 // - II.Step 6 : extra info for QA (ESD)
1323 // miscellaneous pieces of info that may help regarding data quality assessment.
1326 xi->GetPxPyPz( lXiMomX, lXiMomY, lXiMomZ );
1327 lXiTransvMom = TMath::Sqrt( lXiMomX*lXiMomX + lXiMomY*lXiMomY );
1328 lXiTotMom = TMath::Sqrt( lXiMomX*lXiMomX + lXiMomY*lXiMomY + lXiMomZ*lXiMomZ );
1330 xi->GetBPxPyPz( lBachMomX, lBachMomY, lBachMomZ );
1331 //lBachTransvMom = TMath::Sqrt( lBachMomX*lBachMomX + lBachMomY*lBachMomY );
1332 //lBachTotMom = TMath::Sqrt( lBachMomX*lBachMomX + lBachMomY*lBachMomY + lBachMomZ*lBachMomZ );
1334 lChargeXi = xi->Charge();
1336 //lV0toXiCosineOfPointingAngle = xi->GetV0CosineOfPointingAngle( lPosXi[0], lPosXi[1], lPosXi[2] );
1338 lRapXi = xi->RapXi();
1339 lRapOmega = xi->RapOmega();
1341 //lTheta = xi->Theta() *180.0/TMath::Pi();
1342 //lPhi = xi->Phi() *180.0/TMath::Pi();
1343 //lAlphaXi = xi->AlphaXi();
1344 //lPtArmXi = xi->PtArmXi();
1346 //------------------------------------------------
1347 // Set Variables for adding to tree
1348 //------------------------------------------------
1350 /* 1*/ fTreeCascVarCharge = lChargeXi;
1351 /* 2*/ if(lInvMassXiMinus!=0) fTreeCascVarMassAsXi = lInvMassXiMinus;
1352 /* 2*/ if(lInvMassXiPlus!=0) fTreeCascVarMassAsXi = lInvMassXiPlus;
1353 /* 3*/ if(lInvMassOmegaMinus!=0) fTreeCascVarMassAsOmega = lInvMassOmegaMinus;
1354 /* 3*/ if(lInvMassOmegaPlus!=0) fTreeCascVarMassAsOmega = lInvMassOmegaPlus;
1355 /* 4*/ fTreeCascVarPt = lXiTransvMom;
1356 /* 4*/ fTreeCascVarPtMC = lXiTransvMomMC;
1357 /* 5*/ fTreeCascVarRapXi = lRapXi ;
1358 /* 5*/ fTreeCascVarRapMC = lRapMC ;
1359 /* 6*/ fTreeCascVarRapOmega = lRapOmega ;
1360 /* 7*/ fTreeCascVarDCACascDaughters = lDcaXiDaughters;
1361 /* 8*/ fTreeCascVarDCABachToPrimVtx = lDcaBachToPrimVertexXi;
1362 /* 9*/ fTreeCascVarDCAV0Daughters = lDcaV0DaughtersXi;
1363 /*10*/ fTreeCascVarDCAV0ToPrimVtx = lDcaV0ToPrimVertexXi;
1364 /*11*/ fTreeCascVarDCAPosToPrimVtx = lDcaPosToPrimVertexXi;
1365 /*12*/ fTreeCascVarDCANegToPrimVtx = lDcaNegToPrimVertexXi;
1366 /*13*/ fTreeCascVarCascCosPointingAngle = lXiCosineOfPointingAngle;
1367 /*14*/ fTreeCascVarCascRadius = lXiRadius;
1368 /*15*/ fTreeCascVarV0Mass = lInvMassLambdaAsCascDghter;
1369 /*16*/ fTreeCascVarV0CosPointingAngle = lV0CosineOfPointingAngleXi;
1370 /*16*/ fTreeCascVarV0CosPointingAngleSpecial = lV0CosineOfPointingAngleXiSpecial;
1371 /*17*/ fTreeCascVarV0Radius = lV0RadiusXi;
1372 /*20*/ fTreeCascVarLeastNbrClusters = leastnumberofclusters;
1373 /*21*/ fTreeCascVarMultiplicity = lMultiplicity; //multiplicity, whatever that may be
1375 /*23*/ fTreeCascVarDistOverTotMom = TMath::Sqrt(
1376 TMath::Power( lPosXi[0] - lBestPrimaryVtxPos[0] , 2) +
1377 TMath::Power( lPosXi[1] - lBestPrimaryVtxPos[1] , 2) +
1378 TMath::Power( lPosXi[2] - lBestPrimaryVtxPos[2] , 2)
1380 /*23*/ fTreeCascVarDistOverTotMom /= (lXiTotMom+1e-13);
1382 //All vars not specified here: specified elsewhere!
1384 //------------------------------------------------
1386 //------------------------------------------------
1388 // The conditional is meant to decrease excessive
1389 // memory usage! Be careful when loosening the
1392 //Xi Mass window: 150MeV wide
1393 //Omega mass window: 150MeV wide
1395 if( (fTreeCascVarMassAsXi<1.32+0.075&&fTreeCascVarMassAsXi>1.32-0.075) ||
1396 (fTreeCascVarMassAsOmega<1.68+0.075&&fTreeCascVarMassAsOmega>1.68-0.075) ){
1397 fTreeCascade->Fill();
1400 //------------------------------------------------
1402 //------------------------------------------------
1404 }// end of the Cascade loop (ESD or AOD)
1406 // Post output data.
1407 PostData(1, fListHist);
1408 PostData(2, fTreeCascade);
1411 //________________________________________________________________________
1412 void AliAnalysisTaskExtractCascade::Terminate(Option_t *)
1414 // Draw result to the screen
1415 // Called once at the end of the query
1417 TList *cRetrievedList = 0x0;
1418 cRetrievedList = (TList*)GetOutputData(1);
1419 if(!cRetrievedList){
1420 Printf("ERROR - AliAnalysisTaskExtractCascade : ouput data container list not available\n");
1424 fHistV0MultiplicityForTrigEvt = dynamic_cast<TH1F*> ( cRetrievedList->FindObject("fHistV0MultiplicityForTrigEvt") );
1425 if (!fHistV0MultiplicityForTrigEvt) {
1426 Printf("ERROR - AliAnalysisTaskExtractCascade : fHistV0MultiplicityForTrigEvt not available");
1430 TCanvas *canCheck = new TCanvas("AliAnalysisTaskExtractCascade","V0 Multiplicity",10,10,510,510);
1431 canCheck->cd(1)->SetLogy();
1433 fHistV0MultiplicityForTrigEvt->SetMarkerStyle(22);
1434 fHistV0MultiplicityForTrigEvt->DrawCopy("E");
1437 //----------------------------------------------------------------------------
1439 Double_t AliAnalysisTaskExtractCascade::MyRapidity(Double_t rE, Double_t rPz) const
1441 // Local calculation for rapidity
1442 Double_t ReturnValue = -100;
1443 if( (rE-rPz+1.e-13) != 0 && (rE+rPz) != 0 ){
1444 ReturnValue = 0.5*TMath::Log((rE+rPz)/(rE-rPz+1.e-13));