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.
23 // --- Adapted to look for lambdas as well, using code from
24 // AliAnalysisTaskCheckPerformanceStrange.cxx
26 // --- Algorithm Description
27 // 1. Loop over primaries in stack to acquire generated charged Xi
28 // 2. Loop over stack to find V0s, fill TH3Fs "PrimRawPt"s for Efficiency
29 // 3. Perform Physics Selection
30 // 4. Perform Primary Vertex |z|<10cm selection
31 // 5. Perform Primary Vertex NoTPCOnly vertexing selection (>0 contrib.)
32 // 6. Perform Pileup Rejection
34 // 7a. Fill TH3Fs "PrimAnalysisPt" for control purposes only
35 // 7b. Fill TTree object with V0 information, candidates
37 // Please Report Any Bugs!
39 // --- David Dobrigkeit Chinellato
40 // (david.chinellato@gmail.com)
42 // +-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
48 //class AliMCEventHandler;
57 #include <Riostream.h>
63 #include "THnSparse.h"
70 #include "AliESDEvent.h"
71 #include "AliAODEvent.h"
72 #include "AliV0vertexer.h"
73 #include "AliCascadeVertexer.h"
74 #include "AliESDpid.h"
75 #include "AliESDtrack.h"
76 #include "AliESDtrackCuts.h"
77 #include "AliInputEventHandler.h"
78 #include "AliAnalysisManager.h"
79 #include "AliMCEventHandler.h"
80 #include "AliMCEvent.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"
91 #include "AliAnalysisTaskExtractPerformanceV0.h"
96 ClassImp(AliAnalysisTaskExtractPerformanceV0)
98 AliAnalysisTaskExtractPerformanceV0::AliAnalysisTaskExtractPerformanceV0()
99 : AliAnalysisTaskSE(), fListHistV0(0), fTree(0), fPIDResponse(0), fESDtrackCuts(0),
100 fkIsNuclear ( kFALSE ),
101 fkSwitchINT7 ( kFALSE ),
102 fkUseOnTheFly ( kFALSE ),
103 fkTakeAllTracks ( kFALSE ),
104 fpArapidityShift ( 0.465 ),
105 fCentralityEstimator("V0M"),
106 //------------------------------------------------
109 fTreeVariablePrimaryStatus(0),
110 fTreeVariablePrimaryStatusMother(0),
111 fTreeVariableChi2V0(0),
112 fTreeVariableDcaV0Daughters(0),
113 fTreeVariableDcaV0ToPrimVertex(0),
114 fTreeVariableDcaPosToPrimVertex(0),
115 fTreeVariableDcaNegToPrimVertex(0),
116 fTreeVariableV0CosineOfPointingAngle(0),
117 fTreeVariableV0Radius(0),
119 fTreeVariablePtMC(0),
120 fTreeVariableRapK0Short(0),
121 fTreeVariableRapLambda(0),
122 fTreeVariableRapMC(0),
123 fTreeVariableInvMassK0s(0),
124 fTreeVariableInvMassLambda(0),
125 fTreeVariableInvMassAntiLambda(0),
126 fTreeVariableAlphaV0(0),
127 fTreeVariablePtArmV0(0),
128 fTreeVariableNegTotMomentum(0),
129 fTreeVariablePosTotMomentum(0),
130 fTreeVariableNegTransvMomentum(0),
131 fTreeVariablePosTransvMomentum(0),
132 fTreeVariableNegTransvMomentumMC(0),
133 fTreeVariablePosTransvMomentumMC(0),
135 fTreeVariableNSigmasPosProton(0),
136 fTreeVariableNSigmasPosPion(0),
137 fTreeVariableNSigmasNegProton(0),
138 fTreeVariableNSigmasNegPion(0),
140 fTreeVariablePtMother(0),
141 fTreeVariableV0CreationRadius(0),
143 fTreeVariablePIDPositive(0),
144 fTreeVariablePIDNegative(0),
145 fTreeVariablePIDMother(0),
146 fTreeVariableIndexStatus(0),
147 fTreeVariableIndexStatusMother(0),
149 fTreeVariableRunNumber(0),
150 fTreeVariableEventNumber(0),
152 fTreeVariableDistOverTotMom(0),
154 fTreeVariablePosEta(0),
155 fTreeVariableNegEta(0),
157 fTreeVariableVertexZ(0),
159 fTreeVariableLeastNbrCrossedRows(0),
160 fTreeVariableLeastRatioCrossedRowsOverFindable(0),
161 fTreeVariableMultiplicity(0),
162 fTreeVariableMultiplicityMC(0),
168 fTreeVariableV0Px(0),
169 fTreeVariableV0Py(0),
170 fTreeVariableV0Pz(0),
172 fTreeVariableMCV0x(0),
173 fTreeVariableMCV0y(0),
174 fTreeVariableMCV0z(0),
176 fTreeVariableMCV0Px(0),
177 fTreeVariableMCV0Py(0),
178 fTreeVariableMCV0Pz(0),
184 fTreeVariableMCPVx(0),
185 fTreeVariableMCPVy(0),
186 fTreeVariableMCPVz(0),
188 fTreeVariableIsNonInjected(0),
190 //------------------------------------------------
192 // --- Filled on an Event-by-event basis
193 //------------------------------------------------
194 fHistV0MultiplicityBeforeTrigSel(0),
195 fHistV0MultiplicityForTrigEvt(0),
196 fHistV0MultiplicityForSelEvt(0),
197 fHistV0MultiplicityForSelEvtNoTPCOnly(0),
198 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup(0),
199 fHistMultiplicityBeforeTrigSel(0),
200 fHistMultiplicityForTrigEvt(0),
201 fHistMultiplicity(0),
202 fHistMultiplicityNoTPCOnly(0),
203 fHistMultiplicityNoTPCOnlyNoPileup(0),
205 f2dHistMultiplicityVsTrueBeforeTrigSel(0),
206 f2dHistMultiplicityVsTrueForTrigEvt(0),
207 f2dHistMultiplicityVsTrue(0),
208 f2dHistMultiplicityVsTrueNoTPCOnly(0),
209 f2dHistMultiplicityVsTrueNoTPCOnlyNoPileup(0),
211 //Raw Data for Vertex Z position estimator change
212 f2dHistMultiplicityVsVertexZBeforeTrigSel(0),
213 f2dHistMultiplicityVsVertexZForTrigEvt(0),
214 f2dHistMultiplicityVsVertexZ(0),
215 f2dHistMultiplicityVsVertexZNoTPCOnly(0),
216 f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup(0),
218 fHistGenVertexZBeforeTrigSel(0),
219 fHistGenVertexZForTrigEvt(0),
221 fHistGenVertexZNoTPCOnly(0),
222 fHistGenVertexZNoTPCOnlyNoPileup(0),
224 //------------------------------------------------
225 // PARTICLE HISTOGRAMS
226 // --- Filled on a Particle-by-Particle basis
227 //------------------------------------------------
228 f3dHistPrimAnalysisPtVsYVsMultLambda(0),
229 f3dHistPrimAnalysisPtVsYVsMultAntiLambda(0),
230 f3dHistPrimAnalysisPtVsYVsMultK0Short(0),
231 f3dHistPrimAnalysisPtVsYCMSVsMultLambda(0),
232 f3dHistPrimAnalysisPtVsYCMSVsMultAntiLambda(0),
233 f3dHistPrimAnalysisPtVsYCMSVsMultK0Short(0),
234 f3dHistPrimRawPtVsYVsMultLambda(0),
235 f3dHistPrimRawPtVsYVsMultAntiLambda(0),
236 f3dHistPrimRawPtVsYVsMultK0Short(0),
237 f3dHistPrimRawPtVsYCMSVsMultLambda(0),
238 f3dHistPrimRawPtVsYCMSVsMultAntiLambda(0),
239 f3dHistPrimRawPtVsYCMSVsMultK0Short(0),
240 f3dHistPrimRawPtVsYVsMultNonInjLambda(0),
241 f3dHistPrimRawPtVsYVsMultNonInjAntiLambda(0),
242 f3dHistPrimRawPtVsYVsMultNonInjK0Short(0),
243 f3dHistPrimRawPtVsYVsMultMCLambda(0),
244 f3dHistPrimRawPtVsYVsMultMCAntiLambda(0),
245 f3dHistPrimRawPtVsYVsMultMCK0Short(0),
246 f3dHistPrimRawPtVsYVsVertexZLambda(0),
247 f3dHistPrimRawPtVsYVsVertexZAntiLambda(0),
248 f3dHistPrimRawPtVsYVsVertexZK0Short(0),
249 f3dHistPrimCloseToPVPtVsYVsMultLambda(0),
250 f3dHistPrimCloseToPVPtVsYVsMultAntiLambda(0),
251 f3dHistPrimCloseToPVPtVsYVsMultK0Short(0),
252 f3dHistPrimRawPtVsYVsDecayLengthLambda(0),
253 f3dHistPrimRawPtVsYVsDecayLengthAntiLambda(0),
254 f3dHistPrimRawPtVsYVsDecayLengthK0Short(0),
255 f3dHistGenPtVsYVsMultXiMinus(0),
256 f3dHistGenPtVsYVsMultXiPlus(0),
257 f3dHistGenPtVsYVsMultOmegaMinus(0),
258 f3dHistGenPtVsYVsMultOmegaPlus(0),
259 f3dHistGenSelectedPtVsYVsMultXiMinus(0),
260 f3dHistGenSelectedPtVsYVsMultXiPlus(0),
261 f3dHistGenSelectedPtVsYVsMultOmegaMinus(0),
262 f3dHistGenSelectedPtVsYVsMultOmegaPlus(0),
263 f3dHistGenPtVsYCMSVsMultXiMinus(0),
264 f3dHistGenPtVsYCMSVsMultXiPlus(0),
265 f3dHistGenPtVsYCMSVsMultOmegaMinus(0),
266 f3dHistGenPtVsYCMSVsMultOmegaPlus(0),
267 f3dHistGenSelectedPtVsYCMSVsMultXiMinus(0),
268 f3dHistGenSelectedPtVsYCMSVsMultXiPlus(0),
269 f3dHistGenSelectedPtVsYCMSVsMultOmegaMinus(0),
270 f3dHistGenSelectedPtVsYCMSVsMultOmegaPlus(0),
277 fHistPVxAnalysisHasHighPtLambda(0),
278 fHistPVyAnalysisHasHighPtLambda(0),
279 fHistPVzAnalysisHasHighPtLambda(0),
280 fHistSwappedV0Counter(0)
285 AliAnalysisTaskExtractPerformanceV0::AliAnalysisTaskExtractPerformanceV0(const char *name)
286 : AliAnalysisTaskSE(name), fListHistV0(0), fTree(0), fPIDResponse(0), fESDtrackCuts(0),
287 fkIsNuclear ( kFALSE ),
288 fkSwitchINT7 ( kFALSE ),
289 fkUseOnTheFly ( kFALSE ),
290 fkTakeAllTracks ( kFALSE ),
291 fpArapidityShift ( 0.465 ),
292 fCentralityEstimator("V0M"),
293 //------------------------------------------------
296 fTreeVariablePrimaryStatus(0),
297 fTreeVariablePrimaryStatusMother(0),
298 fTreeVariableChi2V0(0),
299 fTreeVariableDcaV0Daughters(0),
300 fTreeVariableDcaV0ToPrimVertex(0),
301 fTreeVariableDcaPosToPrimVertex(0),
302 fTreeVariableDcaNegToPrimVertex(0),
303 fTreeVariableV0CosineOfPointingAngle(0),
304 fTreeVariableV0Radius(0),
306 fTreeVariablePtMC(0),
307 fTreeVariableRapK0Short(0),
308 fTreeVariableRapLambda(0),
309 fTreeVariableRapMC(0),
310 fTreeVariableInvMassK0s(0),
311 fTreeVariableInvMassLambda(0),
312 fTreeVariableInvMassAntiLambda(0),
313 fTreeVariableAlphaV0(0),
314 fTreeVariablePtArmV0(0),
315 fTreeVariableNegTotMomentum(0),
316 fTreeVariablePosTotMomentum(0),
317 fTreeVariableNegTransvMomentum(0),
318 fTreeVariablePosTransvMomentum(0),
319 fTreeVariableNegTransvMomentumMC(0),
320 fTreeVariablePosTransvMomentumMC(0),
322 fTreeVariableNSigmasPosProton(0),
323 fTreeVariableNSigmasPosPion(0),
324 fTreeVariableNSigmasNegProton(0),
325 fTreeVariableNSigmasNegPion(0),
327 fTreeVariablePtMother(0),
328 fTreeVariableV0CreationRadius(0),
330 fTreeVariablePIDPositive(0),
331 fTreeVariablePIDNegative(0),
332 fTreeVariablePIDMother(0),
333 fTreeVariableIndexStatus(0),
334 fTreeVariableIndexStatusMother(0),
336 fTreeVariableRunNumber(0),
337 fTreeVariableEventNumber(0),
339 fTreeVariableDistOverTotMom(0),
341 fTreeVariablePosEta(0),
342 fTreeVariableNegEta(0),
344 fTreeVariableVertexZ(0),
346 fTreeVariableLeastNbrCrossedRows(0),
347 fTreeVariableLeastRatioCrossedRowsOverFindable(0),
348 fTreeVariableMultiplicity(0),
349 fTreeVariableMultiplicityMC(0),
355 fTreeVariableV0Px(0),
356 fTreeVariableV0Py(0),
357 fTreeVariableV0Pz(0),
359 fTreeVariableMCV0x(0),
360 fTreeVariableMCV0y(0),
361 fTreeVariableMCV0z(0),
363 fTreeVariableMCV0Px(0),
364 fTreeVariableMCV0Py(0),
365 fTreeVariableMCV0Pz(0),
371 fTreeVariableMCPVx(0),
372 fTreeVariableMCPVy(0),
373 fTreeVariableMCPVz(0),
375 fTreeVariableIsNonInjected(0),
378 //------------------------------------------------
380 // --- Filled on an Event-by-event basis
381 //------------------------------------------------
382 fHistV0MultiplicityBeforeTrigSel(0),
383 fHistV0MultiplicityForTrigEvt(0),
384 fHistV0MultiplicityForSelEvt(0),
385 fHistV0MultiplicityForSelEvtNoTPCOnly(0),
386 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup(0),
387 fHistMultiplicityBeforeTrigSel(0),
388 fHistMultiplicityForTrigEvt(0),
389 fHistMultiplicity(0),
390 fHistMultiplicityNoTPCOnly(0),
391 fHistMultiplicityNoTPCOnlyNoPileup(0),
393 f2dHistMultiplicityVsTrueBeforeTrigSel(0),
394 f2dHistMultiplicityVsTrueForTrigEvt(0),
395 f2dHistMultiplicityVsTrue(0),
396 f2dHistMultiplicityVsTrueNoTPCOnly(0),
397 f2dHistMultiplicityVsTrueNoTPCOnlyNoPileup(0),
399 //Raw Data for Vertex Z position estimator change
400 f2dHistMultiplicityVsVertexZBeforeTrigSel(0),
401 f2dHistMultiplicityVsVertexZForTrigEvt(0),
402 f2dHistMultiplicityVsVertexZ(0),
403 f2dHistMultiplicityVsVertexZNoTPCOnly(0),
404 f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup(0),
406 fHistGenVertexZBeforeTrigSel(0),
407 fHistGenVertexZForTrigEvt(0),
409 fHistGenVertexZNoTPCOnly(0),
410 fHistGenVertexZNoTPCOnlyNoPileup(0),
412 //------------------------------------------------
413 // PARTICLE HISTOGRAMS
414 // --- Filled on a Particle-by-Particle basis
415 //------------------------------------------------
416 f3dHistPrimAnalysisPtVsYVsMultLambda(0),
417 f3dHistPrimAnalysisPtVsYVsMultAntiLambda(0),
418 f3dHistPrimAnalysisPtVsYVsMultK0Short(0),
419 f3dHistPrimAnalysisPtVsYCMSVsMultLambda(0),
420 f3dHistPrimAnalysisPtVsYCMSVsMultAntiLambda(0),
421 f3dHistPrimAnalysisPtVsYCMSVsMultK0Short(0),
422 f3dHistPrimRawPtVsYVsMultLambda(0),
423 f3dHistPrimRawPtVsYVsMultAntiLambda(0),
424 f3dHistPrimRawPtVsYVsMultK0Short(0),
425 f3dHistPrimRawPtVsYCMSVsMultLambda(0),
426 f3dHistPrimRawPtVsYCMSVsMultAntiLambda(0),
427 f3dHistPrimRawPtVsYCMSVsMultK0Short(0),
428 f3dHistPrimRawPtVsYVsMultNonInjLambda(0),
429 f3dHistPrimRawPtVsYVsMultNonInjAntiLambda(0),
430 f3dHistPrimRawPtVsYVsMultNonInjK0Short(0),
431 f3dHistPrimRawPtVsYVsMultMCLambda(0),
432 f3dHistPrimRawPtVsYVsMultMCAntiLambda(0),
433 f3dHistPrimRawPtVsYVsMultMCK0Short(0),
434 f3dHistPrimRawPtVsYVsVertexZLambda(0),
435 f3dHistPrimRawPtVsYVsVertexZAntiLambda(0),
436 f3dHistPrimRawPtVsYVsVertexZK0Short(0),
437 f3dHistPrimCloseToPVPtVsYVsMultLambda(0),
438 f3dHistPrimCloseToPVPtVsYVsMultAntiLambda(0),
439 f3dHistPrimCloseToPVPtVsYVsMultK0Short(0),
440 f3dHistPrimRawPtVsYVsDecayLengthLambda(0),
441 f3dHistPrimRawPtVsYVsDecayLengthAntiLambda(0),
442 f3dHistPrimRawPtVsYVsDecayLengthK0Short(0),
443 f3dHistGenPtVsYVsMultXiMinus(0),
444 f3dHistGenPtVsYVsMultXiPlus(0),
445 f3dHistGenPtVsYVsMultOmegaMinus(0),
446 f3dHistGenPtVsYVsMultOmegaPlus(0),
447 f3dHistGenSelectedPtVsYVsMultXiMinus(0),
448 f3dHistGenSelectedPtVsYVsMultXiPlus(0),
449 f3dHistGenSelectedPtVsYVsMultOmegaMinus(0),
450 f3dHistGenSelectedPtVsYVsMultOmegaPlus(0),
451 f3dHistGenPtVsYCMSVsMultXiMinus(0),
452 f3dHistGenPtVsYCMSVsMultXiPlus(0),
453 f3dHistGenPtVsYCMSVsMultOmegaMinus(0),
454 f3dHistGenPtVsYCMSVsMultOmegaPlus(0),
455 f3dHistGenSelectedPtVsYCMSVsMultXiMinus(0),
456 f3dHistGenSelectedPtVsYCMSVsMultXiPlus(0),
457 f3dHistGenSelectedPtVsYCMSVsMultOmegaMinus(0),
458 f3dHistGenSelectedPtVsYCMSVsMultOmegaPlus(0),
465 fHistPVxAnalysisHasHighPtLambda(0),
466 fHistPVyAnalysisHasHighPtLambda(0),
467 fHistPVzAnalysisHasHighPtLambda(0),
468 fHistSwappedV0Counter(0)
471 // Output slot #0 writes into a TList container (Cascade)
472 DefineOutput(1, TList::Class());
473 DefineOutput(2, TTree::Class());
477 AliAnalysisTaskExtractPerformanceV0::~AliAnalysisTaskExtractPerformanceV0()
479 //------------------------------------------------
481 //------------------------------------------------
491 //cleanup esd track cuts object too...
493 delete fESDtrackCuts;
498 //________________________________________________________________________
499 void AliAnalysisTaskExtractPerformanceV0::UserCreateOutputObjects()
505 //------------------------------------------------
507 fTree = new TTree("fTree","V0Candidates");
509 //------------------------------------------------
510 // fTree Branch definitions - V0 Tree
511 //------------------------------------------------
513 //-----------BASIC-INFO---------------------------
514 /* 1*/ fTree->Branch("fTreeVariablePrimaryStatus",&fTreeVariablePrimaryStatus,"fTreeVariablePrimaryStatus/I");
515 /* 1*/ fTree->Branch("fTreeVariablePrimaryStatusMother",&fTreeVariablePrimaryStatusMother,"fTreeVariablePrimaryStatusMother/I");
516 /* 2*/ fTree->Branch("fTreeVariableChi2V0",&fTreeVariableChi2V0,"Chi2V0/F");
517 /* 3*/ fTree->Branch("fTreeVariableDcaV0Daughters",&fTreeVariableDcaV0Daughters,"fTreeVariableDcaV0Daughters/F");
518 /* 4*/ fTree->Branch("fTreeVariableDcaPosToPrimVertex",&fTreeVariableDcaPosToPrimVertex,"fTreeVariableDcaPosToPrimVertex/F");
519 /* 5*/ fTree->Branch("fTreeVariableDcaNegToPrimVertex",&fTreeVariableDcaNegToPrimVertex,"fTreeVariableDcaNegToPrimVertex/F");
520 /* 6*/ fTree->Branch("fTreeVariableV0Radius",&fTreeVariableV0Radius,"fTreeVariableV0Radius/F");
521 /* 7*/ fTree->Branch("fTreeVariablePt",&fTreeVariablePt,"fTreeVariablePt/F");
522 /* 7*/ fTree->Branch("fTreeVariablePtMC",&fTreeVariablePtMC,"fTreeVariablePtMC/F");
523 /* 8*/ fTree->Branch("fTreeVariableRapK0Short",&fTreeVariableRapK0Short,"fTreeVariableRapK0Short/F");
524 /* 9*/ fTree->Branch("fTreeVariableRapLambda",&fTreeVariableRapLambda,"fTreeVariableRapLambda/F");
525 /*10*/ fTree->Branch("fTreeVariableRapMC",&fTreeVariableRapMC,"fTreeVariableRapMC/F");
526 /*11*/ fTree->Branch("fTreeVariableInvMassK0s",&fTreeVariableInvMassK0s,"fTreeVariableInvMassK0s/F");
527 /*12*/ fTree->Branch("fTreeVariableInvMassLambda",&fTreeVariableInvMassLambda,"fTreeVariableInvMassLambda/F");
528 /*13*/ fTree->Branch("fTreeVariableInvMassAntiLambda",&fTreeVariableInvMassAntiLambda,"fTreeVariableInvMassAntiLambda/F");
529 /*14*/ fTree->Branch("fTreeVariableAlphaV0",&fTreeVariableAlphaV0,"fTreeVariableAlphaV0/F");
530 /*15*/ fTree->Branch("fTreeVariablePtArmV0",&fTreeVariablePtArmV0,"fTreeVariablePtArmV0/F");
531 /*16*/ fTree->Branch("fTreeVariableNegTransvMomentum",&fTreeVariableNegTransvMomentum,"fTreeVariableNegTransvMomentum/F");
532 /*17*/ fTree->Branch("fTreeVariablePosTransvMomentum",&fTreeVariablePosTransvMomentum,"fTreeVariablePosTransvMomentum/F");
533 /*18*/ fTree->Branch("fTreeVariableNegTransvMomentumMC",&fTreeVariableNegTransvMomentumMC,"fTreeVariableNegTransvMomentumMC/F");
534 /*19*/ fTree->Branch("fTreeVariablePosTransvMomentumMC",&fTreeVariablePosTransvMomentumMC,"fTreeVariablePosTransvMomentumMC/F");
535 /*20*/ fTree->Branch("fTreeVariableLeastNbrCrossedRows",&fTreeVariableLeastNbrCrossedRows,"fTreeVariableLeastNbrCrossedRows/I");
536 /*21*/ fTree->Branch("fTreeVariableLeastRatioCrossedRowsOverFindable",&fTreeVariableLeastRatioCrossedRowsOverFindable,"fTreeVariableLeastRatioCrossedRowsOverFindable/F");
537 /*22*/ fTree->Branch("fTreeVariablePID",&fTreeVariablePID,"fTreeVariablePID/I");
538 /*23*/ fTree->Branch("fTreeVariablePIDPositive",&fTreeVariablePIDPositive,"fTreeVariablePIDPositive/I");
539 /*24*/ fTree->Branch("fTreeVariablePIDNegative",&fTreeVariablePIDNegative,"fTreeVariablePIDNegative/I");
540 /*25*/ fTree->Branch("fTreeVariablePIDMother",&fTreeVariablePIDMother,"fTreeVariablePIDMother/I");
541 /*26*/ fTree->Branch("fTreeVariablePtXiMother",&fTreeVariablePtMother,"fTreeVariablePtMother/F");
542 /*27*/ fTree->Branch("fTreeVariableV0CosineOfPointingAngle",&fTreeVariableV0CosineOfPointingAngle,"fTreeVariableV0CosineOfPointingAngle/F");
543 //-----------MULTIPLICITY-INFO--------------------
544 /*28*/ fTree->Branch("fTreeVariableMultiplicity",&fTreeVariableMultiplicity,"fTreeVariableMultiplicity/I");
545 /*28*/ fTree->Branch("fTreeVariableMultiplicityMC",&fTreeVariableMultiplicityMC,"fTreeVariableMultiplicityMC/I");
546 //------------------------------------------------
547 /*29*/ fTree->Branch("fTreeVariableDistOverTotMom",&fTreeVariableDistOverTotMom,"fTreeVariableDistOverTotMom/F");
548 /*30*/ fTree->Branch("fTreeVariableNSigmasPosProton",&fTreeVariableNSigmasPosProton,"fTreeVariableNSigmasPosProton/F");
549 /*31*/ fTree->Branch("fTreeVariableNSigmasPosPion",&fTreeVariableNSigmasPosPion,"fTreeVariableNSigmasPosPion/F");
550 /*32*/ fTree->Branch("fTreeVariableNSigmasNegProton",&fTreeVariableNSigmasNegProton,"fTreeVariableNSigmasNegProton/F");
551 /*33*/ fTree->Branch("fTreeVariableNSigmasNegPion",&fTreeVariableNSigmasNegPion,"fTreeVariableNSigmasNegPion/F");
552 //------------------------------------------------
553 /*34*/ fTree->Branch("fTreeVariableNegEta",&fTreeVariableNegEta,"fTreeVariableNegEta/F");
554 /*35*/ fTree->Branch("fTreeVariablePosEta",&fTreeVariablePosEta,"fTreeVariablePosEta/F");
555 /*36*/ fTree->Branch("fTreeVariableV0CreationRadius",&fTreeVariableV0CreationRadius,"fTreeVariableV0CreationRadius/F");
556 /*37*/ fTree->Branch("fTreeVariableIndexStatus",&fTreeVariableIndexStatus,"fTreeVariableIndexStatus/I");
557 /*38*/ fTree->Branch("fTreeVariableIndexStatusMother",&fTreeVariableIndexStatusMother,"fTreeVariableIndexStatusMother/I");
559 /*39*/ fTree->Branch("fTreeVariableRunNumber",&fTreeVariableRunNumber,"fTreeVariableRunNumber/I");
560 /*40*/ fTree->Branch("fTreeVariableEventNumber",&fTreeVariableEventNumber,"fTreeVariableEventNumber/l");
562 /*34*/ fTree->Branch("fTreeVariableVertexZ",&fTreeVariableVertexZ,"fTreeVariableVertexZ/F");
564 //-----------FOR CTAU DEBUGGING: Full Phase Space + Decay Position Info
565 fTree->Branch("fTreeVariableV0x",&fTreeVariableV0x,"fTreeVariableV0x/F");
566 fTree->Branch("fTreeVariableV0y",&fTreeVariableV0y,"fTreeVariableV0y/F");
567 fTree->Branch("fTreeVariableV0z",&fTreeVariableV0z,"fTreeVariableV0z/F");
569 fTree->Branch("fTreeVariableV0Px",&fTreeVariableV0Px,"fTreeVariableV0Px/F");
570 fTree->Branch("fTreeVariableV0Py",&fTreeVariableV0Py,"fTreeVariableV0Py/F");
571 fTree->Branch("fTreeVariableV0Pz",&fTreeVariableV0Pz,"fTreeVariableV0Pz/F");
573 //-----------FOR CTAU DEBUGGING: Full Phase Space + Decay Position Info, perfect info from MC
574 fTree->Branch("fTreeVariableMCV0x",&fTreeVariableMCV0x,"fTreeVariableMCV0x/F");
575 fTree->Branch("fTreeVariableMCV0y",&fTreeVariableMCV0y,"fTreeVariableMCV0y/F");
576 fTree->Branch("fTreeVariableMCV0z",&fTreeVariableMCV0z,"fTreeVariableMCV0z/F");
578 fTree->Branch("fTreeVariableMCV0Px",&fTreeVariableMCV0Px,"fTreeVariableMCV0Px/F");
579 fTree->Branch("fTreeVariableMCV0Py",&fTreeVariableMCV0Py,"fTreeVariableMCV0Py/F");
580 fTree->Branch("fTreeVariableMCV0Pz",&fTreeVariableMCV0Pz,"fTreeVariableMCV0Pz/F");
582 //-----------FOR CTAU DEBUGGING: Primary vertex info
583 fTree->Branch("fTreeVariablePVx",&fTreeVariablePVx,"fTreeVariablePVx/F");
584 fTree->Branch("fTreeVariablePVy",&fTreeVariablePVy,"fTreeVariablePVy/F");
585 fTree->Branch("fTreeVariablePVz",&fTreeVariablePVz,"fTreeVariablePVz/F");
587 fTree->Branch("fTreeVariableMCPVx",&fTreeVariableMCPVx,"fTreeVariableMCPVx/F");
588 fTree->Branch("fTreeVariableMCPVy",&fTreeVariableMCPVy,"fTreeVariableMCPVy/F");
589 fTree->Branch("fTreeVariableMCPVz",&fTreeVariableMCPVz,"fTreeVariableMCPVz/F");
591 fTree->Branch("fTreeVariableIsNonInjected",&fTreeVariableIsNonInjected,"fTreeVariableIsNonInjected/O"); //O for bOOlean...
592 //------------------------------------------------
593 // Particle Identification Setup
594 //------------------------------------------------
596 AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
597 AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
598 fPIDResponse = inputHandler->GetPIDResponse();
602 if(! fESDtrackCuts ){
603 fESDtrackCuts = new AliESDtrackCuts();
606 //------------------------------------------------
607 // V0 Multiplicity Histograms
608 //------------------------------------------------
612 fListHistV0 = new TList();
613 fListHistV0->SetOwner(); // See http://root.cern.ch/root/html/TCollection.html#TCollection:SetOwner
616 if(! fHistV0MultiplicityBeforeTrigSel) {
617 fHistV0MultiplicityBeforeTrigSel = new TH1F("fHistV0MultiplicityBeforeTrigSel",
618 "V0s per event (before Trig. Sel.);Nbr of V0s/Evt;Events",
620 fListHistV0->Add(fHistV0MultiplicityBeforeTrigSel);
623 if(! fHistV0MultiplicityForTrigEvt) {
624 fHistV0MultiplicityForTrigEvt = new TH1F("fHistV0MultiplicityForTrigEvt",
625 "V0s per event (for triggered evt);Nbr of V0s/Evt;Events",
627 fListHistV0->Add(fHistV0MultiplicityForTrigEvt);
630 if(! fHistV0MultiplicityForSelEvt) {
631 fHistV0MultiplicityForSelEvt = new TH1F("fHistV0MultiplicityForSelEvt",
632 "V0s per event;Nbr of V0s/Evt;Events",
634 fListHistV0->Add(fHistV0MultiplicityForSelEvt);
637 if(! fHistV0MultiplicityForSelEvtNoTPCOnly) {
638 fHistV0MultiplicityForSelEvtNoTPCOnly = new TH1F("fHistV0MultiplicityForSelEvtNoTPCOnly",
639 "V0s per event;Nbr of V0s/Evt;Events",
641 fListHistV0->Add(fHistV0MultiplicityForSelEvtNoTPCOnly);
643 if(! fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup) {
644 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup = new TH1F("fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup",
645 "V0s per event;Nbr of V0s/Evt;Events",
647 fListHistV0->Add(fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup);
650 //------------------------------------------------
651 // Track Multiplicity Histograms
652 //------------------------------------------------
654 if(! fHistMultiplicityBeforeTrigSel) {
655 fHistMultiplicityBeforeTrigSel = new TH1F("fHistMultiplicityBeforeTrigSel",
656 "Tracks per event;Nbr of Tracks;Events",
658 fListHistV0->Add(fHistMultiplicityBeforeTrigSel);
660 if(! fHistMultiplicityForTrigEvt) {
661 fHistMultiplicityForTrigEvt = new TH1F("fHistMultiplicityForTrigEvt",
662 "Tracks per event;Nbr of Tracks;Events",
664 fListHistV0->Add(fHistMultiplicityForTrigEvt);
666 if(! fHistMultiplicity) {
667 fHistMultiplicity = new TH1F("fHistMultiplicity",
668 "Tracks per event;Nbr of Tracks;Events",
670 fListHistV0->Add(fHistMultiplicity);
672 if(! fHistMultiplicityNoTPCOnly) {
673 fHistMultiplicityNoTPCOnly = new TH1F("fHistMultiplicityNoTPCOnly",
674 "Tracks per event;Nbr of Tracks;Events",
676 fListHistV0->Add(fHistMultiplicityNoTPCOnly);
678 if(! fHistMultiplicityNoTPCOnlyNoPileup) {
679 fHistMultiplicityNoTPCOnlyNoPileup = new TH1F("fHistMultiplicityNoTPCOnlyNoPileup",
680 "Tracks per event;Nbr of Tracks;Events",
682 fListHistV0->Add(fHistMultiplicityNoTPCOnlyNoPileup);
685 //Raw Data for J/Psi paper Technique
686 //TH2F *f2dHistMultiplicityVsTrueBeforeTrigSel; //! multiplicity distribution
687 //TH2F *f2dHistMultiplicityVsTrueForTrigEvt; //! multiplicity distribution
688 //TH2F *f2dHistMultiplicityVsTrue; //! multiplicity distribution
689 //TH2F *f2dHistMultiplicityVsTrueNoTPCOnly; //! multiplicity distribution
690 //TH2F *f2dHistMultiplicityVsTrueNoTPCOnlyNoPileup; //! multiplicity distribution
692 if(! f2dHistMultiplicityVsTrueBeforeTrigSel) {
693 f2dHistMultiplicityVsTrueBeforeTrigSel = new TH2F("f2dHistMultiplicityVsTrueBeforeTrigSel",
694 "Tracks per event", 200, 0, 200, 200, 0, 200);
695 fListHistV0->Add(f2dHistMultiplicityVsTrueBeforeTrigSel);
697 if(! f2dHistMultiplicityVsTrueForTrigEvt) {
698 f2dHistMultiplicityVsTrueForTrigEvt = new TH2F("f2dHistMultiplicityVsTrueForTrigEvt",
699 "Tracks per event", 200, 0, 200, 200, 0, 200);
700 fListHistV0->Add(f2dHistMultiplicityVsTrueForTrigEvt);
702 if(! f2dHistMultiplicityVsTrue) {
703 f2dHistMultiplicityVsTrue = new TH2F("f2dHistMultiplicityVsTrue",
704 "Tracks per event", 200, 0, 200, 200, 0, 200);
705 fListHistV0->Add(f2dHistMultiplicityVsTrue);
707 if(! f2dHistMultiplicityVsTrueNoTPCOnly) {
708 f2dHistMultiplicityVsTrueNoTPCOnly = new TH2F("f2dHistMultiplicityVsTrueNoTPCOnly",
709 "Tracks per event", 200, 0, 200, 200, 0, 200);
710 fListHistV0->Add(f2dHistMultiplicityVsTrueNoTPCOnly);
712 if(! f2dHistMultiplicityVsTrueNoTPCOnlyNoPileup) {
713 f2dHistMultiplicityVsTrueNoTPCOnlyNoPileup = new TH2F("f2dHistMultiplicityVsTrueNoTPCOnlyNoPileup",
714 "Tracks per event", 200, 0, 200, 200, 0, 200);
715 fListHistV0->Add(f2dHistMultiplicityVsTrueNoTPCOnlyNoPileup);
719 //Raw Data for Vertex Z position estimator change
720 //TH2F *f2dHistMultiplicityVsVertexZBeforeTrigSel; //! multiplicity distribution
721 //TH2F *f2dHistMultiplicityVsVertexZForTrigEvt; //! multiplicity distribution
722 //TH2F *f2dHistMultiplicityVsVertexZ; //! multiplicity distribution
723 //TH2F *f2dHistMultiplicityVsVertexZNoTPCOnly; //! multiplicity distribution
724 //TH2F *f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup; //! multiplicity distribution
726 if(! f2dHistMultiplicityVsVertexZBeforeTrigSel) {
727 f2dHistMultiplicityVsVertexZBeforeTrigSel = new TH2F("f2dHistMultiplicityVsVertexZBeforeTrigSel",
728 "Tracks per event", 200, 0, 200,400, -20, 20);
729 fListHistV0->Add(f2dHistMultiplicityVsVertexZBeforeTrigSel);
731 if(! f2dHistMultiplicityVsVertexZForTrigEvt) {
732 f2dHistMultiplicityVsVertexZForTrigEvt = new TH2F("f2dHistMultiplicityVsVertexZForTrigEvt",
733 "Tracks per event", 200, 0, 200, 400, -20, 20);
734 fListHistV0->Add(f2dHistMultiplicityVsVertexZForTrigEvt);
736 if(! f2dHistMultiplicityVsVertexZ) {
737 f2dHistMultiplicityVsVertexZ = new TH2F("f2dHistMultiplicityVsVertexZ",
738 "Tracks per event", 200, 0, 200, 400, -20, 20);
739 fListHistV0->Add(f2dHistMultiplicityVsVertexZ);
741 if(! f2dHistMultiplicityVsVertexZNoTPCOnly) {
742 f2dHistMultiplicityVsVertexZNoTPCOnly = new TH2F("f2dHistMultiplicityVsVertexZNoTPCOnly",
743 "Tracks per event", 200, 0, 200, 400, -20, 20);
744 fListHistV0->Add(f2dHistMultiplicityVsVertexZNoTPCOnly);
746 if(! f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup) {
747 f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup = new TH2F("f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup",
748 "Tracks per event", 200, 0, 200, 400, -20, 20);
749 fListHistV0->Add(f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup);
754 // TH1F *fHistGenVertexZBeforeTrigSel; //! multiplicity distribution
755 // TH1F *fHistGenVertexZForTrigEvt; //! multiplicity distribution
756 // TH1F *fHistGenVertexZ; //! multiplicity distribution
757 // TH1F *fHistGenVertexZNoTPCOnly; //! multiplicity distribution
758 // TH1F *fHistGenVertexZNoTPCOnlyNoPileup; //! multiplicity distribution
760 if(! fHistGenVertexZBeforeTrigSel) {
761 fHistGenVertexZBeforeTrigSel = new TH1F("fHistGenVertexZBeforeTrigSel",
762 "PV z position;Nbr of Evts;z",
764 fListHistV0->Add(fHistGenVertexZBeforeTrigSel);
766 if(! fHistGenVertexZForTrigEvt) {
767 fHistGenVertexZForTrigEvt = new TH1F("fHistGenVertexZForTrigEvt",
768 "PV z position;Nbr of Evts;z",
770 fListHistV0->Add(fHistGenVertexZForTrigEvt);
772 if(! fHistGenVertexZ) {
773 fHistGenVertexZ = new TH1F("fHistGenVertexZ",
774 "PV z position;Nbr of Evts;z",
776 fListHistV0->Add(fHistGenVertexZ);
778 if(! fHistGenVertexZNoTPCOnly) {
779 fHistGenVertexZNoTPCOnly = new TH1F("fHistGenVertexZNoTPCOnly",
780 "PV z position;Nbr of Evts;z",
782 fListHistV0->Add(fHistGenVertexZNoTPCOnly);
784 if(! fHistGenVertexZNoTPCOnlyNoPileup) {
785 fHistGenVertexZNoTPCOnlyNoPileup = new TH1F("fHistGenVertexZNoTPCOnlyNoPileup",
786 "PV z position;Nbr of Evts;z",
788 fListHistV0->Add(fHistGenVertexZNoTPCOnlyNoPileup);
792 //------------------------------------------------
793 // Generated Particle Histograms
794 //------------------------------------------------
796 Int_t lCustomNBins = 200;
797 Double_t lCustomPtUpperLimit = 20;
798 Int_t lCustomNBinsMultiplicity = 100;
800 //----------------------------------
801 // Raw Generated (Pre-physics-selection)
802 //----------------------------------
804 //--- 3D Histo (Pt, Y, Multiplicity)
806 if(! f3dHistPrimRawPtVsYVsMultLambda) {
807 f3dHistPrimRawPtVsYVsMultLambda = new TH3F( "f3dHistPrimRawPtVsYVsMultLambda", "Pt_{lambda} Vs Y_{#Lambda} Vs Multiplicity; Pt_{lambda} (GeV/c); Y_{#Lambda} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
808 fListHistV0->Add(f3dHistPrimRawPtVsYVsMultLambda);
810 if(! f3dHistPrimRawPtVsYVsMultAntiLambda) {
811 f3dHistPrimRawPtVsYVsMultAntiLambda = new TH3F( "f3dHistPrimRawPtVsYVsMultAntiLambda", "Pt_{antilambda} Vs Y_{#Lambda} Vs Multiplicity; Pt_{antilambda} (GeV/c); Y_{#Lambda} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
812 fListHistV0->Add(f3dHistPrimRawPtVsYVsMultAntiLambda);
814 if(! f3dHistPrimRawPtVsYVsMultK0Short) {
815 f3dHistPrimRawPtVsYVsMultK0Short = new TH3F( "f3dHistPrimRawPtVsYVsMultK0Short", "Pt_{K0S} Vs Y_{K0S} Vs Multiplicity; Pt_{K0S} (GeV/c); Y_{K0S} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
816 fListHistV0->Add(f3dHistPrimRawPtVsYVsMultK0Short);
819 if(! f3dHistPrimRawPtVsYCMSVsMultLambda) {
820 f3dHistPrimRawPtVsYCMSVsMultLambda = new TH3F( "f3dHistPrimRawPtVsYCMSVsMultLambda", "Pt_{lambda} Vs Y_{#Lambda} Vs Multiplicity; Pt_{lambda} (GeV/c); Y_{#Lambda} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
821 fListHistV0->Add(f3dHistPrimRawPtVsYCMSVsMultLambda);
823 if(! f3dHistPrimRawPtVsYCMSVsMultAntiLambda) {
824 f3dHistPrimRawPtVsYCMSVsMultAntiLambda = new TH3F( "f3dHistPrimRawPtVsYCMSVsMultAntiLambda", "Pt_{antilambda} Vs Y_{#Lambda} Vs Multiplicity; Pt_{antilambda} (GeV/c); Y_{#Lambda} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
825 fListHistV0->Add(f3dHistPrimRawPtVsYCMSVsMultAntiLambda);
827 if(! f3dHistPrimRawPtVsYCMSVsMultK0Short) {
828 f3dHistPrimRawPtVsYCMSVsMultK0Short = new TH3F( "f3dHistPrimRawPtVsYCMSVsMultK0Short", "Pt_{K0S} Vs Y_{K0S} Vs Multiplicity; Pt_{K0S} (GeV/c); Y_{K0S} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
829 fListHistV0->Add(f3dHistPrimRawPtVsYCMSVsMultK0Short);
832 //---> Non-injected particles
834 if(! f3dHistPrimRawPtVsYVsMultNonInjLambda) {
835 f3dHistPrimRawPtVsYVsMultNonInjLambda = new TH3F( "f3dHistPrimRawPtVsYVsMultNonInjLambda", "Pt_{lambda} Vs Y_{#Lambda} Vs Multiplicity; Pt_{lambda} (GeV/c); Y_{#Lambda} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
836 fListHistV0->Add(f3dHistPrimRawPtVsYVsMultNonInjLambda);
838 if(! f3dHistPrimRawPtVsYVsMultNonInjAntiLambda) {
839 f3dHistPrimRawPtVsYVsMultNonInjAntiLambda = new TH3F( "f3dHistPrimRawPtVsYVsMultNonInjAntiLambda", "Pt_{antilambda} Vs Y_{#Lambda} Vs Multiplicity; Pt_{antilambda} (GeV/c); Y_{#Lambda} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
840 fListHistV0->Add(f3dHistPrimRawPtVsYVsMultNonInjAntiLambda);
842 if(! f3dHistPrimRawPtVsYVsMultNonInjK0Short) {
843 f3dHistPrimRawPtVsYVsMultNonInjK0Short = new TH3F( "f3dHistPrimRawPtVsYVsMultNonInjK0Short", "Pt_{K0S} Vs Y_{K0S} Vs Multiplicity; Pt_{K0S} (GeV/c); Y_{K0S} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
844 fListHistV0->Add(f3dHistPrimRawPtVsYVsMultNonInjK0Short);
847 //--- 3D Histo (Pt, Y, MultiplicityMC)
849 if(! f3dHistPrimRawPtVsYVsMultMCLambda) {
850 f3dHistPrimRawPtVsYVsMultMCLambda = new TH3F( "f3dHistPrimRawPtVsYVsMultMCLambda", "Pt_{lambda} Vs Y_{#Lambda} Vs Multiplicity; Pt_{lambda} (GeV/c); Y_{#Lambda} ; MultMC", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
851 fListHistV0->Add(f3dHistPrimRawPtVsYVsMultMCLambda);
853 if(! f3dHistPrimRawPtVsYVsMultMCAntiLambda) {
854 f3dHistPrimRawPtVsYVsMultMCAntiLambda = new TH3F( "f3dHistPrimRawPtVsYVsMultMCAntiLambda", "Pt_{antilambda} Vs Y_{#Lambda} Vs Multiplicity; Pt_{antilambda} (GeV/c); Y_{#Lambda} ; MultMC", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
855 fListHistV0->Add(f3dHistPrimRawPtVsYVsMultMCAntiLambda);
857 if(! f3dHistPrimRawPtVsYVsMultMCK0Short) {
858 f3dHistPrimRawPtVsYVsMultMCK0Short = new TH3F( "f3dHistPrimRawPtVsYVsMultMCK0Short", "Pt_{K0S} Vs Y_{K0S} Vs Multiplicity; Pt_{K0S} (GeV/c); Y_{K0S} ; MultMC", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
859 fListHistV0->Add(f3dHistPrimRawPtVsYVsMultMCK0Short);
862 //--- 3D Histo (Pt, Y, VertexZ)
864 if(! f3dHistPrimRawPtVsYVsVertexZLambda) {
865 f3dHistPrimRawPtVsYVsVertexZLambda = new TH3F( "f3dHistPrimRawPtVsYVsVertexZLambda", "Pt_{lambda} Vs Y_{#Lambda} Vs VertexZiplicity; Pt_{lambda} (GeV/c); Y_{#Lambda} ; VertexZ", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,40,-10,10);
866 fListHistV0->Add(f3dHistPrimRawPtVsYVsVertexZLambda);
868 if(! f3dHistPrimRawPtVsYVsVertexZAntiLambda) {
869 f3dHistPrimRawPtVsYVsVertexZAntiLambda = new TH3F( "f3dHistPrimRawPtVsYVsVertexZAntiLambda", "Pt_{antilambda} Vs Y_{#Lambda} Vs VertexZiplicity; Pt_{antilambda} (GeV/c); Y_{#Lambda} ; VertexZ", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,40,-10,10);
870 fListHistV0->Add(f3dHistPrimRawPtVsYVsVertexZAntiLambda);
872 if(! f3dHistPrimRawPtVsYVsVertexZK0Short) {
873 f3dHistPrimRawPtVsYVsVertexZK0Short = new TH3F( "f3dHistPrimRawPtVsYVsVertexZK0Short", "Pt_{K0S} Vs Y_{K0S} Vs VertexZiplicity; Pt_{K0S} (GeV/c); Y_{K0S} ; VertexZ", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,40,-10,10);
874 fListHistV0->Add(f3dHistPrimRawPtVsYVsVertexZK0Short);
877 //--- 3D Histo (Pt, Y, Multiplicity), close to PV criterion
879 if(! f3dHistPrimCloseToPVPtVsYVsMultLambda) {
880 f3dHistPrimCloseToPVPtVsYVsMultLambda = new TH3F( "f3dHistPrimCloseToPVPtVsYVsMultLambda", "Pt_{lambda} Vs Y_{#Lambda} Vs Multiplicity; Pt_{lambda} (GeV/c); Y_{#Lambda} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
881 fListHistV0->Add(f3dHistPrimCloseToPVPtVsYVsMultLambda);
883 if(! f3dHistPrimCloseToPVPtVsYVsMultAntiLambda) {
884 f3dHistPrimCloseToPVPtVsYVsMultAntiLambda = new TH3F( "f3dHistPrimCloseToPVPtVsYVsMultAntiLambda", "Pt_{antilambda} Vs Y_{#Lambda} Vs Multiplicity; Pt_{antilambda} (GeV/c); Y_{#Lambda} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
885 fListHistV0->Add(f3dHistPrimCloseToPVPtVsYVsMultAntiLambda);
887 if(! f3dHistPrimCloseToPVPtVsYVsMultK0Short) {
888 f3dHistPrimCloseToPVPtVsYVsMultK0Short = new TH3F( "f3dHistPrimCloseToPVPtVsYVsMultK0Short", "Pt_{K0S} Vs Y_{K0S} Vs Multiplicity; Pt_{K0S} (GeV/c); Y_{K0S} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
889 fListHistV0->Add(f3dHistPrimCloseToPVPtVsYVsMultK0Short);
893 //--- 3D Histo (Pt, Y, Proper Decay Length)
895 if(! f3dHistPrimRawPtVsYVsDecayLengthLambda) {
896 f3dHistPrimRawPtVsYVsDecayLengthLambda = new TH3F( "f3dHistPrimRawPtVsYVsDecayLengthLambda", "Pt_{lambda} Vs Y_{#Lambda} Vs DecayLength; Pt_{lambda} (GeV/c); Y_{#Lambda} ; DecayLength", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,200,0,50);
897 fListHistV0->Add(f3dHistPrimRawPtVsYVsDecayLengthLambda);
899 if(! f3dHistPrimRawPtVsYVsDecayLengthAntiLambda) {
900 f3dHistPrimRawPtVsYVsDecayLengthAntiLambda = new TH3F( "f3dHistPrimRawPtVsYVsDecayLengthAntiLambda", "Pt_{antilambda} Vs Y_{#Lambda} Vs DecayLength; Pt_{antilambda} (GeV/c); Y_{#Lambda} ; DecayLength", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,200,0,50);
901 fListHistV0->Add(f3dHistPrimRawPtVsYVsDecayLengthAntiLambda);
903 if(! f3dHistPrimRawPtVsYVsDecayLengthK0Short) {
904 f3dHistPrimRawPtVsYVsDecayLengthK0Short = new TH3F( "f3dHistPrimRawPtVsYVsDecayLengthK0Short", "Pt_{K0S} Vs Y_{K0S} Vs DecayLength; Pt_{K0S} (GeV/c); Y_{K0S} ; DecayLength", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,200,0,50);
905 fListHistV0->Add(f3dHistPrimRawPtVsYVsDecayLengthK0Short);
908 //--------------------------------------------------------------------------------------
909 //--- 3D Histo (Pt, Y, Multiplicity) for generated XiMinus/Plus, all generated
911 if(! f3dHistGenPtVsYVsMultXiMinus) {
912 f3dHistGenPtVsYVsMultXiMinus = new TH3F( "f3dHistGenPtVsYVsMultXiMinus", "Pt_{#Xi} Vs Y_{#Xi} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Xi} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
913 fListHistV0->Add(f3dHistGenPtVsYVsMultXiMinus);
915 if(! f3dHistGenPtVsYVsMultXiPlus) {
916 f3dHistGenPtVsYVsMultXiPlus = new TH3F( "f3dHistGenPtVsYVsMultXiPlus", "Pt_{#Xi} Vs Y_{#Xi} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Xi} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
917 fListHistV0->Add(f3dHistGenPtVsYVsMultXiPlus);
919 //--- 3D Histo (Pt, Y, Multiplicity) for generated OmegaMinus/Plus
921 if(! f3dHistGenPtVsYVsMultOmegaMinus) {
922 f3dHistGenPtVsYVsMultOmegaMinus = new TH3F( "f3dHistGenPtVsYVsMultOmegaMinus", "Pt_{#Omega} Vs Y_{#Omega} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Omega} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
923 fListHistV0->Add(f3dHistGenPtVsYVsMultOmegaMinus);
925 if(! f3dHistGenPtVsYVsMultOmegaPlus) {
926 f3dHistGenPtVsYVsMultOmegaPlus = new TH3F( "f3dHistGenPtVsYVsMultOmegaPlus", "Pt_{#Omega} Vs Y_{#Omega} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Omega} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
927 fListHistV0->Add(f3dHistGenPtVsYVsMultOmegaPlus);
932 if(! f3dHistGenPtVsYCMSVsMultXiMinus) {
933 f3dHistGenPtVsYCMSVsMultXiMinus = new TH3F( "f3dHistGenPtVsYCMSVsMultXiMinus", "Pt_{#Xi} Vs Y_{#Xi} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Xi} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
934 fListHistV0->Add(f3dHistGenPtVsYCMSVsMultXiMinus);
936 if(! f3dHistGenPtVsYCMSVsMultXiPlus) {
937 f3dHistGenPtVsYCMSVsMultXiPlus = new TH3F( "f3dHistGenPtVsYCMSVsMultXiPlus", "Pt_{#Xi} Vs Y_{#Xi} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Xi} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
938 fListHistV0->Add(f3dHistGenPtVsYCMSVsMultXiPlus);
940 //--- 3D Histo (Pt, Y, Multiplicity) for generated OmegaMinus/Plus
942 if(! f3dHistGenPtVsYCMSVsMultOmegaMinus) {
943 f3dHistGenPtVsYCMSVsMultOmegaMinus = new TH3F( "f3dHistGenPtVsYCMSVsMultOmegaMinus", "Pt_{#Omega} Vs Y_{#Omega} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Omega} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
944 fListHistV0->Add(f3dHistGenPtVsYCMSVsMultOmegaMinus);
946 if(! f3dHistGenPtVsYCMSVsMultOmegaPlus) {
947 f3dHistGenPtVsYCMSVsMultOmegaPlus = new TH3F( "f3dHistGenPtVsYCMSVsMultOmegaPlus", "Pt_{#Omega} Vs Y_{#Omega} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Omega} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
948 fListHistV0->Add(f3dHistGenPtVsYCMSVsMultOmegaPlus);
953 //--------------------------------------------------------------------------------------
954 //--- 3D Histo (Pt, Y, Multiplicity) for generated XiMinus/Plus, at selected analysis evts
956 if(! f3dHistGenSelectedPtVsYVsMultXiMinus) {
957 f3dHistGenSelectedPtVsYVsMultXiMinus = new TH3F( "f3dHistGenSelectedPtVsYVsMultXiMinus", "Pt_{#Xi} Vs Y_{#Xi} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Xi} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
958 fListHistV0->Add(f3dHistGenSelectedPtVsYVsMultXiMinus);
960 if(! f3dHistGenSelectedPtVsYVsMultXiPlus) {
961 f3dHistGenSelectedPtVsYVsMultXiPlus = new TH3F( "f3dHistGenSelectedPtVsYVsMultXiPlus", "Pt_{#Xi} Vs Y_{#Xi} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Xi} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
962 fListHistV0->Add(f3dHistGenSelectedPtVsYVsMultXiPlus);
964 //--- 3D Histo (Pt, Y, Multiplicity) for generated OmegaMinus/Plus
966 if(! f3dHistGenSelectedPtVsYVsMultOmegaMinus) {
967 f3dHistGenSelectedPtVsYVsMultOmegaMinus = new TH3F( "f3dHistGenSelectedPtVsYVsMultOmegaMinus", "Pt_{#Omega} Vs Y_{#Omega} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Omega} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
968 fListHistV0->Add(f3dHistGenSelectedPtVsYVsMultOmegaMinus);
970 if(! f3dHistGenSelectedPtVsYVsMultOmegaPlus) {
971 f3dHistGenSelectedPtVsYVsMultOmegaPlus = new TH3F( "f3dHistGenSelectedPtVsYVsMultOmegaPlus", "Pt_{#Omega} Vs Y_{#Omega} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Omega} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
972 fListHistV0->Add(f3dHistGenSelectedPtVsYVsMultOmegaPlus);
975 //CASCADES, analysis level, y CMS
977 //--------------------------------------------------------------------------------------
978 //--- 3D Histo (Pt, Y, Multiplicity) for generated XiMinus/Plus, at selected analysis evts
980 if(! f3dHistGenSelectedPtVsYCMSVsMultXiMinus) {
981 f3dHistGenSelectedPtVsYCMSVsMultXiMinus = new TH3F( "f3dHistGenSelectedPtVsYCMSVsMultXiMinus", "Pt_{#Xi} Vs Y_{#Xi} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Xi} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
982 fListHistV0->Add(f3dHistGenSelectedPtVsYCMSVsMultXiMinus);
984 if(! f3dHistGenSelectedPtVsYCMSVsMultXiPlus) {
985 f3dHistGenSelectedPtVsYCMSVsMultXiPlus = new TH3F( "f3dHistGenSelectedPtVsYCMSVsMultXiPlus", "Pt_{#Xi} Vs Y_{#Xi} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Xi} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
986 fListHistV0->Add(f3dHistGenSelectedPtVsYCMSVsMultXiPlus);
988 //--- 3D Histo (Pt, Y, Multiplicity) for generated OmegaMinus/Plus
990 if(! f3dHistGenSelectedPtVsYCMSVsMultOmegaMinus) {
991 f3dHistGenSelectedPtVsYCMSVsMultOmegaMinus = new TH3F( "f3dHistGenSelectedPtVsYCMSVsMultOmegaMinus", "Pt_{#Omega} Vs Y_{#Omega} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Omega} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
992 fListHistV0->Add(f3dHistGenSelectedPtVsYCMSVsMultOmegaMinus);
994 if(! f3dHistGenSelectedPtVsYCMSVsMultOmegaPlus) {
995 f3dHistGenSelectedPtVsYCMSVsMultOmegaPlus = new TH3F( "f3dHistGenSelectedPtVsYCMSVsMultOmegaPlus", "Pt_{#Omega} Vs Y_{#Omega} Vs Multiplicity; Pt_{cascade} (GeV/c); Y_{#Omega} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
996 fListHistV0->Add(f3dHistGenSelectedPtVsYCMSVsMultOmegaPlus);
1000 //----------------------------------
1001 // Histos at analysis level
1002 //----------------------------------
1004 if(! f3dHistPrimAnalysisPtVsYVsMultLambda) {
1005 f3dHistPrimAnalysisPtVsYVsMultLambda = new TH3F( "f3dHistPrimAnalysisPtVsYVsMultLambda", "Pt_{lambda} Vs Y_{#Lambda} Vs Multiplicity; Pt_{lambda} (GeV/c); Y_{#Lambda} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
1006 fListHistV0->Add(f3dHistPrimAnalysisPtVsYVsMultLambda);
1008 if(! f3dHistPrimAnalysisPtVsYVsMultAntiLambda) {
1009 f3dHistPrimAnalysisPtVsYVsMultAntiLambda = new TH3F( "f3dHistPrimAnalysisPtVsYVsMultAntiLambda", "Pt_{antilambda} Vs Y_{#Lambda} Vs Multiplicity; Pt_{antilambda} (GeV/c); Y_{#Lambda} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
1010 fListHistV0->Add(f3dHistPrimAnalysisPtVsYVsMultAntiLambda);
1012 if(! f3dHistPrimAnalysisPtVsYVsMultK0Short) {
1013 f3dHistPrimAnalysisPtVsYVsMultK0Short = new TH3F( "f3dHistPrimAnalysisPtVsYVsMultK0Short", "Pt_{K0S} Vs Y_{K0S} Vs Multiplicity; Pt_{K0S} (GeV/c); Y_{K0S} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
1014 fListHistV0->Add(f3dHistPrimAnalysisPtVsYVsMultK0Short);
1017 if(! f3dHistPrimAnalysisPtVsYCMSVsMultLambda) {
1018 f3dHistPrimAnalysisPtVsYCMSVsMultLambda = new TH3F( "f3dHistPrimAnalysisPtVsYCMSVsMultLambda", "Pt_{lambda} Vs Y_{#Lambda} Vs Multiplicity; Pt_{lambda} (GeV/c); Y_{#Lambda} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
1019 fListHistV0->Add(f3dHistPrimAnalysisPtVsYCMSVsMultLambda);
1021 if(! f3dHistPrimAnalysisPtVsYCMSVsMultAntiLambda) {
1022 f3dHistPrimAnalysisPtVsYCMSVsMultAntiLambda = new TH3F( "f3dHistPrimAnalysisPtVsYCMSVsMultAntiLambda", "Pt_{antilambda} Vs Y_{#Lambda} Vs Multiplicity; Pt_{antilambda} (GeV/c); Y_{#Lambda} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
1023 fListHistV0->Add(f3dHistPrimAnalysisPtVsYCMSVsMultAntiLambda);
1025 if(! f3dHistPrimAnalysisPtVsYCMSVsMultK0Short) {
1026 f3dHistPrimAnalysisPtVsYCMSVsMultK0Short = new TH3F( "f3dHistPrimAnalysisPtVsYCMSVsMultK0Short", "Pt_{K0S} Vs Y_{K0S} Vs Multiplicity; Pt_{K0S} (GeV/c); Y_{K0S} ; Mult", lCustomNBins, 0., lCustomPtUpperLimit, 48, -1.2,1.2,lCustomNBinsMultiplicity,0,lCustomNBinsMultiplicity);
1027 fListHistV0->Add(f3dHistPrimAnalysisPtVsYCMSVsMultK0Short);
1030 //----------------------------------
1031 // Primary Vertex Position Histos
1032 //----------------------------------
1035 fHistPVx = new TH1F("fHistPVx",
1036 "PV x position;Nbr of Evts;x",
1038 fListHistV0->Add(fHistPVx);
1041 fHistPVy = new TH1F("fHistPVy",
1042 "PV y position;Nbr of Evts;y",
1044 fListHistV0->Add(fHistPVy);
1047 fHistPVz = new TH1F("fHistPVz",
1048 "PV z position;Nbr of Evts;z",
1050 fListHistV0->Add(fHistPVz);
1053 if(! fHistPVxAnalysis) {
1054 fHistPVxAnalysis = new TH1F("fHistPVxAnalysis",
1055 "PV x position;Nbr of Evts;x",
1057 fListHistV0->Add(fHistPVxAnalysis);
1059 if(! fHistPVyAnalysis) {
1060 fHistPVyAnalysis = new TH1F("fHistPVyAnalysis",
1061 "PV y position;Nbr of Evts;y",
1063 fListHistV0->Add(fHistPVyAnalysis);
1065 if(! fHistPVzAnalysis) {
1066 fHistPVzAnalysis = new TH1F("fHistPVzAnalysis",
1067 "PV z position;Nbr of Evts;z",
1069 fListHistV0->Add(fHistPVzAnalysis);
1072 if(! fHistPVxAnalysisHasHighPtLambda) {
1073 fHistPVxAnalysisHasHighPtLambda = new TH1F("fHistPVxAnalysisHasHighPtLambda",
1074 "PV x position;Nbr of Evts;x",
1076 fListHistV0->Add(fHistPVxAnalysisHasHighPtLambda);
1078 if(! fHistPVyAnalysisHasHighPtLambda) {
1079 fHistPVyAnalysisHasHighPtLambda = new TH1F("fHistPVyAnalysisHasHighPtLambda",
1080 "PV y position;Nbr of Evts;y",
1082 fListHistV0->Add(fHistPVyAnalysisHasHighPtLambda);
1084 if(! fHistPVzAnalysisHasHighPtLambda) {
1085 fHistPVzAnalysisHasHighPtLambda = new TH1F("fHistPVzAnalysisHasHighPtLambda",
1086 "PV z position;Nbr of Evts;z",
1088 fListHistV0->Add(fHistPVzAnalysisHasHighPtLambda);
1090 if(! fHistSwappedV0Counter) {
1091 fHistSwappedV0Counter = new TH1F("fHistSwappedV0Counter",
1092 "Swap or not histo;Swapped (1) or not (0); count",
1094 fListHistV0->Add(fHistSwappedV0Counter);
1097 //List of Histograms: Normal
1098 PostData(1, fListHistV0);
1100 //TTree Object: Saved to base directory. Should cache to disk while saving.
1101 //(Important to avoid excessive memory usage, particularly when merging)
1104 }// end UserCreateOutputObjects
1107 //________________________________________________________________________
1108 void AliAnalysisTaskExtractPerformanceV0::UserExec(Option_t *)
1111 // Called for each event
1113 AliESDEvent *lESDevent = 0x0;
1114 AliMCEvent *lMCevent = 0x0;
1115 AliStack *lMCstack = 0x0;
1117 Int_t lNumberOfV0s = -1;
1118 Double_t lTrkgPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
1119 Double_t lBestPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
1120 Double_t lMagneticField = -10.;
1122 // Connect to the InputEvent
1123 // After these lines, we should have an ESD/AOD event + the number of V0s in it.
1125 // Appropriate for ESD analysis!
1127 lESDevent = dynamic_cast<AliESDEvent*>( InputEvent() );
1129 AliWarning("ERROR: lESDevent not available \n");
1133 fTreeVariableRunNumber = lESDevent->GetRunNumber();
1134 fTreeVariableEventNumber =
1135 ( ( ((ULong64_t)lESDevent->GetPeriodNumber() ) << 36 ) |
1136 ( ((ULong64_t)lESDevent->GetOrbitNumber () ) << 12 ) |
1137 ((ULong64_t)lESDevent->GetBunchCrossNumber() ) );
1139 lMCevent = MCEvent();
1141 Printf("ERROR: Could not retrieve MC event \n");
1142 cout << "Name of the file with pb :" << fInputHandler->GetTree()->GetCurrentFile()->GetName() << endl;
1146 lMCstack = lMCevent->Stack();
1148 Printf("ERROR: Could not retrieve MC stack \n");
1149 cout << "Name of the file with pb :" << fInputHandler->GetTree()->GetCurrentFile()->GetName() << endl;
1152 TArrayF mcPrimaryVtx;
1153 AliGenEventHeader* mcHeader=lMCevent->GenEventHeader();
1154 if(!mcHeader) return;
1155 mcHeader->PrimaryVertex(mcPrimaryVtx);
1157 //------------------------------------------------
1158 // Multiplicity Information Acquistion
1159 //------------------------------------------------
1161 //REVISED multiplicity estimator after 'multiplicity day' (2011)
1162 Int_t lMultiplicity = -100;
1165 if(fkIsNuclear == kFALSE) lMultiplicity = fESDtrackCuts->GetReferenceMultiplicity(lESDevent, AliESDtrackCuts::kTrackletsITSTPC,0.5);
1167 //---> If this is a nuclear collision, then go nuclear on "multiplicity" variable...
1168 //---> Warning: Experimental
1169 if(fkIsNuclear == kTRUE){
1170 AliCentrality* centrality;
1171 centrality = lESDevent->GetCentrality();
1172 lMultiplicity = ( ( Int_t ) ( centrality->GetCentralityPercentile( fCentralityEstimator.Data() ) ) );
1173 if (centrality->GetQuality()>1) {
1174 PostData(1, fListHistV0);
1180 //Set variable for filling tree afterwards!
1181 //---> pp case......: GetReferenceMultiplicity
1182 //---> Pb-Pb case...: Centrality by V0M
1183 fTreeVariableMultiplicity = lMultiplicity;
1185 fHistV0MultiplicityBeforeTrigSel->Fill ( lESDevent->GetNumberOfV0s() );
1186 fHistMultiplicityBeforeTrigSel->Fill ( lMultiplicity );
1188 //------------------------------------------------
1189 // MC Information Acquistion
1190 //------------------------------------------------
1192 Int_t iNumberOfPrimaries = -1;
1193 iNumberOfPrimaries = lMCstack->GetNprimary();
1194 if(iNumberOfPrimaries < 1) return;
1195 Bool_t lHasHighPtLambda = kFALSE;
1197 //------------------------------------------------
1198 // Variable Definition
1199 //------------------------------------------------
1201 Int_t lNbMCPrimary = 0;
1203 Int_t lPdgcodeCurrentPart = 0;
1204 Double_t lRapCurrentPart = 0;
1205 Double_t lPtCurrentPart = 0;
1207 //Int_t lComeFromSigma = 0;
1209 // current mc particle 's mother
1210 //Int_t iCurrentMother = 0;
1211 lNbMCPrimary = lMCstack->GetNprimary();
1213 //------------------------------------------------
1214 // Pre-Physics Selection
1215 //------------------------------------------------
1217 //----- Loop on primary Xi, Omega --------------------------------------------------------------
1218 for (Int_t iCurrentLabelStack = 0; iCurrentLabelStack < lNbMCPrimary; iCurrentLabelStack++)
1219 {// This is the begining of the loop on primaries
1221 TParticle* lCurrentParticlePrimary = 0x0;
1222 lCurrentParticlePrimary = lMCstack->Particle( iCurrentLabelStack );
1223 if(!lCurrentParticlePrimary){
1224 Printf("Cascade loop %d - MC TParticle pointer to current stack particle = 0x0 ! Skip ...\n", iCurrentLabelStack );
1227 if ( TMath::Abs(lCurrentParticlePrimary->GetPdgCode()) == 3312 || TMath::Abs(lCurrentParticlePrimary->GetPdgCode()) == 3334 ) {
1228 Double_t lRapXiMCPrimary = -100;
1229 if( (lCurrentParticlePrimary->Energy() - lCurrentParticlePrimary->Pz() +1.e-13) != 0 ) {
1230 if ( (lCurrentParticlePrimary->Energy() + lCurrentParticlePrimary->Pz()) / (lCurrentParticlePrimary->Energy() - lCurrentParticlePrimary->Pz() +1.e-13) !=0 ){
1231 lRapXiMCPrimary = 0.5*TMath::Log( (lCurrentParticlePrimary->Energy() + lCurrentParticlePrimary->Pz()) / (lCurrentParticlePrimary->Energy() - lCurrentParticlePrimary->Pz() +1.e-13) );
1235 //=================================================================================
1237 if( lCurrentParticlePrimary->GetPdgCode() == 3312 ){
1238 lPtCurrentPart = lCurrentParticlePrimary->Pt();
1239 f3dHistGenPtVsYVsMultXiMinus->Fill(lPtCurrentPart, lRapXiMCPrimary, lMultiplicity);
1240 f3dHistGenPtVsYCMSVsMultXiMinus->Fill(lPtCurrentPart, lRapXiMCPrimary+fpArapidityShift, lMultiplicity);
1242 if( lCurrentParticlePrimary->GetPdgCode() == -3312 ){
1243 lPtCurrentPart = lCurrentParticlePrimary->Pt();
1244 f3dHistGenPtVsYVsMultXiPlus->Fill(lPtCurrentPart, lRapXiMCPrimary, lMultiplicity);
1245 f3dHistGenPtVsYCMSVsMultXiPlus->Fill(lPtCurrentPart, lRapXiMCPrimary+fpArapidityShift, lMultiplicity);
1248 if( lCurrentParticlePrimary->GetPdgCode() == 3334 ){
1249 lPtCurrentPart = lCurrentParticlePrimary->Pt();
1250 f3dHistGenPtVsYVsMultOmegaMinus->Fill(lPtCurrentPart, lRapXiMCPrimary, lMultiplicity);
1251 f3dHistGenPtVsYCMSVsMultOmegaMinus->Fill(lPtCurrentPart, lRapXiMCPrimary+fpArapidityShift, lMultiplicity);
1253 if( lCurrentParticlePrimary->GetPdgCode() == -3334 ){
1254 lPtCurrentPart = lCurrentParticlePrimary->Pt();
1255 f3dHistGenPtVsYVsMultOmegaPlus->Fill(lPtCurrentPart, lRapXiMCPrimary, lMultiplicity);
1256 f3dHistGenPtVsYCMSVsMultOmegaPlus->Fill(lPtCurrentPart, lRapXiMCPrimary+fpArapidityShift, lMultiplicity);
1260 //----- End Loop on primary Xi, Omega ----------------------------------------------------------
1262 //--------- GENERATED NUMBER OF CHARGED PARTICLES
1263 // ---> Set Variables to Zero again
1264 // ---> Variable Definition
1266 Long_t lNumberOfCharged = 0;
1268 //----- Loop on Stack ----------------------------------------------------------------
1269 for (Int_t iCurrentLabelStack = 0; iCurrentLabelStack < (lMCstack->GetNtrack()); iCurrentLabelStack++)
1270 {// This is the begining of the loop on tracks
1271 TParticle* particleOne = lMCstack->Particle(iCurrentLabelStack);
1272 if(!particleOne) continue;
1273 if(!particleOne->GetPDG()) continue;
1274 Double_t lThisCharge = particleOne->GetPDG()->Charge()/3.;
1275 if(TMath::Abs(lThisCharge)<0.001) continue;
1276 if(! (lMCstack->IsPhysicalPrimary(iCurrentLabelStack)) ) continue;
1278 Double_t gpt = particleOne -> Pt();
1279 Double_t geta = particleOne -> Eta();
1281 if( TMath::Abs(geta) < 0.5) lNumberOfCharged++;
1282 }//End of loop on tracks
1283 //----- End Loop on Stack ------------------------------------------------------------
1285 //Double_t lpArapidityShift = 0.465;
1286 Bool_t lStackNatural = kTRUE;
1287 //----- Loop on Lambda, K0Short ----------------------------------------------------------------
1288 for (Int_t iCurrentLabelStack = 0; iCurrentLabelStack < (lMCstack->GetNtrack()); iCurrentLabelStack++)
1289 {// This is the begining of the loop on tracks
1291 TParticle* lCurrentParticleForLambdaCheck = 0x0;
1292 lCurrentParticleForLambdaCheck = lMCstack->Particle( iCurrentLabelStack );
1293 if(!lCurrentParticleForLambdaCheck){
1294 Printf("V0s loop %d - MC TParticle pointer to current stack particle = 0x0 ! Skip ...\n", iCurrentLabelStack );
1298 //=================================================================================
1299 //Single-Strange checks
1300 // Keep only K0s, Lambda and AntiLambda:
1301 lPdgcodeCurrentPart = lCurrentParticleForLambdaCheck->GetPdgCode();
1303 if ( (lCurrentParticleForLambdaCheck->GetPdgCode() == 310 ) ||
1304 (lCurrentParticleForLambdaCheck->GetPdgCode() == 3122 ) ||
1305 (lCurrentParticleForLambdaCheck->GetPdgCode() == -3122 ) )
1307 lRapCurrentPart = MyRapidity(lCurrentParticleForLambdaCheck->Energy(),lCurrentParticleForLambdaCheck->Pz());
1308 lPtCurrentPart = lCurrentParticleForLambdaCheck->Pt();
1310 //Use Close to PV for filling CloseToPV histograms!
1311 Double_t dx, dy, dz;
1313 dx = ( (mcPrimaryVtx.At(0)) - (lCurrentParticleForLambdaCheck->Vx()) );
1314 dy = ( (mcPrimaryVtx.At(1)) - (lCurrentParticleForLambdaCheck->Vy()) );
1315 dz = ( (mcPrimaryVtx.At(2)) - (lCurrentParticleForLambdaCheck->Vz()) );
1316 Double_t lDistToPV = TMath::Sqrt(dx*dx + dy*dy + dz*dz);
1317 if( lDistToPV <= 0.001){
1318 if( lPdgcodeCurrentPart == 3122 ){
1319 f3dHistPrimCloseToPVPtVsYVsMultLambda->Fill(lPtCurrentPart, lRapCurrentPart, lMultiplicity);
1321 if( lPdgcodeCurrentPart == -3122 ){
1322 f3dHistPrimCloseToPVPtVsYVsMultAntiLambda->Fill(lPtCurrentPart, lRapCurrentPart, lMultiplicity);
1324 if( lPdgcodeCurrentPart == 310 ){
1325 f3dHistPrimCloseToPVPtVsYVsMultK0Short->Fill(lPtCurrentPart, lRapCurrentPart, lMultiplicity);
1329 //Use Physical Primaries only for filling PrimRaw Histograms!
1330 if ( lMCstack->IsPhysicalPrimary(iCurrentLabelStack)!=kTRUE ) continue;
1332 lStackNatural = lMCevent->IsFromBGEvent(iCurrentLabelStack); //Is it?
1333 if (!lStackNatural){
1334 if (!(lCurrentParticleForLambdaCheck->GetFirstMother()<0))
1335 {lStackNatural = kTRUE;} // because there are primaries (ALICE definition) not produced in the collision
1338 if( lPdgcodeCurrentPart == 3122 ){
1339 f3dHistPrimRawPtVsYVsMultLambda->Fill(lPtCurrentPart, lRapCurrentPart, lMultiplicity);
1340 f3dHistPrimRawPtVsYCMSVsMultLambda->Fill(lPtCurrentPart, lRapCurrentPart+fpArapidityShift, lMultiplicity);
1341 if(lStackNatural){f3dHistPrimRawPtVsYVsMultNonInjLambda->Fill(lPtCurrentPart, lRapCurrentPart, lMultiplicity);}
1342 f3dHistPrimRawPtVsYVsMultMCLambda->Fill(lPtCurrentPart, lRapCurrentPart, lNumberOfCharged);
1343 f3dHistPrimRawPtVsYVsVertexZLambda->Fill(lPtCurrentPart, lRapCurrentPart, mcPrimaryVtx.At(2));
1344 if( TMath::Abs( lCurrentParticleForLambdaCheck->Eta() )<1.2 && lPtCurrentPart>2 ){
1345 lHasHighPtLambda = kTRUE; //Keep track of events with Lambda within |eta|<1.2 and pt>2
1348 if( lPdgcodeCurrentPart == -3122 ){
1349 f3dHistPrimRawPtVsYVsMultAntiLambda->Fill(lPtCurrentPart, lRapCurrentPart, lMultiplicity);
1350 f3dHistPrimRawPtVsYCMSVsMultAntiLambda->Fill(lPtCurrentPart, lRapCurrentPart+fpArapidityShift, lMultiplicity);
1351 if(lStackNatural){f3dHistPrimRawPtVsYVsMultNonInjAntiLambda->Fill(lPtCurrentPart, lRapCurrentPart, lMultiplicity);}
1352 f3dHistPrimRawPtVsYVsMultMCAntiLambda->Fill(lPtCurrentPart, lRapCurrentPart, lNumberOfCharged);
1353 f3dHistPrimRawPtVsYVsVertexZAntiLambda->Fill(lPtCurrentPart, lRapCurrentPart, mcPrimaryVtx.At(2));
1355 if( lPdgcodeCurrentPart == 310 ){
1356 f3dHistPrimRawPtVsYVsMultK0Short->Fill(lPtCurrentPart, lRapCurrentPart, lMultiplicity);
1357 f3dHistPrimRawPtVsYCMSVsMultK0Short->Fill(lPtCurrentPart, lRapCurrentPart+fpArapidityShift, lMultiplicity);
1358 if(lStackNatural){f3dHistPrimRawPtVsYVsMultNonInjK0Short->Fill(lPtCurrentPart, lRapCurrentPart, lMultiplicity);}
1359 f3dHistPrimRawPtVsYVsMultMCK0Short->Fill(lPtCurrentPart, lRapCurrentPart, lNumberOfCharged);
1360 f3dHistPrimRawPtVsYVsVertexZK0Short->Fill(lPtCurrentPart, lRapCurrentPart, mcPrimaryVtx.At(2));
1362 //Decay Length Acquisition=====================================================
1363 Double_t decaylength = -1;
1364 Double_t lV0Mass = -1;
1366 if( !(lCurrentParticleForLambdaCheck->GetDaughter(0) < 0) ) {
1367 TParticle* lDght0ofV0 = lMCstack->Particle( lCurrentParticleForLambdaCheck->GetDaughter(0) ); //get first daughter
1368 if(lDght0ofV0){ // skip if not defined.
1369 decaylength = TMath::Sqrt(
1370 TMath::Power( lCurrentParticleForLambdaCheck->Vx() - lDght0ofV0->Vx() , 2) +
1371 TMath::Power( lCurrentParticleForLambdaCheck->Vy() - lDght0ofV0->Vy() , 2) +
1372 TMath::Power( lCurrentParticleForLambdaCheck->Vz() - lDght0ofV0->Vz() , 2)
1374 //Need to correct for relativitity! Involves multiplying by mass and dividing by momentum.
1375 if(TMath::Abs( lPdgcodeCurrentPart ) == 3122 ) { lV0Mass = 1.115683; }
1376 if(TMath::Abs( lPdgcodeCurrentPart ) == 310 ) { lV0Mass = 0.497614; }
1377 if( lCurrentParticleForLambdaCheck->P() + 1e-10 != 0 ) decaylength = ( lV0Mass * decaylength ) / ( lCurrentParticleForLambdaCheck->P() + 1e-10 );
1378 if( lCurrentParticleForLambdaCheck->P() + 1e-10 == 0 ) decaylength = 1e+5;
1381 if( lPdgcodeCurrentPart == 3122) f3dHistPrimRawPtVsYVsDecayLengthLambda ->Fill( lPtCurrentPart, lRapCurrentPart , decaylength );
1382 if( lPdgcodeCurrentPart == -3122) f3dHistPrimRawPtVsYVsDecayLengthAntiLambda ->Fill( lPtCurrentPart, lRapCurrentPart , decaylength );
1383 if( lPdgcodeCurrentPart == 310) f3dHistPrimRawPtVsYVsDecayLengthK0Short ->Fill( lPtCurrentPart, lRapCurrentPart , decaylength );
1385 }//End of loop on tracks
1386 //----- End Loop on Lambda, K0Short ------------------------------------------------------------
1389 f2dHistMultiplicityVsTrueBeforeTrigSel->Fill ( lMultiplicity , lNumberOfCharged );
1391 fTreeVariableMultiplicityMC = lNumberOfCharged;
1393 fHistGenVertexZBeforeTrigSel->Fill( (mcPrimaryVtx.At(2)) );
1395 lPdgcodeCurrentPart = 0;
1396 lRapCurrentPart = 0;
1399 //------------------------------------------------
1400 // Physics Selection
1401 //------------------------------------------------
1403 UInt_t maskIsSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
1404 Bool_t isSelected = 0;
1405 isSelected = (maskIsSelected & AliVEvent::kMB) == AliVEvent::kMB;
1407 //pA triggering: CINT7
1408 if( fkSwitchINT7 ) isSelected = (maskIsSelected & AliVEvent::kINT7) == AliVEvent::kINT7;
1410 //Standard Min-Bias Selection
1411 if ( ! isSelected ) {
1412 PostData(1, fListHistV0);
1417 f2dHistMultiplicityVsTrueForTrigEvt->Fill ( lMultiplicity , lNumberOfCharged );
1418 fHistGenVertexZForTrigEvt->Fill( mcPrimaryVtx.At(2) );
1419 //------------------------------------------------
1420 // After Trigger Selection
1421 //------------------------------------------------
1423 lNumberOfV0s = lESDevent->GetNumberOfV0s();
1425 //Set variable for filling tree afterwards!
1426 fHistV0MultiplicityForTrigEvt->Fill(lNumberOfV0s);
1427 fHistMultiplicityForTrigEvt->Fill ( lMultiplicity );
1429 //------------------------------------------------
1430 // Getting: Primary Vertex + MagField Info
1431 //------------------------------------------------
1433 const AliESDVertex *lPrimaryTrackingESDVtx = lESDevent->GetPrimaryVertexTracks();
1434 // get the vtx stored in ESD found with tracks
1435 lPrimaryTrackingESDVtx->GetXYZ( lTrkgPrimaryVtxPos );
1437 const AliESDVertex *lPrimaryBestESDVtx = lESDevent->GetPrimaryVertex();
1438 // get the best primary vertex available for the event
1439 // As done in AliCascadeVertexer, we keep the one which is the best one available.
1440 // between : Tracking vertex > SPD vertex > TPC vertex > default SPD vertex
1441 // This one will be used for next calculations (DCA essentially)
1442 lPrimaryBestESDVtx->GetXYZ( lBestPrimaryVtxPos );
1444 Double_t lPrimaryVtxPosition[3];
1445 const AliVVertex *primaryVtx = lESDevent->GetPrimaryVertex();
1446 lPrimaryVtxPosition[0] = primaryVtx->GetX();
1447 lPrimaryVtxPosition[1] = primaryVtx->GetY();
1448 lPrimaryVtxPosition[2] = primaryVtx->GetZ();
1449 fHistPVx->Fill( lPrimaryVtxPosition[0] );
1450 fHistPVy->Fill( lPrimaryVtxPosition[1] );
1451 fHistPVz->Fill( lPrimaryVtxPosition[2] );
1453 f2dHistMultiplicityVsVertexZForTrigEvt->Fill( lMultiplicity, lPrimaryVtxPosition[2] );
1455 //------------------------------------------------
1456 // Primary Vertex Z position: SKIP
1457 //------------------------------------------------
1459 if(TMath::Abs(lBestPrimaryVtxPos[2]) > 10.0 ) {
1460 AliWarning("Pb / | Z position of Best Prim Vtx | > 10.0 cm ... return !");
1461 PostData(1, fListHistV0);
1466 f2dHistMultiplicityVsVertexZ->Fill( lMultiplicity, lPrimaryVtxPosition[2] );
1468 lMagneticField = lESDevent->GetMagneticField( );
1469 fHistV0MultiplicityForSelEvt ->Fill( lNumberOfV0s );
1470 fHistMultiplicity->Fill(lMultiplicity);
1471 f2dHistMultiplicityVsTrue->Fill ( lMultiplicity , lNumberOfCharged );
1472 fHistGenVertexZ->Fill( (mcPrimaryVtx.At(2)) );
1473 //------------------------------------------------
1474 // SKIP: Events with well-established PVtx
1475 //------------------------------------------------
1477 const AliESDVertex *lPrimaryTrackingESDVtxCheck = lESDevent->GetPrimaryVertexTracks();
1478 const AliESDVertex *lPrimarySPDVtx = lESDevent->GetPrimaryVertexSPD();
1479 if (!lPrimarySPDVtx->GetStatus() && !lPrimaryTrackingESDVtxCheck->GetStatus() ){
1480 AliWarning("Pb / No SPD prim. vertex nor prim. Tracking vertex ... return !");
1481 PostData(1, fListHistV0);
1486 f2dHistMultiplicityVsVertexZNoTPCOnly->Fill( lMultiplicity, lPrimaryVtxPosition[2] );
1487 fHistV0MultiplicityForSelEvtNoTPCOnly ->Fill( lNumberOfV0s );
1488 fHistMultiplicityNoTPCOnly->Fill(lMultiplicity);
1489 f2dHistMultiplicityVsTrueNoTPCOnly->Fill ( lMultiplicity , lNumberOfCharged );
1490 fHistGenVertexZNoTPCOnly->Fill( (mcPrimaryVtx.At(2)) );
1491 //------------------------------------------------
1492 // Pileup Rejection Studies
1493 //------------------------------------------------
1495 // FIXME : quality selection regarding pile-up rejection
1496 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
1497 AliWarning("Pb / This is tagged as Pileup from SPD... return !");
1498 PostData(1, fListHistV0);
1502 f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup->Fill( lMultiplicity, lPrimaryVtxPosition[2] );
1503 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup ->Fill( lNumberOfV0s );
1504 fHistMultiplicityNoTPCOnlyNoPileup->Fill(lMultiplicity);
1505 f2dHistMultiplicityVsTrueNoTPCOnlyNoPileup->Fill ( lMultiplicity , lNumberOfCharged );
1506 fHistGenVertexZNoTPCOnlyNoPileup->Fill( (mcPrimaryVtx.At(2)) );
1507 //Do control histograms without the IsFromVertexerZ events, but consider them in analysis...
1508 if( ! (lESDevent->GetPrimaryVertex()->IsFromVertexerZ() ) ){
1509 fHistPVxAnalysis->Fill( lPrimaryVtxPosition[0] );
1510 fHistPVyAnalysis->Fill( lPrimaryVtxPosition[1] );
1511 fHistPVzAnalysis->Fill( lPrimaryVtxPosition[2] );
1512 if ( lHasHighPtLambda == kTRUE ){
1513 fHistPVxAnalysisHasHighPtLambda->Fill( lPrimaryVtxPosition[0] );
1514 fHistPVyAnalysisHasHighPtLambda->Fill( lPrimaryVtxPosition[1] );
1515 fHistPVzAnalysisHasHighPtLambda->Fill( lPrimaryVtxPosition[2] );
1519 fTreeVariableVertexZ = lPrimaryVtxPosition[2];
1521 fTreeVariablePVx = lPrimaryVtxPosition[0];
1522 fTreeVariablePVy = lPrimaryVtxPosition[1];
1523 fTreeVariablePVz = lPrimaryVtxPosition[2];
1525 fTreeVariableMCPVx = (mcPrimaryVtx.At(0));
1526 fTreeVariableMCPVy = (mcPrimaryVtx.At(1));
1527 fTreeVariableMCPVz = (mcPrimaryVtx.At(2));
1529 //------------------------------------------------
1530 // stack loop starts here
1531 //------------------------------------------------
1533 //---> Loop over ALL PARTICLES
1535 for (Int_t iMc = 0; iMc < (lMCstack->GetNtrack()); iMc++) {
1536 TParticle *p0 = lMCstack->Particle(iMc);
1538 //Printf("ERROR: particle with label %d not found in lMCstack (mc loop)", iMc);
1541 lPdgcodeCurrentPart = p0->GetPdgCode();
1543 // Keep only K0s, Lambda and AntiLambda:
1544 if ( (lPdgcodeCurrentPart != 310 ) && (lPdgcodeCurrentPart != 3122 ) && (lPdgcodeCurrentPart != -3122 ) ) continue;
1546 lRapCurrentPart = MyRapidity(p0->Energy(),p0->Pz());
1547 lPtCurrentPart = p0->Pt();
1549 //Use Physical Primaries only for filling PrimRaw Histograms!
1550 if ( lMCstack->IsPhysicalPrimary(iMc)!=kTRUE ) continue;
1552 if( lPdgcodeCurrentPart == 3122 ){
1553 f3dHistPrimAnalysisPtVsYVsMultLambda->Fill(lPtCurrentPart, lRapCurrentPart, lMultiplicity);
1554 f3dHistPrimAnalysisPtVsYCMSVsMultLambda->Fill(lPtCurrentPart, lRapCurrentPart+fpArapidityShift, lMultiplicity);
1556 if( lPdgcodeCurrentPart == -3122 ){
1557 f3dHistPrimAnalysisPtVsYVsMultAntiLambda->Fill(lPtCurrentPart, lRapCurrentPart, lMultiplicity);
1558 f3dHistPrimAnalysisPtVsYCMSVsMultAntiLambda->Fill(lPtCurrentPart, lRapCurrentPart+fpArapidityShift, lMultiplicity);
1560 if( lPdgcodeCurrentPart == 310 ){
1561 f3dHistPrimAnalysisPtVsYVsMultK0Short->Fill(lPtCurrentPart, lRapCurrentPart, lMultiplicity);
1562 f3dHistPrimAnalysisPtVsYCMSVsMultK0Short->Fill(lPtCurrentPart, lRapCurrentPart+fpArapidityShift, lMultiplicity);
1566 //----- Loop on primary Xi, Omega --------------------------------------------------------------
1567 for (Int_t iCurrentLabelStack = 0; iCurrentLabelStack < lNbMCPrimary; iCurrentLabelStack++)
1568 {// This is the begining of the loop on primaries
1570 TParticle* lCurrentParticlePrimary = 0x0;
1571 lCurrentParticlePrimary = lMCstack->Particle( iCurrentLabelStack );
1572 if(!lCurrentParticlePrimary){
1573 Printf("Cascade loop %d - MC TParticle pointer to current stack particle = 0x0 ! Skip ...\n", iCurrentLabelStack );
1576 if ( TMath::Abs(lCurrentParticlePrimary->GetPdgCode()) == 3312 || TMath::Abs(lCurrentParticlePrimary->GetPdgCode()) == 3334 ) {
1577 Double_t lRapXiMCPrimary = -100;
1578 if( (lCurrentParticlePrimary->Energy() - lCurrentParticlePrimary->Pz() +1.e-13) != 0 ) {
1579 if ( (lCurrentParticlePrimary->Energy() + lCurrentParticlePrimary->Pz()) / (lCurrentParticlePrimary->Energy() - lCurrentParticlePrimary->Pz() +1.e-13) !=0 ){
1580 lRapXiMCPrimary = 0.5*TMath::Log( (lCurrentParticlePrimary->Energy() + lCurrentParticlePrimary->Pz()) / (lCurrentParticlePrimary->Energy() - lCurrentParticlePrimary->Pz() +1.e-13) );
1584 //=================================================================================
1586 if( lCurrentParticlePrimary->GetPdgCode() == 3312 ){
1587 lPtCurrentPart = lCurrentParticlePrimary->Pt();
1588 f3dHistGenSelectedPtVsYVsMultXiMinus->Fill(lPtCurrentPart, lRapXiMCPrimary, lMultiplicity);
1589 f3dHistGenSelectedPtVsYCMSVsMultXiMinus->Fill(lPtCurrentPart, lRapXiMCPrimary+fpArapidityShift, lMultiplicity);
1591 if( lCurrentParticlePrimary->GetPdgCode() == -3312 ){
1592 lPtCurrentPart = lCurrentParticlePrimary->Pt();
1593 f3dHistGenSelectedPtVsYVsMultXiPlus->Fill(lPtCurrentPart, lRapXiMCPrimary, lMultiplicity);
1594 f3dHistGenSelectedPtVsYCMSVsMultXiPlus->Fill(lPtCurrentPart, lRapXiMCPrimary+fpArapidityShift, lMultiplicity);
1597 if( lCurrentParticlePrimary->GetPdgCode() == 3334 ){
1598 lPtCurrentPart = lCurrentParticlePrimary->Pt();
1599 f3dHistGenSelectedPtVsYVsMultOmegaMinus->Fill(lPtCurrentPart, lRapXiMCPrimary, lMultiplicity);
1600 f3dHistGenSelectedPtVsYCMSVsMultOmegaMinus->Fill(lPtCurrentPart, lRapXiMCPrimary+fpArapidityShift, lMultiplicity);
1602 if( lCurrentParticlePrimary->GetPdgCode() == -3334 ){
1603 lPtCurrentPart = lCurrentParticlePrimary->Pt();
1604 f3dHistGenSelectedPtVsYVsMultOmegaPlus->Fill(lPtCurrentPart, lRapXiMCPrimary, lMultiplicity);
1605 f3dHistGenSelectedPtVsYCMSVsMultOmegaPlus->Fill(lPtCurrentPart, lRapXiMCPrimary+fpArapidityShift, lMultiplicity);
1609 //----- End Loop on primary Xi, Omega ----------------------------------------------------------
1611 //------------------------------------------------
1612 // MAIN LAMBDA LOOP STARTS HERE
1613 //------------------------------------------------
1615 //Variable definition
1616 Int_t lOnFlyStatus = 0;
1617 Double_t lChi2V0 = 0;
1618 Double_t lDcaV0Daughters = 0, lDcaV0ToPrimVertex = 0;
1619 Double_t lDcaPosToPrimVertex = 0, lDcaNegToPrimVertex = 0;
1620 Double_t lV0CosineOfPointingAngle = 0;
1621 Double_t lV0Radius = 0, lPt = 0;
1622 Double_t lRapK0Short = 0, lRapLambda = 0;
1623 Double_t lInvMassK0s = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
1624 Double_t lAlphaV0 = 0, lPtArmV0 = 0;
1625 Double_t fMinV0Pt = 0;
1626 Double_t fMaxV0Pt = 100;
1629 nv0s = lESDevent->GetNumberOfV0s();
1631 for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
1632 {// This is the begining of the V0 loop
1633 AliESDv0 *v0 = ((AliESDEvent*)lESDevent)->GetV0(iV0);
1636 //---> Fix On-the-Fly candidates
1637 if( v0->GetParamN()->Charge() > 0 && v0->GetParamP()->Charge() < 0 ){
1638 fHistSwappedV0Counter -> Fill( 1 );
1640 fHistSwappedV0Counter -> Fill( 0 );
1642 if ( fkUseOnTheFly ) CheckChargeV0(v0);
1646 v0->GetPxPyPz( tV0mom[0],tV0mom[1],tV0mom[2] );
1647 Double_t lV0TotalMomentum = TMath::Sqrt(
1648 tV0mom[0]*tV0mom[0]+tV0mom[1]*tV0mom[1]+tV0mom[2]*tV0mom[2] );
1650 Double_t tDecayVertexV0[3]; v0->GetXYZ(tDecayVertexV0[0],tDecayVertexV0[1],tDecayVertexV0[2]);
1651 lV0Radius = TMath::Sqrt(tDecayVertexV0[0]*tDecayVertexV0[0]+tDecayVertexV0[1]*tDecayVertexV0[1]);
1653 lRapK0Short = v0->RapK0Short();
1654 lRapLambda = v0->RapLambda();
1656 //Set Variables for later filling
1657 fTreeVariableV0x = tDecayVertexV0[0];
1658 fTreeVariableV0y = tDecayVertexV0[1];
1659 fTreeVariableV0z = tDecayVertexV0[2];
1661 //Set Variables for later filling
1662 fTreeVariableV0Px = tV0mom[0];
1663 fTreeVariableV0Py = tV0mom[1];
1664 fTreeVariableV0Pz = tV0mom[2];
1666 if ((lPt<fMinV0Pt)||(fMaxV0Pt<lPt)) continue;
1668 UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
1669 UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
1671 Double_t lMomPos[3]; v0->GetPPxPyPz(lMomPos[0],lMomPos[1],lMomPos[2]);
1672 Double_t lMomNeg[3]; v0->GetNPxPyPz(lMomNeg[0],lMomNeg[1],lMomNeg[2]);
1674 AliESDtrack *pTrack=((AliESDEvent*)lESDevent)->GetTrack(lKeyPos);
1675 AliESDtrack *nTrack=((AliESDEvent*)lESDevent)->GetTrack(lKeyNeg);
1676 if (!pTrack || !nTrack) {
1677 Printf("ERROR: Could not retreive one of the daughter track");
1681 fTreeVariableNegEta = nTrack->Eta();
1682 fTreeVariablePosEta = pTrack->Eta();
1684 // Filter like-sign V0 (next: add counter and distribution)
1685 if ( pTrack->GetSign() == nTrack->GetSign()){
1689 //________________________________________________________________________
1690 // Track quality cuts
1691 Float_t lPosTrackCrossedRows = pTrack->GetTPCClusterInfo(2,1);
1692 Float_t lNegTrackCrossedRows = nTrack->GetTPCClusterInfo(2,1);
1693 fTreeVariableLeastNbrCrossedRows = (Int_t) lPosTrackCrossedRows;
1694 if( lNegTrackCrossedRows < fTreeVariableLeastNbrCrossedRows )
1695 fTreeVariableLeastNbrCrossedRows = (Int_t) lNegTrackCrossedRows;
1697 // TPC refit condition (done during reconstruction for Offline but not for On-the-fly)
1698 if( !(pTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
1699 if( !(nTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
1701 if ( ( ( ( pTrack->GetTPCClusterInfo(2,1) ) < 70 ) || ( ( nTrack->GetTPCClusterInfo(2,1) ) < 70 ) )&&(fkTakeAllTracks==kFALSE) ) continue;
1703 //GetKinkIndex condition
1704 if( pTrack->GetKinkIndex(0)>0 || nTrack->GetKinkIndex(0)>0 ) continue;
1706 //Findable clusters > 0 condition
1707 if( pTrack->GetTPCNclsF()<=0 || nTrack->GetTPCNclsF()<=0 ) continue;
1709 //Compute ratio Crossed Rows / Findable clusters
1710 //Note: above test avoids division by zero!
1711 Float_t lPosTrackCrossedRowsOverFindable = -1;
1712 Float_t lNegTrackCrossedRowsOverFindable = -1;
1713 if ( ((double)(pTrack->GetTPCNclsF()) ) != 0 ) lPosTrackCrossedRowsOverFindable = lPosTrackCrossedRows / ((double)(pTrack->GetTPCNclsF()));
1714 if ( ((double)(nTrack->GetTPCNclsF()) ) != 0 ) lNegTrackCrossedRowsOverFindable = lNegTrackCrossedRows / ((double)(nTrack->GetTPCNclsF()));
1716 fTreeVariableLeastRatioCrossedRowsOverFindable = lPosTrackCrossedRowsOverFindable;
1717 if( lNegTrackCrossedRowsOverFindable < fTreeVariableLeastRatioCrossedRowsOverFindable )
1718 fTreeVariableLeastRatioCrossedRowsOverFindable = lNegTrackCrossedRowsOverFindable;
1720 //Lowest Cut Level for Ratio Crossed Rows / Findable = 0.8, set here
1721 if ( (fTreeVariableLeastRatioCrossedRowsOverFindable < 0.8)&&(fkTakeAllTracks==kFALSE) ) continue;
1723 //End track Quality Cuts
1724 //________________________________________________________________________
1726 lDcaPosToPrimVertex = TMath::Abs(pTrack->GetD(lPrimaryVtxPosition[0],
1727 lPrimaryVtxPosition[1],
1730 lDcaNegToPrimVertex = TMath::Abs(nTrack->GetD(lPrimaryVtxPosition[0],
1731 lPrimaryVtxPosition[1],
1734 lOnFlyStatus = v0->GetOnFlyStatus();
1735 lChi2V0 = v0->GetChi2V0();
1736 lDcaV0Daughters = v0->GetDcaV0Daughters();
1737 lDcaV0ToPrimVertex = v0->GetD(lPrimaryVtxPosition[0],lPrimaryVtxPosition[1],lPrimaryVtxPosition[2]);
1738 lV0CosineOfPointingAngle = v0->GetV0CosineOfPointingAngle(lPrimaryVtxPosition[0],lPrimaryVtxPosition[1],lPrimaryVtxPosition[2]);
1739 fTreeVariableV0CosineOfPointingAngle=lV0CosineOfPointingAngle;
1741 // Getting invariant mass infos directly from ESD
1742 v0->ChangeMassHypothesis(310);
1743 lInvMassK0s = v0->GetEffMass();
1744 v0->ChangeMassHypothesis(3122);
1745 lInvMassLambda = v0->GetEffMass();
1746 v0->ChangeMassHypothesis(-3122);
1747 lInvMassAntiLambda = v0->GetEffMass();
1748 lAlphaV0 = v0->AlphaV0();
1749 lPtArmV0 = v0->PtArmV0();
1751 //fTreeVariableOnFlyStatus = lOnFlyStatus;
1752 //fHistV0OnFlyStatus->Fill(lOnFlyStatus);
1754 //===============================================
1755 // Monte Carlo Association starts here
1756 //===============================================
1758 //---> Set Everything to "I don't know" before starting
1760 fTreeVariablePIDPositive = 0;
1761 fTreeVariablePIDNegative = 0;
1763 fTreeVariableIndexStatus = 0;
1764 fTreeVariableIndexStatusMother = 0;
1766 fTreeVariablePtMother = -1;
1767 fTreeVariablePtMC = -1;
1768 fTreeVariableRapMC = -100;
1770 fTreeVariablePID = -1;
1771 fTreeVariablePIDMother = -1;
1773 fTreeVariablePrimaryStatus = 0;
1774 fTreeVariablePrimaryStatusMother = 0;
1775 fTreeVariableV0CreationRadius = -1;
1777 Int_t lblPosV0Dghter = (Int_t) TMath::Abs( pTrack->GetLabel() );
1778 Int_t lblNegV0Dghter = (Int_t) TMath::Abs( nTrack->GetLabel() );
1780 TParticle* mcPosV0Dghter = lMCstack->Particle( lblPosV0Dghter );
1781 TParticle* mcNegV0Dghter = lMCstack->Particle( lblNegV0Dghter );
1783 fTreeVariablePosTransvMomentumMC = mcPosV0Dghter->Pt();
1784 fTreeVariableNegTransvMomentumMC = mcNegV0Dghter->Pt();
1786 Int_t lPIDPositive = mcPosV0Dghter -> GetPdgCode();
1787 Int_t lPIDNegative = mcNegV0Dghter -> GetPdgCode();
1789 fTreeVariablePIDPositive = lPIDPositive;
1790 fTreeVariablePIDNegative = lPIDNegative;
1792 Int_t lblMotherPosV0Dghter = mcPosV0Dghter->GetFirstMother() ;
1793 Int_t lblMotherNegV0Dghter = mcNegV0Dghter->GetFirstMother();
1795 if( lblMotherPosV0Dghter == lblMotherNegV0Dghter && lblMotherPosV0Dghter > -1 ){
1796 //either label is fine, they're equal at this stage
1797 TParticle* pThisV0 = lMCstack->Particle( lblMotherPosV0Dghter );
1798 //Set tree variables
1799 fTreeVariablePID = pThisV0->GetPdgCode(); //PDG Code
1800 fTreeVariablePtMC = pThisV0->Pt(); //Perfect Pt
1802 fTreeVariableIsNonInjected = lMCevent->IsFromBGEvent(lblMotherPosV0Dghter); //Is it?
1803 if (!fTreeVariableIsNonInjected){
1804 if (!(pThisV0->GetFirstMother()<0))
1805 {fTreeVariableIsNonInjected = kTRUE;} // because there are primaries (ALICE definition) not produced in the collision
1808 //Set Variables for later filling
1809 //Be careful: Vx, Vy, Vz: Creation vertex. So decay position is the
1810 //Creation vertex of any one of the daughters!
1811 fTreeVariableMCV0x = mcPosV0Dghter->Vx();
1812 fTreeVariableMCV0y = mcPosV0Dghter->Vy();
1813 fTreeVariableMCV0z = mcPosV0Dghter->Vz();
1815 //Set Variables for later filling
1816 fTreeVariableMCV0Px = pThisV0->Px();
1817 fTreeVariableMCV0Py = pThisV0->Py();
1818 fTreeVariableMCV0Pz = pThisV0->Pz();
1820 //Only Interested if it's a Lambda, AntiLambda or K0s
1821 //Avoid the Junction Bug! PYTHIA has particles with Px=Py=Pz=E=0 occasionally,
1822 //having particle code 88 (unrecognized by PDG), for documentation purposes.
1823 //Even ROOT's TParticle::Y() is not prepared to deal with that exception!
1824 //Note that TParticle::Pt() is immune (that would just return 0)...
1825 //Though granted that that should be extremely rare in this precise condition...
1826 if( TMath::Abs(fTreeVariablePID) == 3122 || fTreeVariablePID==310 ){
1827 fTreeVariableRapMC = pThisV0->Y(); //Perfect Y
1829 fTreeVariableV0CreationRadius = TMath::Sqrt(
1830 TMath::Power( ( (mcPrimaryVtx.At(0)) - (pThisV0->Vx()) ) , 2) +
1831 TMath::Power( ( (mcPrimaryVtx.At(1)) - (pThisV0->Vy()) ) , 2) +
1832 TMath::Power( ( (mcPrimaryVtx.At(2)) - (pThisV0->Vz()) ) , 2)
1834 if( lblMotherPosV0Dghter < lNbMCPrimary ) fTreeVariableIndexStatus = 1; //looks primary
1835 if( lblMotherPosV0Dghter >= lNbMCPrimary ) fTreeVariableIndexStatus = 2; //looks secondary
1836 if( lMCstack->IsPhysicalPrimary (lblMotherPosV0Dghter) ) fTreeVariablePrimaryStatus = 1; //Is Primary!
1837 if( lMCstack->IsSecondaryFromWeakDecay(lblMotherPosV0Dghter) ) fTreeVariablePrimaryStatus = 2; //Weak Decay!
1838 if( lMCstack->IsSecondaryFromMaterial (lblMotherPosV0Dghter) ) fTreeVariablePrimaryStatus = 3; //Material Int!
1840 //Now we try to acquire the V0 parent particle, if possible
1841 Int_t lblThisV0Parent = pThisV0->GetFirstMother();
1842 if ( lblThisV0Parent > -1 ){ //if it has a parent, get it and store specs
1843 TParticle* pThisV0Parent = lMCstack->Particle( lblThisV0Parent );
1844 fTreeVariablePIDMother = pThisV0Parent->GetPdgCode(); //V0 Mother PDG
1845 fTreeVariablePtMother = pThisV0Parent->Pt(); //V0 Mother Pt
1846 //Primary Status for the V0 Mother particle
1847 if( lblThisV0Parent < lNbMCPrimary ) fTreeVariableIndexStatusMother = 1; //looks primary
1848 if( lblThisV0Parent >= lNbMCPrimary ) fTreeVariableIndexStatusMother = 2; //looks secondary
1849 if( lMCstack->IsPhysicalPrimary (lblThisV0Parent) ) fTreeVariablePrimaryStatusMother = 1; //Is Primary!
1850 if( lMCstack->IsSecondaryFromWeakDecay(lblThisV0Parent) ) fTreeVariablePrimaryStatusMother = 2; //Weak Decay!
1851 if( lMCstack->IsSecondaryFromMaterial (lblThisV0Parent) ) fTreeVariablePrimaryStatusMother = 3; //Material Int!
1855 fTreeVariablePt = v0->Pt();
1856 fTreeVariableChi2V0 = lChi2V0;
1857 fTreeVariableDcaV0ToPrimVertex = lDcaV0ToPrimVertex;
1858 fTreeVariableDcaV0Daughters = lDcaV0Daughters;
1859 fTreeVariableV0CosineOfPointingAngle = lV0CosineOfPointingAngle;
1860 fTreeVariableV0Radius = lV0Radius;
1861 fTreeVariableDcaPosToPrimVertex = lDcaPosToPrimVertex;
1862 fTreeVariableDcaNegToPrimVertex = lDcaNegToPrimVertex;
1863 fTreeVariableInvMassK0s = lInvMassK0s;
1864 fTreeVariableInvMassLambda = lInvMassLambda;
1865 fTreeVariableInvMassAntiLambda = lInvMassAntiLambda;
1866 fTreeVariableRapK0Short = lRapK0Short;
1868 fTreeVariableRapLambda = lRapLambda;
1869 fTreeVariableAlphaV0 = lAlphaV0;
1870 fTreeVariablePtArmV0 = lPtArmV0;
1872 //Official means of acquiring N-sigmas
1873 fTreeVariableNSigmasPosProton = fPIDResponse->NumberOfSigmasTPC( pTrack, AliPID::kProton );
1874 fTreeVariableNSigmasPosPion = fPIDResponse->NumberOfSigmasTPC( pTrack, AliPID::kPion );
1875 fTreeVariableNSigmasNegProton = fPIDResponse->NumberOfSigmasTPC( nTrack, AliPID::kProton );
1876 fTreeVariableNSigmasNegPion = fPIDResponse->NumberOfSigmasTPC( nTrack, AliPID::kPion );
1878 //tDecayVertexV0[0],tDecayVertexV0[1],tDecayVertexV0[2]
1879 Double_t lDistanceTravelled = TMath::Sqrt(
1880 TMath::Power( tDecayVertexV0[0] - lBestPrimaryVtxPos[0] , 2) +
1881 TMath::Power( tDecayVertexV0[1] - lBestPrimaryVtxPos[1] , 2) +
1882 TMath::Power( tDecayVertexV0[2] - lBestPrimaryVtxPos[2] , 2)
1884 fTreeVariableDistOverTotMom = 1e+5;
1885 if( lV0TotalMomentum + 1e-10 != 0 ) fTreeVariableDistOverTotMom = lDistanceTravelled / (lV0TotalMomentum + 1e-10); //avoid division by zero, to be sure
1887 Double_t lMomentumPosTemp[3];
1888 pTrack->GetPxPyPz(lMomentumPosTemp);
1889 Double_t lPtPosTemporary = sqrt(pow(lMomentumPosTemp[0],2) + pow(lMomentumPosTemp[1],2));
1891 Double_t lMomentumNegTemp[3];
1892 nTrack->GetPxPyPz(lMomentumNegTemp);
1893 Double_t lPtNegTemporary = sqrt(pow(lMomentumNegTemp[0],2) + pow(lMomentumNegTemp[1],2));
1895 fTreeVariablePosTransvMomentum = lPtPosTemporary;
1896 fTreeVariableNegTransvMomentum = lPtNegTemporary;
1899 //------------------------------------------------
1901 //------------------------------------------------
1903 // The conditionals are meant to decrease excessive
1906 //Modified version: Keep only OnFlyStatus == 0
1907 //Keep only if included in a parametric InvMass Region 20 sigmas away from peak
1909 //First Selection: Reject OnFly
1910 if( (lOnFlyStatus == 0 && fkUseOnTheFly == kFALSE) || (lOnFlyStatus != 0 && fkUseOnTheFly == kTRUE ) ){
1911 //Second Selection: rough 20-sigma band, parametric.
1912 //K0Short: Enough to parametrize peak broadening with linear function.
1913 Double_t lUpperLimitK0Short = (5.63707e-01) + (1.14979e-02)*fTreeVariablePt;
1914 Double_t lLowerLimitK0Short = (4.30006e-01) - (1.10029e-02)*fTreeVariablePt;
1915 //Lambda: Linear (for higher pt) plus exponential (for low-pt broadening)
1916 //[0]+[1]*x+[2]*TMath::Exp(-[3]*x)
1917 Double_t lUpperLimitLambda = (1.13688e+00) + (5.27838e-03)*fTreeVariablePt + (8.42220e-02)*TMath::Exp(-(3.80595e+00)*fTreeVariablePt);
1918 Double_t lLowerLimitLambda = (1.09501e+00) - (5.23272e-03)*fTreeVariablePt - (7.52690e-02)*TMath::Exp(-(3.46339e+00)*fTreeVariablePt);
1920 if( (fTreeVariableInvMassLambda < lUpperLimitLambda && fTreeVariableInvMassLambda > lLowerLimitLambda ) ||
1921 (fTreeVariableInvMassAntiLambda < lUpperLimitLambda && fTreeVariableInvMassAntiLambda > lLowerLimitLambda ) ||
1922 (fTreeVariableInvMassK0s < lUpperLimitK0Short && fTreeVariableInvMassK0s > lLowerLimitK0Short ) ){
1923 //Pre-selection in case this is AA...
1924 if( fkIsNuclear == kFALSE ) fTree->Fill();
1925 if( fkIsNuclear == kTRUE){
1926 //If this is a nuclear collision___________________
1927 // ... pre-filter with daughter eta selection only (not TPC)
1928 if ( TMath::Abs(fTreeVariableNegEta)<0.8 && TMath::Abs(fTreeVariablePosEta)<0.8 ) fTree->Fill();
1929 }//end nuclear_____________________________________
1933 //------------------------------------------------
1935 //------------------------------------------------
1938 }// This is the end of the V0 loop
1940 //------------------------------------------------
1942 // Post output data.
1943 PostData(1, fListHistV0);
1947 //________________________________________________________________________
1948 void AliAnalysisTaskExtractPerformanceV0::Terminate(Option_t *)
1950 // Draw result to the screen
1951 // Called once at the end of the query
1953 TList *cRetrievedList = 0x0;
1954 cRetrievedList = (TList*)GetOutputData(1);
1955 if(!cRetrievedList){
1956 Printf("ERROR - AliAnalysisTaskExtractV0 : ouput data container list not available\n");
1960 fHistV0MultiplicityForTrigEvt = dynamic_cast<TH1F*> ( cRetrievedList->FindObject("fHistV0MultiplicityForTrigEvt") );
1961 if (!fHistV0MultiplicityForTrigEvt) {
1962 Printf("ERROR - AliAnalysisTaskExtractV0 : fHistV0MultiplicityForTrigEvt not available");
1966 TCanvas *canCheck = new TCanvas("AliAnalysisTaskExtractV0","V0 Multiplicity",10,10,510,510);
1967 canCheck->cd(1)->SetLogy();
1969 fHistV0MultiplicityForTrigEvt->SetMarkerStyle(22);
1970 fHistV0MultiplicityForTrigEvt->DrawCopy("E");
1973 //----------------------------------------------------------------------------
1975 Double_t AliAnalysisTaskExtractPerformanceV0::MyRapidity(Double_t rE, Double_t rPz) const
1977 // Local calculation for rapidity
1978 Double_t ReturnValue = -100;
1979 if( (rE-rPz+1.e-13) != 0 && (rE+rPz) != 0 ){
1980 ReturnValue = 0.5*TMath::Log((rE+rPz)/(rE-rPz+1.e-13));
1985 //________________________________________________________________________
1986 void AliAnalysisTaskExtractPerformanceV0::CheckChargeV0(AliESDv0 *v0)
1988 // This function checks charge of negative and positive daughter tracks.
1989 // If incorrectly defined (onfly vertexer), swaps out.
1990 if( v0->GetParamN()->Charge() > 0 && v0->GetParamP()->Charge() < 0 ){
1991 //V0 daughter track swapping is required! Note: everything is swapped here... P->N, N->P
1992 Long_t lCorrectNidx = v0->GetPindex();
1993 Long_t lCorrectPidx = v0->GetNindex();
1994 Double32_t lCorrectNmom[3];
1995 Double32_t lCorrectPmom[3];
1996 v0->GetPPxPyPz( lCorrectNmom[0], lCorrectNmom[1], lCorrectNmom[2] );
1997 v0->GetNPxPyPz( lCorrectPmom[0], lCorrectPmom[1], lCorrectPmom[2] );
1999 AliExternalTrackParam lCorrectParamN(
2000 v0->GetParamP()->GetX() ,
2001 v0->GetParamP()->GetAlpha() ,
2002 v0->GetParamP()->GetParameter() ,
2003 v0->GetParamP()->GetCovariance()
2005 AliExternalTrackParam lCorrectParamP(
2006 v0->GetParamN()->GetX() ,
2007 v0->GetParamN()->GetAlpha() ,
2008 v0->GetParamN()->GetParameter() ,
2009 v0->GetParamN()->GetCovariance()
2011 lCorrectParamN.SetMostProbablePt( v0->GetParamP()->GetMostProbablePt() );
2012 lCorrectParamP.SetMostProbablePt( v0->GetParamN()->GetMostProbablePt() );
2014 //Get Variables___________________________________________________
2015 Double_t lDcaV0Daughters = v0 -> GetDcaV0Daughters();
2016 Double_t lCosPALocal = v0 -> GetV0CosineOfPointingAngle();
2017 Bool_t lOnFlyStatusLocal = v0 -> GetOnFlyStatus();
2019 //Create Replacement Object_______________________________________
2020 AliESDv0 *v0correct = new AliESDv0(lCorrectParamN,lCorrectNidx,lCorrectParamP,lCorrectPidx);
2021 v0correct->SetDcaV0Daughters ( lDcaV0Daughters );
2022 v0correct->SetV0CosineOfPointingAngle ( lCosPALocal );
2023 v0correct->ChangeMassHypothesis ( kK0Short );
2024 v0correct->SetOnFlyStatus ( lOnFlyStatusLocal );
2026 //Reverse Cluster info..._________________________________________
2027 v0correct->SetClusters( v0->GetClusters( 1 ), v0->GetClusters ( 0 ) );
2030 //Proper cleanup..._______________________________________________
2031 v0correct->Delete();
2034 //Just another cross-check and output_____________________________
2035 if( v0->GetParamN()->Charge() > 0 && v0->GetParamP()->Charge() < 0 ) {
2036 AliWarning("Found Swapped Charges, tried to correct but something FAILED!");
2038 //AliWarning("Found Swapped Charges and fixed.");
2040 //________________________________________________________________
2042 //Don't touch it! ---
2043 //Printf("Ah, nice. Charges are already ordered...");