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 // This task is meant to explore the possibility of using a VZERO amplitude
19 // based multiplicity estimator for proton-proton collisions. For this, two
20 // main operation methods for this task are foreseen:
22 // 1) (under development) it should act as an auxiliary task and provide a
23 // calibrated estimator
25 // 2) "Debug mode" which will also create a ROOT TTree object with event
26 // by event info potentially used for exploration / calibration. This
27 // includes the following info:
29 // --- All VZERO Amplitudes (saved as Float_t)
30 // --- (optional) time for each channel
31 // --- (optional) time width for each channel
32 // --- GetReferenceMultiplicity Estimator, |eta|<0.5
33 // --- GetReferenceMultiplicity Estimator, |eta|<0.8
34 // --- (if MC) True Multiplicity, |eta|<0.5
35 // --- (if MC) True Multiplicity, 2.8 < eta < 5.8 (VZEROA region)
36 // --- (if MC) True Multiplicity, -3.7 < eta <-1.7 (VZEROC region)
39 // +-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
45 //class AliMCEventHandler;
54 #include <Riostream.h>
60 #include "THnSparse.h"
67 #include "AliESDEvent.h"
68 #include "AliAODEvent.h"
69 #include "AliV0vertexer.h"
70 #include "AliCascadeVertexer.h"
71 #include "AliESDpid.h"
72 #include "AliESDtrack.h"
73 #include "AliESDtrackCuts.h"
74 #include "AliInputEventHandler.h"
75 #include "AliAnalysisManager.h"
76 #include "AliMCEventHandler.h"
77 #include "AliMCEvent.h"
79 #include "AliCentrality.h"
80 #include "AliPPVsMultUtils.h"
82 #include "AliCFContainer.h"
83 #include "AliMultiplicity.h"
84 #include "AliAODMCParticle.h"
85 #include "AliESDcascade.h"
86 #include "AliAODcascade.h"
87 #include "AliESDUtils.h"
88 #include "AliGenEventHeader.h"
89 #include "AliAnalysisTaskSE.h"
90 #include "AliAnalysisUtils.h"
91 #include "AliAnalysisTaskStrangenessVsMultiplicity.h"
96 ClassImp(AliAnalysisTaskStrangenessVsMultiplicity)
98 AliAnalysisTaskStrangenessVsMultiplicity::AliAnalysisTaskStrangenessVsMultiplicity()
99 : AliAnalysisTaskSE(), fListHist(0), fTreeEvent(0), fTreeV0(0), fTreeCascade(0), fPIDResponse(0), fESDtrackCuts(0), fPPVsMultUtils(0), fUtils(0),
100 fkSaveV0Tree ( kFALSE ),
101 fkSaveCascadeTree ( kTRUE ),
102 fkRunVertexers ( kTRUE ),
103 fkSkipEventSelection( kFALSE ),
104 //---> Variables for fTreeEvent
107 fAmplitude_V0AEq (0),
108 fAmplitude_V0CEq (0),
112 fCentrality_V0AEq(0),
113 fCentrality_V0CEq(0),
114 fCentrality_V0MEq(0),
115 fCustomCentrality_V0M(0),
116 fCustomCentrality_V0MEq(0),
120 fEvSel_HasAtLeastSPDVertex(0),
122 fEvSel_IsNotPileup(0),
123 fEvSel_IsNotPileupMV(0),
124 fEvSel_IsNotPileupInMultBins(0),
126 fEvSel_nTracklets(0),
127 fEvSel_nSPDClusters(0),
129 fEvSel_nSPDPrimVertices(0),
131 fEvSel_nContributors(0),
132 fEvSel_nContributorsPileup(0),
134 //---> Variables for fTreeV0
135 fTreeVariableChi2V0(0),
136 fTreeVariableDcaV0Daughters(0),
137 fTreeVariableDcaV0ToPrimVertex(0),
138 fTreeVariableDcaPosToPrimVertex(0),
139 fTreeVariableDcaNegToPrimVertex(0),
140 fTreeVariableV0CosineOfPointingAngle(0),
141 fTreeVariableV0Radius(0),
143 fTreeVariableRapK0Short(0),
144 fTreeVariableRapLambda(0),
145 fTreeVariableInvMassK0s(0),
146 fTreeVariableInvMassLambda(0),
147 fTreeVariableInvMassAntiLambda(0),
148 fTreeVariableAlphaV0(0),
149 fTreeVariablePtArmV0(0),
150 fTreeVariableNegEta(0),
151 fTreeVariablePosEta(0),
153 fTreeVariableNSigmasPosProton(0),
154 fTreeVariableNSigmasPosPion(0),
155 fTreeVariableNSigmasNegProton(0),
156 fTreeVariableNSigmasNegPion(0),
158 fTreeVariableDistOverTotMom(0),
159 fTreeVariableLeastNbrCrossedRows(0),
160 fTreeVariableLeastRatioCrossedRowsOverFindable(0),
162 fTreeVariableCentV0A(0),
163 fTreeVariableCentV0C(0),
164 fTreeVariableCentV0M(0),
165 fTreeVariableCentV0AEq(0),
166 fTreeVariableCentV0CEq(0),
167 fTreeVariableCentV0MEq(0),
168 fTreeVariableAmpV0A(0),
169 fTreeVariableAmpV0C(0),
170 fTreeVariableAmpV0AEq(0),
171 fTreeVariableAmpV0CEq(0),
172 fTreeVariableRefMultEta8(0),
173 fTreeVariableRefMultEta5(0),
174 fTreeVariableRunNumber(0),
175 //---> Variables for fTreeCascade
176 fTreeCascVarCharge(0),
177 fTreeCascVarMassAsXi(0),
178 fTreeCascVarMassAsOmega(0),
180 fTreeCascVarRapXi(0),
181 fTreeCascVarRapOmega(0),
182 fTreeCascVarNegEta(0),
183 fTreeCascVarPosEta(0),
184 fTreeCascVarBachEta(0),
185 fTreeCascVarDCACascDaughters(0),
186 fTreeCascVarDCABachToPrimVtx(0),
187 fTreeCascVarDCAV0Daughters(0),
188 fTreeCascVarDCAV0ToPrimVtx(0),
189 fTreeCascVarDCAPosToPrimVtx(0),
190 fTreeCascVarDCANegToPrimVtx(0),
191 fTreeCascVarCascCosPointingAngle(0),
192 fTreeCascVarCascRadius(0),
193 fTreeCascVarV0Mass(0),
194 fTreeCascVarV0CosPointingAngle(0),
195 fTreeCascVarV0CosPointingAngleSpecial(0),
196 fTreeCascVarV0Radius(0),
197 fTreeCascVarLeastNbrClusters(0),
198 fTreeCascVarDistOverTotMom(0),
199 fTreeCascVarNegNSigmaPion(0),
200 fTreeCascVarNegNSigmaProton(0),
201 fTreeCascVarPosNSigmaPion(0),
202 fTreeCascVarPosNSigmaProton(0),
203 fTreeCascVarBachNSigmaPion(0),
204 fTreeCascVarBachNSigmaKaon(0),
205 fTreeCascVarCentV0A(0),
206 fTreeCascVarCentV0C(0),
207 fTreeCascVarCentV0M(0),
208 fTreeCascVarCentV0AEq(0),
209 fTreeCascVarCentV0CEq(0),
210 fTreeCascVarCentV0MEq(0),
211 fTreeCascVarAmpV0A(0),
212 fTreeCascVarAmpV0C(0),
213 fTreeCascVarAmpV0AEq(0),
214 fTreeCascVarAmpV0CEq(0),
215 fTreeCascVarRefMultEta8(0),
216 fTreeCascVarRefMultEta5(0),
217 fTreeCascVarRunNumber(0),
220 //------------------------------------------------
226 AliAnalysisTaskStrangenessVsMultiplicity::AliAnalysisTaskStrangenessVsMultiplicity(const char *name)
227 : AliAnalysisTaskSE(name), fListHist(0), fTreeEvent(0), fTreeV0(0), fTreeCascade(0), fPIDResponse(0), fESDtrackCuts(0), fPPVsMultUtils(0), fUtils(0),
228 fkSaveV0Tree ( kFALSE ),
229 fkSaveCascadeTree ( kTRUE ),
230 fkRunVertexers ( kTRUE ),
231 fkSkipEventSelection( kFALSE ),
232 //---> Variables for fTreeEvent
235 fAmplitude_V0AEq (0),
236 fAmplitude_V0CEq (0),
240 fCentrality_V0AEq(0),
241 fCentrality_V0CEq(0),
242 fCentrality_V0MEq(0),
243 fCustomCentrality_V0M(0),
244 fCustomCentrality_V0MEq(0),
248 fEvSel_HasAtLeastSPDVertex(0),
250 fEvSel_IsNotPileup(0),
251 fEvSel_IsNotPileupMV(0),
252 fEvSel_IsNotPileupInMultBins(0),
254 fEvSel_nTracklets(0),
255 fEvSel_nSPDClusters(0),
257 fEvSel_nSPDPrimVertices(0),
259 fEvSel_nContributors(0),
260 fEvSel_nContributorsPileup(0),
262 //---> Variables for fTreeV0
263 fTreeVariableChi2V0(0),
264 fTreeVariableDcaV0Daughters(0),
265 fTreeVariableDcaV0ToPrimVertex(0),
266 fTreeVariableDcaPosToPrimVertex(0),
267 fTreeVariableDcaNegToPrimVertex(0),
268 fTreeVariableV0CosineOfPointingAngle(0),
269 fTreeVariableV0Radius(0),
271 fTreeVariableRapK0Short(0),
272 fTreeVariableRapLambda(0),
273 fTreeVariableInvMassK0s(0),
274 fTreeVariableInvMassLambda(0),
275 fTreeVariableInvMassAntiLambda(0),
276 fTreeVariableAlphaV0(0),
277 fTreeVariablePtArmV0(0),
278 fTreeVariableNegEta(0),
279 fTreeVariablePosEta(0),
281 fTreeVariableNSigmasPosProton(0),
282 fTreeVariableNSigmasPosPion(0),
283 fTreeVariableNSigmasNegProton(0),
284 fTreeVariableNSigmasNegPion(0),
286 fTreeVariableDistOverTotMom(0),
287 fTreeVariableLeastNbrCrossedRows(0),
288 fTreeVariableLeastRatioCrossedRowsOverFindable(0),
290 fTreeVariableCentV0A(0),
291 fTreeVariableCentV0C(0),
292 fTreeVariableCentV0M(0),
293 fTreeVariableCentV0AEq(0),
294 fTreeVariableCentV0CEq(0),
295 fTreeVariableCentV0MEq(0),
296 fTreeVariableAmpV0A(0),
297 fTreeVariableAmpV0C(0),
298 fTreeVariableAmpV0AEq(0),
299 fTreeVariableAmpV0CEq(0),
300 fTreeVariableRefMultEta8(0),
301 fTreeVariableRefMultEta5(0),
302 fTreeVariableRunNumber(0),
303 //---> Variables for fTreeCascade
304 fTreeCascVarCharge(0),
305 fTreeCascVarMassAsXi(0),
306 fTreeCascVarMassAsOmega(0),
308 fTreeCascVarRapXi(0),
309 fTreeCascVarRapOmega(0),
310 fTreeCascVarNegEta(0),
311 fTreeCascVarPosEta(0),
312 fTreeCascVarBachEta(0),
313 fTreeCascVarDCACascDaughters(0),
314 fTreeCascVarDCABachToPrimVtx(0),
315 fTreeCascVarDCAV0Daughters(0),
316 fTreeCascVarDCAV0ToPrimVtx(0),
317 fTreeCascVarDCAPosToPrimVtx(0),
318 fTreeCascVarDCANegToPrimVtx(0),
319 fTreeCascVarCascCosPointingAngle(0),
320 fTreeCascVarCascRadius(0),
321 fTreeCascVarV0Mass(0),
322 fTreeCascVarV0CosPointingAngle(0),
323 fTreeCascVarV0CosPointingAngleSpecial(0),
324 fTreeCascVarV0Radius(0),
325 fTreeCascVarLeastNbrClusters(0),
326 fTreeCascVarDistOverTotMom(0),
327 fTreeCascVarNegNSigmaPion(0),
328 fTreeCascVarNegNSigmaProton(0),
329 fTreeCascVarPosNSigmaPion(0),
330 fTreeCascVarPosNSigmaProton(0),
331 fTreeCascVarBachNSigmaPion(0),
332 fTreeCascVarBachNSigmaKaon(0),
333 fTreeCascVarCentV0A(0),
334 fTreeCascVarCentV0C(0),
335 fTreeCascVarCentV0M(0),
336 fTreeCascVarCentV0AEq(0),
337 fTreeCascVarCentV0CEq(0),
338 fTreeCascVarCentV0MEq(0),
339 fTreeCascVarAmpV0A(0),
340 fTreeCascVarAmpV0C(0),
341 fTreeCascVarAmpV0AEq(0),
342 fTreeCascVarAmpV0CEq(0),
343 fTreeCascVarRefMultEta8(0),
344 fTreeCascVarRefMultEta5(0),
345 fTreeCascVarRunNumber(0),
350 //Re-vertex: Will only apply for cascade candidates
352 fV0VertexerSels[0] = 33. ; // max allowed chi2
353 fV0VertexerSels[1] = 0.02; // min allowed impact parameter for the 1st daughter (LHC09a4 : 0.05)
354 fV0VertexerSels[2] = 0.02; // min allowed impact parameter for the 2nd daughter (LHC09a4 : 0.05)
355 fV0VertexerSels[4] = 0.95; // min allowed cosine of V0's pointing angle (LHC09a4 : 0.99)
356 fV0VertexerSels[5] = 1.0 ; // min radius of the fiducial volume (LHC09a4 : 0.2)
357 fV0VertexerSels[6] = 200. ; // max radius of the fiducial volume (LHC09a4 : 100.0)
359 fCascadeVertexerSels[0] = 33. ; // max allowed chi2 (same as PDC07)
360 fCascadeVertexerSels[1] = 0.05 ; // min allowed V0 impact parameter (PDC07 : 0.05 / LHC09a4 : 0.025 )
361 fCascadeVertexerSels[2] = 0.010; // "window" around the Lambda mass (PDC07 : 0.008 / LHC09a4 : 0.010 )
362 fCascadeVertexerSels[3] = 0.03 ; // min allowed bachelor's impact parameter (PDC07 : 0.035 / LHC09a4 : 0.025 )
363 fCascadeVertexerSels[4] = 2.0 ; // max allowed DCA between the V0 and the bachelor (PDC07 : 0.1 / LHC09a4 : 0.2 )
364 fCascadeVertexerSels[5] = 0.95 ; // min allowed cosine of the cascade pointing angle (PDC07 : 0.9985 / LHC09a4 : 0.998 )
365 fCascadeVertexerSels[6] = 0.4 ; // min radius of the fiducial volume (PDC07 : 0.9 / LHC09a4 : 0.2 )
366 fCascadeVertexerSels[7] = 100. ; // max radius of the fiducial volume (PDC07 : 100 / LHC09a4 : 100 )
369 DefineOutput(1, TList::Class()); // Event Counter Histo
370 DefineOutput(2, TTree::Class()); // Event Tree
371 DefineOutput(3, TTree::Class()); // V0 Tree
372 DefineOutput(4, TTree::Class()); // Cascade Tree
376 AliAnalysisTaskStrangenessVsMultiplicity::~AliAnalysisTaskStrangenessVsMultiplicity()
378 //------------------------------------------------
380 //------------------------------------------------
399 delete fPPVsMultUtils;
400 fPPVsMultUtils = 0x0;
408 //________________________________________________________________________
409 void AliAnalysisTaskStrangenessVsMultiplicity::UserCreateOutputObjects()
415 //------------------------------------------------
417 fTreeEvent = new TTree("fTreeEvent","Event");
419 //------------------------------------------------
420 // fTree Branch definitions - Event by Event info
421 //------------------------------------------------
423 //-----------BASIC-INFO---------------------------
425 //--- VZERO Data (Integrated)
426 fTreeEvent->Branch("fAmplitude_V0A",&fAmplitude_V0A,"fAmplitude_V0A/F");
427 fTreeEvent->Branch("fAmplitude_V0C",&fAmplitude_V0C,"fAmplitude_V0C/F");
428 fTreeEvent->Branch("fAmplitude_V0AEq",&fAmplitude_V0AEq,"fAmplitude_V0AEq/F");
429 fTreeEvent->Branch("fAmplitude_V0CEq",&fAmplitude_V0CEq,"fAmplitude_V0CEq/F");
431 //Info from AliCentrality (not necessarily 'centrality' per se)
432 fTreeEvent->Branch("fCentrality_V0A",&fCentrality_V0A,"fCentrality_V0A/F");
433 fTreeEvent->Branch("fCentrality_V0C",&fCentrality_V0C,"fCentrality_V0C/F");
434 fTreeEvent->Branch("fCentrality_V0M",&fCentrality_V0M,"fCentrality_V0M/F");
435 fTreeEvent->Branch("fCentrality_V0AEq",&fCentrality_V0AEq,"fCentrality_V0AEq/F");
436 fTreeEvent->Branch("fCentrality_V0CEq",&fCentrality_V0CEq,"fCentrality_V0CEq/F");
437 fTreeEvent->Branch("fCentrality_V0MEq",&fCentrality_V0MEq,"fCentrality_V0MEq/F");
438 fTreeEvent->Branch("fCustomCentrality_V0M",&fCustomCentrality_V0M,"fCustomCentrality_V0M/F");
439 fTreeEvent->Branch("fCustomCentrality_V0MEq",&fCustomCentrality_V0MEq,"fCustomCentrality_V0MEq/F");
441 //Official GetReferenceMultiplicity
442 fTreeEvent->Branch("fRefMultEta5",&fRefMultEta5,"fRefMultEta5/I");
443 fTreeEvent->Branch("fRefMultEta8",&fRefMultEta8,"fRefMultEta8/I");
446 fTreeEvent->Branch("fRunNumber", &fRunNumber, "fRunNumber/I");
448 //Booleans for Event Selection
449 fTreeEvent->Branch("fEvSel_HasAtLeastSPDVertex", &fEvSel_HasAtLeastSPDVertex, "fEvSel_HasAtLeastSPDVertex/O");
450 fTreeEvent->Branch("fEvSel_VtxZCut", &fEvSel_VtxZCut, "fEvSel_VtxZCut/O");
451 fTreeEvent->Branch("fEvSel_IsNotPileup", &fEvSel_IsNotPileup, "fEvSel_IsNotPileup/O");
452 fTreeEvent->Branch("fEvSel_IsNotPileupMV", &fEvSel_IsNotPileupMV, "fEvSel_IsNotPileupMV/O");
453 fTreeEvent->Branch("fEvSel_IsNotPileupInMultBins", &fEvSel_IsNotPileupInMultBins, "fEvSel_IsNotPileupInMultBins/O");
454 fTreeEvent->Branch("fEvSel_Triggered", &fEvSel_Triggered, "fEvSel_Triggered/O");
456 //Tracklets vs clusters test
457 fTreeEvent->Branch("fEvSel_nTracklets", &fEvSel_nTracklets, "fEvSel_nTracklets/I");
458 fTreeEvent->Branch("fEvSel_nSPDClusters", &fEvSel_nSPDClusters, "fEvSel_nSPDClusters/I");
460 fTreeEvent->Branch("fEvSel_VtxZ", &fEvSel_VtxZ, "fEvSel_VtxZ/F");
461 fTreeEvent->Branch("fEvSel_nSPDPrimVertices", &fEvSel_nSPDPrimVertices, "fEvSel_nSPDPrimVertices/I");
462 fTreeEvent->Branch("fEvSel_distZ", &fEvSel_distZ, "fEvSel_distZ/F");
463 fTreeEvent->Branch("fEvSel_nContributors", &fEvSel_nContributors, "fEvSel_nContributors/I");
464 fTreeEvent->Branch("fEvSel_nContributorsPileup", &fEvSel_nContributorsPileup, "fEvSel_nContributorsPileup/I");
466 //Create Basic V0 Output Tree
467 fTreeV0 = new TTree( "fTreeV0", "V0 Candidates");
469 //------------------------------------------------
470 // fTreeV0 Branch definitions
471 //------------------------------------------------
473 //-----------BASIC-INFO---------------------------
474 fTreeV0->Branch("fTreeVariableChi2V0",&fTreeVariableChi2V0,"fTreeVariableChi2V0/F");
475 fTreeV0->Branch("fTreeVariableDcaV0Daughters",&fTreeVariableDcaV0Daughters,"fTreeVariableDcaV0Daughters/F");
476 fTreeV0->Branch("fTreeVariableDcaV0ToPrimVertex",&fTreeVariableDcaV0ToPrimVertex,"fTreeVariableDcaV0ToPrimVertex/F");
477 fTreeV0->Branch("fTreeVariableDcaPosToPrimVertex",&fTreeVariableDcaPosToPrimVertex,"fTreeVariableDcaPosToPrimVertex/F");
478 fTreeV0->Branch("fTreeVariableDcaNegToPrimVertex",&fTreeVariableDcaNegToPrimVertex,"fTreeVariableDcaNegToPrimVertex/F");
479 fTreeV0->Branch("fTreeVariableV0Radius",&fTreeVariableV0Radius,"fTreeVariableV0Radius/F");
480 fTreeV0->Branch("fTreeVariablePt",&fTreeVariablePt,"fTreeVariablePt/F");
481 fTreeV0->Branch("fTreeVariableRapK0Short",&fTreeVariableRapK0Short,"fTreeVariableRapK0Short/F");
482 fTreeV0->Branch("fTreeVariableRapLambda",&fTreeVariableRapLambda,"fTreeVariableRapLambda/F");
483 fTreeV0->Branch("fTreeVariableInvMassK0s",&fTreeVariableInvMassK0s,"fTreeVariableInvMassK0s/F");
484 fTreeV0->Branch("fTreeVariableInvMassLambda",&fTreeVariableInvMassLambda,"fTreeVariableInvMassLambda/F");
485 fTreeV0->Branch("fTreeVariableInvMassAntiLambda",&fTreeVariableInvMassAntiLambda,"fTreeVariableInvMassAntiLambda/F");
486 fTreeV0->Branch("fTreeVariableV0CosineOfPointingAngle",&fTreeVariableV0CosineOfPointingAngle,"fTreeVariableV0CosineOfPointingAngle/F");
487 fTreeV0->Branch("fTreeVariableAlphaV0",&fTreeVariableAlphaV0,"fTreeVariableAlphaV0/F");
488 fTreeV0->Branch("fTreeVariablePtArmV0",&fTreeVariablePtArmV0,"fTreeVariablePtArmV0/F");
489 fTreeV0->Branch("fTreeVariableLeastNbrCrossedRows",&fTreeVariableLeastNbrCrossedRows,"fTreeVariableLeastNbrCrossedRows/I");
490 fTreeV0->Branch("fTreeVariableLeastRatioCrossedRowsOverFindable",&fTreeVariableLeastRatioCrossedRowsOverFindable,"fTreeVariableLeastRatioCrossedRowsOverFindable/F");
491 fTreeV0->Branch("fTreeVariableDistOverTotMom",&fTreeVariableDistOverTotMom,"fTreeVariableDistOverTotMom/F");
492 fTreeV0->Branch("fTreeVariableNSigmasPosProton",&fTreeVariableNSigmasPosProton,"fTreeVariableNSigmasPosProton/F");
493 fTreeV0->Branch("fTreeVariableNSigmasPosPion",&fTreeVariableNSigmasPosPion,"fTreeVariableNSigmasPosPion/F");
494 fTreeV0->Branch("fTreeVariableNSigmasNegProton",&fTreeVariableNSigmasNegProton,"fTreeVariableNSigmasNegProton/F");
495 fTreeV0->Branch("fTreeVariableNSigmasNegPion",&fTreeVariableNSigmasNegPion,"fTreeVariableNSigmasNegPion/F");
496 fTreeV0->Branch("fTreeVariableNegEta",&fTreeVariableNegEta,"fTreeVariableNegEta/F");
497 fTreeV0->Branch("fTreeVariablePosEta",&fTreeVariablePosEta,"fTreeVariablePosEta/F");
498 //-----------MULTIPLICITY-INFO--------------------
499 fTreeV0->Branch("fTreeVariableCentV0A",&fTreeVariableCentV0A,"fTreeVariableCentV0A/F");
500 fTreeV0->Branch("fTreeVariableCentV0C",&fTreeVariableCentV0C,"fTreeVariableCentV0C/F");
501 fTreeV0->Branch("fTreeVariableCentV0M",&fTreeVariableCentV0M,"fTreeVariableCentV0M/F");
502 fTreeV0->Branch("fTreeVariableCentV0AEq",&fTreeVariableCentV0AEq,"fTreeVariableCentV0AEq/F");
503 fTreeV0->Branch("fTreeVariableCentV0CEq",&fTreeVariableCentV0CEq,"fTreeVariableCentV0CEq/F");
504 fTreeV0->Branch("fTreeVariableCentV0MEq",&fTreeVariableCentV0MEq,"fTreeVariableCentV0MEq/F");
505 fTreeV0->Branch("fTreeVariableAmpV0A",&fTreeVariableAmpV0A,"fTreeVariableAmpV0A/F");
506 fTreeV0->Branch("fTreeVariableAmpV0C",&fTreeVariableAmpV0C,"fTreeVariableAmpV0C/F");
507 fTreeV0->Branch("fTreeVariableAmpV0AEq",&fTreeVariableAmpV0AEq,"fTreeVariableAmpV0AEq/F");
508 fTreeV0->Branch("fTreeVariableAmpV0CEq",&fTreeVariableAmpV0CEq,"fTreeVariableAmpV0CEq/F");
509 fTreeV0->Branch("fTreeVariableRefMultEta8",&fTreeVariableRefMultEta8,"fTreeVariableRefMultEta8/I");
510 fTreeV0->Branch("fTreeVariableRefMultEta5",&fTreeVariableRefMultEta5,"fTreeVariableRefMultEta5/I");
511 fTreeV0->Branch("fTreeVariableRunNumber",&fTreeVariableRunNumber,"fTreeVariableRunNumber/I");
512 //------------------------------------------------
514 //Create Cascade output tree
515 fTreeCascade = new TTree("fTreeCascade","CascadeCandidates");
517 //------------------------------------------------
518 // fTreeCascade Branch definitions - Cascade Tree
519 //------------------------------------------------
521 //-----------BASIC-INFO---------------------------
522 fTreeCascade->Branch("fTreeCascVarCharge",&fTreeCascVarCharge,"fTreeCascVarCharge/I");
523 fTreeCascade->Branch("fTreeCascVarMassAsXi",&fTreeCascVarMassAsXi,"fTreeCascVarMassAsXi/F");
524 fTreeCascade->Branch("fTreeCascVarMassAsOmega",&fTreeCascVarMassAsOmega,"fTreeCascVarMassAsOmega/F");
525 fTreeCascade->Branch("fTreeCascVarPt",&fTreeCascVarPt,"fTreeCascVarPt/F");
526 fTreeCascade->Branch("fTreeCascVarRapXi",&fTreeCascVarRapXi,"fTreeCascVarRapXi/F");
527 fTreeCascade->Branch("fTreeCascVarRapOmega",&fTreeCascVarRapOmega,"fTreeCascVarRapOmega/F");
528 fTreeCascade->Branch("fTreeCascVarNegEta",&fTreeCascVarNegEta,"fTreeCascVarNegEta/F");
529 fTreeCascade->Branch("fTreeCascVarPosEta",&fTreeCascVarPosEta,"fTreeCascVarPosEta/F");
530 fTreeCascade->Branch("fTreeCascVarBachEta",&fTreeCascVarBachEta,"fTreeCascVarBachEta/F");
531 //-----------INFO-FOR-CUTS------------------------
532 fTreeCascade->Branch("fTreeCascVarDCACascDaughters",&fTreeCascVarDCACascDaughters,"fTreeCascVarDCACascDaughters/F");
533 fTreeCascade->Branch("fTreeCascVarDCABachToPrimVtx",&fTreeCascVarDCABachToPrimVtx,"fTreeCascVarDCABachToPrimVtx/F");
534 fTreeCascade->Branch("fTreeCascVarDCAV0Daughters",&fTreeCascVarDCAV0Daughters,"fTreeCascVarDCAV0Daughters/F");
535 fTreeCascade->Branch("fTreeCascVarDCAV0ToPrimVtx",&fTreeCascVarDCAV0ToPrimVtx,"fTreeCascVarDCAV0ToPrimVtx/F");
536 fTreeCascade->Branch("fTreeCascVarDCAPosToPrimVtx",&fTreeCascVarDCAPosToPrimVtx,"fTreeCascVarDCAPosToPrimVtx/F");
537 fTreeCascade->Branch("fTreeCascVarDCANegToPrimVtx",&fTreeCascVarDCANegToPrimVtx,"fTreeCascVarDCANegToPrimVtx/F");
538 fTreeCascade->Branch("fTreeCascVarCascCosPointingAngle",&fTreeCascVarCascCosPointingAngle,"fTreeCascVarCascCosPointingAngle/F");
539 fTreeCascade->Branch("fTreeCascVarCascRadius",&fTreeCascVarCascRadius,"fTreeCascVarCascRadius/F");
540 fTreeCascade->Branch("fTreeCascVarV0Mass",&fTreeCascVarV0Mass,"fTreeCascVarV0Mass/F");
541 fTreeCascade->Branch("fTreeCascVarV0CosPointingAngle",&fTreeCascVarV0CosPointingAngle,"fTreeCascVarV0CosPointingAngle/F");
542 fTreeCascade->Branch("fTreeCascVarV0CosPointingAngleSpecial",&fTreeCascVarV0CosPointingAngleSpecial,"fTreeCascVarV0CosPointingAngleSpecial/F");
543 fTreeCascade->Branch("fTreeCascVarV0Radius",&fTreeCascVarV0Radius,"fTreeCascVarV0Radius/F");
544 fTreeCascade->Branch("fTreeCascVarLeastNbrClusters",&fTreeCascVarLeastNbrClusters,"fTreeCascVarLeastNbrClusters/I");
545 //-----------MULTIPLICITY-INFO--------------------
546 fTreeCascade->Branch("fTreeCascVarCentV0A",&fTreeCascVarCentV0A,"fTreeCascVarCentV0A/F");
547 fTreeCascade->Branch("fTreeCascVarCentV0C",&fTreeCascVarCentV0C,"fTreeCascVarCentV0C/F");
548 fTreeCascade->Branch("fTreeCascVarCentV0M",&fTreeCascVarCentV0M,"fTreeCascVarCentV0M/F");
549 fTreeCascade->Branch("fTreeCascVarCentV0AEq",&fTreeCascVarCentV0AEq,"fTreeCascVarCentV0AEq/F");
550 fTreeCascade->Branch("fTreeCascVarCentV0CEq",&fTreeCascVarCentV0CEq,"fTreeCascVarCentV0CEq/F");
551 fTreeCascade->Branch("fTreeCascVarCentV0MEq",&fTreeCascVarCentV0MEq,"fTreeCascVarCentV0MEq/F");
552 fTreeCascade->Branch("fTreeCascVarAmpV0A",&fTreeCascVarAmpV0A,"fTreeCascVarAmpV0A/F");
553 fTreeCascade->Branch("fTreeCascVarAmpV0C",&fTreeCascVarAmpV0C,"fTreeCascVarAmpV0C/F");
554 fTreeCascade->Branch("fTreeCascVarAmpV0AEq",&fTreeCascVarAmpV0AEq,"fTreeCascVarAmpV0AEq/F");
555 fTreeCascade->Branch("fTreeCascVarAmpV0CEq",&fTreeCascVarAmpV0CEq,"fTreeCascVarAmpV0CEq/F");
556 fTreeCascade->Branch("fTreeCascVarRefMultEta8",&fTreeCascVarRefMultEta8,"fTreeCascVarRefMultEta8/I");
557 fTreeCascade->Branch("fTreeCascVarRefMultEta5",&fTreeCascVarRefMultEta5,"fTreeCascVarRefMultEta5/I");
558 fTreeCascade->Branch("fTreeCascVarRunNumber",&fTreeCascVarRunNumber,"fTreeCascVarRunNumber/I");
559 //-----------DECAY-LENGTH-INFO--------------------
560 fTreeCascade->Branch("fTreeCascVarDistOverTotMom",&fTreeCascVarDistOverTotMom,"fTreeCascVarDistOverTotMom/F");
561 //------------------------------------------------
562 fTreeCascade->Branch("fTreeCascVarNegNSigmaPion",&fTreeCascVarNegNSigmaPion,"fTreeCascVarNegNSigmaPion/F");
563 fTreeCascade->Branch("fTreeCascVarNegNSigmaProton",&fTreeCascVarNegNSigmaProton,"fTreeCascVarNegNSigmaProton/F");
564 fTreeCascade->Branch("fTreeCascVarPosNSigmaPion",&fTreeCascVarPosNSigmaPion,"fTreeCascVarPosNSigmaPion/F");
565 fTreeCascade->Branch("fTreeCascVarPosNSigmaProton",&fTreeCascVarPosNSigmaProton,"fTreeCascVarPosNSigmaProton/F");
566 fTreeCascade->Branch("fTreeCascVarBachNSigmaPion",&fTreeCascVarBachNSigmaPion,"fTreeCascVarBachNSigmaPion/F");
567 fTreeCascade->Branch("fTreeCascVarBachNSigmaKaon",&fTreeCascVarBachNSigmaKaon,"fTreeCascVarBachNSigmaKaon/F");
569 //------------------------------------------------
570 // Particle Identification Setup
571 //------------------------------------------------
573 AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
574 AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
575 fPIDResponse = inputHandler->GetPIDResponse();
578 if(! fESDtrackCuts ){
579 fESDtrackCuts = new AliESDtrackCuts();
582 if(! fPPVsMultUtils ){
583 fPPVsMultUtils = new AliPPVsMultUtils();
587 fUtils = new AliAnalysisUtils();
590 //------------------------------------------------
591 // V0 Multiplicity Histograms
592 //------------------------------------------------
596 fListHist = new TList();
597 fListHist->SetOwner(); // See http://root.cern.ch/root/html/TCollection.html#TCollection:SetOwner
599 if(! fHistEventCounter ) {
600 //Histogram Output: Event-by-Event
601 fHistEventCounter = new TH1D( "fHistEventCounter", ";Evt. Sel. Step;Count",5,0,5);
602 fHistEventCounter->GetXaxis()->SetBinLabel(1, "Processed");
603 fHistEventCounter->GetXaxis()->SetBinLabel(2, "Phys-Sel");
604 fHistEventCounter->GetXaxis()->SetBinLabel(3, "Has Vtx");
605 fHistEventCounter->GetXaxis()->SetBinLabel(4, "Vtx |z|<10cm");
606 fHistEventCounter->GetXaxis()->SetBinLabel(5, "Isn't Pileup");
607 fListHist->Add(fHistEventCounter);
610 //List of Histograms: Normal
611 PostData(1, fListHist);
613 //TTree Object: Saved to base directory. Should cache to disk while saving.
614 //(Important to avoid excessive memory usage, particularly when merging)
615 PostData(2, fTreeEvent);
616 PostData(3, fTreeV0);
617 PostData(4, fTreeCascade);
619 }// end UserCreateOutputObjects
622 //________________________________________________________________________
623 void AliAnalysisTaskStrangenessVsMultiplicity::UserExec(Option_t *)
626 // Called for each event
628 AliESDEvent *lESDevent = 0x0;
630 //Zero all booleans, etc: safe initialization per event
631 fEvSel_HasAtLeastSPDVertex = kFALSE;
632 fEvSel_VtxZCut = kFALSE;
633 fEvSel_IsNotPileup = kFALSE;
634 fEvSel_IsNotPileupMV = kFALSE;
635 fEvSel_IsNotPileupInMultBins = kFALSE;
636 fEvSel_Triggered = kFALSE;
638 fEvSel_nTracklets = -1;
639 fEvSel_nSPDClusters = -1;
640 fEvSel_nContributors = -1;
641 fEvSel_nContributorsPileup = -1;
642 fEvSel_nSPDPrimVertices = -1;
646 // Connect to the InputEvent
647 // After these lines, we should have an ESD/AOD event + the number of V0s in it.
649 // Appropriate for ESD analysis!
651 lESDevent = dynamic_cast<AliESDEvent*>( InputEvent() );
653 AliWarning("ERROR: lESDevent not available \n");
657 //Get VZERO Information for multiplicity later
658 AliVVZERO* esdV0 = lESDevent->GetVZEROData();
660 AliError("AliVVZERO not available");
664 fRunNumber = lESDevent->GetRunNumber();
666 Double_t lMagneticField = -10;
667 lMagneticField = lESDevent->GetMagneticField( );
669 //------------------------------------------------
670 // Variable Definition
671 //------------------------------------------------
673 //------------------------------------------------
675 //------------------------------------------------
677 fHistEventCounter->Fill(0.5);
679 UInt_t maskIsSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
680 Bool_t isSelected = 0;
681 isSelected = (maskIsSelected & AliVEvent::kMB) == AliVEvent::kMB;
682 fEvSel_Triggered = isSelected;
684 //Standard Min-Bias Selection - always do this!
685 if ( (! isSelected) && (! fkSkipEventSelection ) ) {
686 PostData(1, fListHist);
687 PostData(2, fTreeEvent);
688 PostData(3, fTreeV0);
689 PostData(4, fTreeCascade);
693 fHistEventCounter->Fill(1.5);
695 //------------------------------------------------
696 // Primary Vertex Requirements Section:
697 // ---> pp: has vertex, |z|<10cm
698 //------------------------------------------------
700 //classical Proton-proton like selection
701 const AliESDVertex *lPrimaryBestESDVtx = lESDevent->GetPrimaryVertex();
702 const AliESDVertex *lPrimaryTrackingESDVtx = lESDevent->GetPrimaryVertexTracks();
703 const AliESDVertex *lPrimarySPDVtx = lESDevent->GetPrimaryVertexSPD();
705 Double_t lBestPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
706 lPrimaryBestESDVtx->GetXYZ( lBestPrimaryVtxPos );
708 //Only accept if Tracking or SPD vertex is fine
709 if (!lPrimarySPDVtx->GetStatus() && !lPrimaryTrackingESDVtx->GetStatus() && !fkSkipEventSelection ){
710 AliWarning("Pb / No SPD prim. vertex nor prim. Tracking vertex ... return !");
711 PostData(1, fListHist);
712 PostData(2, fTreeEvent);
713 PostData(3, fTreeV0);
714 PostData(4, fTreeCascade);
718 if(! (!lPrimarySPDVtx->GetStatus() && !lPrimaryTrackingESDVtx->GetStatus()) ){
720 fEvSel_HasAtLeastSPDVertex = kTRUE;
723 //Has SPD or Tracking Vertex
724 fHistEventCounter -> Fill(2.5);
726 //Always do Primary Vertex Selection
727 if(TMath::Abs(lBestPrimaryVtxPos[2]) > 10.0 && !fkSkipEventSelection ) {
728 AliWarning("Pb / | Z position of Best Prim Vtx | > 10.0 cm ... return !");
729 PostData(1, fListHist);
730 PostData(2, fTreeEvent);
731 PostData(3, fTreeV0);
732 PostData(4, fTreeCascade);
736 if(TMath::Abs(lBestPrimaryVtxPos[2]) <= 10.0 ){
738 fEvSel_VtxZCut = kTRUE;
740 fEvSel_VtxZ = lBestPrimaryVtxPos[2] ; //Set
742 //Fill Event selected counter
743 fHistEventCounter -> Fill(3.5);
745 //------------------------------------------------
746 // Check if this isn't pileup
747 //------------------------------------------------
749 if(lESDevent->IsPileupFromSPDInMultBins() && !fkSkipEventSelection ){
750 // minContributors=3, minZdist=0.8, nSigmaZdist=3., nSigmaDiamXY=2., nSigmaDiamZ=5.
751 //-> see http://alisoft.cern.ch/viewvc/trunk/STEER/AliESDEvent.h?root=AliRoot&r1=41914&r2=42199&pathrev=42199
752 AliWarning("Pb / Event tagged as pile-up by SPD... return !");
753 PostData(1, fListHist);
754 PostData(2, fTreeEvent);
755 PostData(3, fTreeV0);
756 PostData(4, fTreeCascade);
760 if( !lESDevent->IsPileupFromSPD() ) fEvSel_IsNotPileup = kTRUE;
761 if( !lESDevent->IsPileupFromSPDInMultBins() ) fEvSel_IsNotPileupInMultBins = kTRUE;
763 //Acquire information to compute residual pileup
764 fEvSel_nSPDPrimVertices = lESDevent->GetNumberOfPileupVerticesSPD();
766 //Long_t lNcontributorsSPDvtx = lPrimarySPDVtx -> GetNContributors();
767 Long_t lNcontributorsSecondLargest = -1;
768 Long_t lIndexSecondLargest = -1;
769 //Look for the two events with the largest numbers of contributors...
770 for(Int_t i=0; i<fEvSel_nSPDPrimVertices;i++){
771 const AliESDVertex* pv=lESDevent -> GetPileupVertexSPD(i);
772 if( pv->GetNContributors() > lNcontributorsSecondLargest ){
773 lNcontributorsSecondLargest = pv->GetNContributors();
774 lIndexSecondLargest = i;
777 fEvSel_nContributors = lPrimaryBestESDVtx -> GetNContributors();
778 if( fEvSel_nSPDPrimVertices > 0 && lIndexSecondLargest > -1){
779 const AliESDVertex* largestpv=lESDevent ->GetPileupVertexSPD(lIndexSecondLargest);
780 fEvSel_distZ = lPrimarySPDVtx->GetZ() - largestpv->GetZ();
781 fEvSel_nContributorsPileup = largestpv -> GetNContributors();
784 //First implementation of pileup from multi-vertexer (simple use of analysis utils)
785 if ( !fUtils->IsPileUpMV( lESDevent ) ) fEvSel_IsNotPileupMV = kTRUE;
787 //Fill Event isn't pileup counter
788 fHistEventCounter -> Fill(4.5);
790 //------------------------------------------------
791 // Multiplicity Information Acquistion
792 //------------------------------------------------
794 //Standard GetReferenceMultiplicity Estimator (0.5 and 0.8)
795 fRefMultEta5 = fESDtrackCuts->GetReferenceMultiplicity(lESDevent, AliESDtrackCuts::kTrackletsITSTPC,0.5);
796 fRefMultEta8 = fESDtrackCuts->GetReferenceMultiplicity(lESDevent, AliESDtrackCuts::kTrackletsITSTPC,0.8);
799 Float_t multV0A = 0; // multiplicity from V0 reco side A
800 Float_t multV0C = 0; // multiplicity from V0 reco side C
801 Float_t multV0AEq = 0; // multiplicity from V0 reco side A
802 Float_t multV0CEq = 0; // multiplicity from V0 reco side C
803 Float_t multV0ACorr = 0; // multiplicity from V0 reco side A
804 Float_t multV0CCorr = 0; // multiplicity from V0 reco side C
806 //Non-Equalized Signal: copy of multV0ACorr and multV0CCorr from AliCentralitySelectionTask
807 //Getters for uncorrected multiplicity
808 multV0A=esdV0->GetMTotV0A();
809 multV0C=esdV0->GetMTotV0C();
811 //Get Z vertex position of SPD vertex (why not Tracking if available?)
812 Float_t zvtx = lPrimarySPDVtx->GetZ();
814 //Acquire Corrected multV0A
815 multV0ACorr = AliESDUtils::GetCorrV0A(multV0A,zvtx);
816 multV0CCorr = AliESDUtils::GetCorrV0C(multV0C,zvtx);
818 //Copy to Event Tree for extra information
819 fAmplitude_V0A = multV0ACorr;
820 fAmplitude_V0C = multV0CCorr;
822 // Equalized signals // From AliCentralitySelectionTask // Updated
823 for(Int_t iCh = 32; iCh < 64; ++iCh) {
824 Double_t mult = lESDevent->GetVZEROEqMultiplicity(iCh);
827 for(Int_t iCh = 0; iCh < 32; ++iCh) {
828 Double_t mult = lESDevent->GetVZEROEqMultiplicity(iCh);
831 fAmplitude_V0AEq = multV0AEq;
832 fAmplitude_V0CEq = multV0CEq;
834 fCentrality_V0A = -100;
835 fCentrality_V0C = -100;
836 fCentrality_V0M = -100;
837 fCentrality_V0AEq = -100;
838 fCentrality_V0CEq = -100;
839 fCentrality_V0MEq = -100;
841 //AliCentrality... Check if working?
842 AliCentrality* centrality;
843 centrality = lESDevent->GetCentrality();
844 if ( !(centrality->GetQuality()>1) ){
845 fCentrality_V0A = centrality->GetCentralityPercentile( "V0A" );
846 fCentrality_V0C = centrality->GetCentralityPercentile( "V0C" );
847 fCentrality_V0M = centrality->GetCentralityPercentile( "V0M" );
848 fCentrality_V0AEq = centrality->GetCentralityPercentile( "V0AEq" );
849 fCentrality_V0CEq = centrality->GetCentralityPercentile( "V0CEq" );
850 fCentrality_V0MEq = centrality->GetCentralityPercentile( "V0MEq" );
854 fCustomCentrality_V0M = fPPVsMultUtils->GetMultiplicityPercentile(lESDevent,"V0M");
855 fCustomCentrality_V0MEq = fPPVsMultUtils->GetMultiplicityPercentile(lESDevent,"V0MEq");
857 //Tracklets vs Clusters Exploratory data
858 fEvSel_nTracklets = lESDevent->GetMultiplicity()->GetNumberOfTracklets();
859 fEvSel_nSPDClusters = lESDevent->GetNumberOfITSClusters(0) + lESDevent->GetNumberOfITSClusters(1);
864 //STOP HERE if skipping event selections (no point in doing the rest...)
865 if( fkSkipEventSelection ){
866 PostData(1, fListHist);
867 PostData(2, fTreeEvent);
868 PostData(3, fTreeV0);
869 PostData(4, fTreeCascade);
873 //------------------------------------------------
875 //------------------------------------------------
876 // Fill V0 Tree as needed
877 //------------------------------------------------
879 //Variable definition
880 Int_t lOnFlyStatus = 0;// nv0sOn = 0, nv0sOff = 0;
881 Double_t lChi2V0 = 0;
882 Double_t lDcaV0Daughters = 0, lDcaV0ToPrimVertex = 0;
883 Double_t lDcaPosToPrimVertex = 0, lDcaNegToPrimVertex = 0;
884 Double_t lV0CosineOfPointingAngle = 0;
885 Double_t lV0Radius = 0, lPt = 0;
886 Double_t lRapK0Short = 0, lRapLambda = 0;
887 Double_t lInvMassK0s = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
888 Double_t lAlphaV0 = 0, lPtArmV0 = 0;
890 Double_t fMinV0Pt = 0;
891 Double_t fMaxV0Pt = 100;
894 nv0s = lESDevent->GetNumberOfV0s();
896 for (Int_t iV0 = 0; iV0 < nv0s; iV0++) //extra-crazy test
897 {// This is the begining of the V0 loop
898 AliESDv0 *v0 = ((AliESDEvent*)lESDevent)->GetV0(iV0);
901 Double_t tDecayVertexV0[3]; v0->GetXYZ(tDecayVertexV0[0],tDecayVertexV0[1],tDecayVertexV0[2]);
904 v0->GetPxPyPz( tV0mom[0],tV0mom[1],tV0mom[2] );
905 Double_t lV0TotalMomentum = TMath::Sqrt(
906 tV0mom[0]*tV0mom[0]+tV0mom[1]*tV0mom[1]+tV0mom[2]*tV0mom[2] );
908 lV0Radius = TMath::Sqrt(tDecayVertexV0[0]*tDecayVertexV0[0]+tDecayVertexV0[1]*tDecayVertexV0[1]);
911 lRapK0Short = v0->RapK0Short();
912 lRapLambda = v0->RapLambda();
913 if ((lPt<fMinV0Pt)||(fMaxV0Pt<lPt)) continue;
915 UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
916 UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
918 Double_t lMomPos[3]; v0->GetPPxPyPz(lMomPos[0],lMomPos[1],lMomPos[2]);
919 Double_t lMomNeg[3]; v0->GetNPxPyPz(lMomNeg[0],lMomNeg[1],lMomNeg[2]);
921 AliESDtrack *pTrack=((AliESDEvent*)lESDevent)->GetTrack(lKeyPos);
922 AliESDtrack *nTrack=((AliESDEvent*)lESDevent)->GetTrack(lKeyNeg);
923 if (!pTrack || !nTrack) {
924 Printf("ERROR: Could not retreive one of the daughter track");
928 //Daughter Eta for Eta selection, afterwards
929 fTreeVariableNegEta = nTrack->Eta();
930 fTreeVariablePosEta = pTrack->Eta();
932 // Filter like-sign V0 (next: add counter and distribution)
933 if ( pTrack->GetSign() == nTrack->GetSign()){
937 //________________________________________________________________________
938 // Track quality cuts
939 Float_t lPosTrackCrossedRows = pTrack->GetTPCClusterInfo(2,1);
940 Float_t lNegTrackCrossedRows = nTrack->GetTPCClusterInfo(2,1);
941 fTreeVariableLeastNbrCrossedRows = (Int_t) lPosTrackCrossedRows;
942 if( lNegTrackCrossedRows < fTreeVariableLeastNbrCrossedRows )
943 fTreeVariableLeastNbrCrossedRows = (Int_t) lNegTrackCrossedRows;
945 // TPC refit condition (done during reconstruction for Offline but not for On-the-fly)
946 if( !(pTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
947 if( !(nTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
950 if ( ( ( pTrack->GetTPCClusterInfo(2,1) ) < 70 ) || ( ( nTrack->GetTPCClusterInfo(2,1) ) < 70 ) ) continue;
952 //GetKinkIndex condition
953 if( pTrack->GetKinkIndex(0)>0 || nTrack->GetKinkIndex(0)>0 ) continue;
955 //Findable clusters > 0 condition
956 if( pTrack->GetTPCNclsF()<=0 || nTrack->GetTPCNclsF()<=0 ) continue;
958 //Compute ratio Crossed Rows / Findable clusters
959 //Note: above test avoids division by zero!
960 Float_t lPosTrackCrossedRowsOverFindable = lPosTrackCrossedRows / ((double)(pTrack->GetTPCNclsF()));
961 Float_t lNegTrackCrossedRowsOverFindable = lNegTrackCrossedRows / ((double)(nTrack->GetTPCNclsF()));
963 fTreeVariableLeastRatioCrossedRowsOverFindable = lPosTrackCrossedRowsOverFindable;
964 if( lNegTrackCrossedRowsOverFindable < fTreeVariableLeastRatioCrossedRowsOverFindable )
965 fTreeVariableLeastRatioCrossedRowsOverFindable = lNegTrackCrossedRowsOverFindable;
967 //Lowest Cut Level for Ratio Crossed Rows / Findable = 0.8, set here
968 if ( fTreeVariableLeastRatioCrossedRowsOverFindable < 0.8 ) continue;
970 //End track Quality Cuts
971 //________________________________________________________________________
973 lDcaPosToPrimVertex = TMath::Abs(pTrack->GetD(lBestPrimaryVtxPos[0],
974 lBestPrimaryVtxPos[1],
977 lDcaNegToPrimVertex = TMath::Abs(nTrack->GetD(lBestPrimaryVtxPos[0],
978 lBestPrimaryVtxPos[1],
981 lOnFlyStatus = v0->GetOnFlyStatus();
982 lChi2V0 = v0->GetChi2V0();
983 lDcaV0Daughters = v0->GetDcaV0Daughters();
984 lDcaV0ToPrimVertex = v0->GetD(lBestPrimaryVtxPos[0],lBestPrimaryVtxPos[1],lBestPrimaryVtxPos[2]);
985 lV0CosineOfPointingAngle = v0->GetV0CosineOfPointingAngle(lBestPrimaryVtxPos[0],lBestPrimaryVtxPos[1],lBestPrimaryVtxPos[2]);
986 fTreeVariableV0CosineOfPointingAngle=lV0CosineOfPointingAngle;
988 // Getting invariant mass infos directly from ESD
989 v0->ChangeMassHypothesis(310);
990 lInvMassK0s = v0->GetEffMass();
991 v0->ChangeMassHypothesis(3122);
992 lInvMassLambda = v0->GetEffMass();
993 v0->ChangeMassHypothesis(-3122);
994 lInvMassAntiLambda = v0->GetEffMass();
995 lAlphaV0 = v0->AlphaV0();
996 lPtArmV0 = v0->PtArmV0();
998 fTreeVariablePt = v0->Pt();
999 fTreeVariableChi2V0 = lChi2V0;
1000 fTreeVariableDcaV0ToPrimVertex = lDcaV0ToPrimVertex;
1001 fTreeVariableDcaV0Daughters = lDcaV0Daughters;
1002 fTreeVariableV0CosineOfPointingAngle = lV0CosineOfPointingAngle;
1003 fTreeVariableV0Radius = lV0Radius;
1004 fTreeVariableDcaPosToPrimVertex = lDcaPosToPrimVertex;
1005 fTreeVariableDcaNegToPrimVertex = lDcaNegToPrimVertex;
1006 fTreeVariableInvMassK0s = lInvMassK0s;
1007 fTreeVariableInvMassLambda = lInvMassLambda;
1008 fTreeVariableInvMassAntiLambda = lInvMassAntiLambda;
1009 fTreeVariableRapK0Short = lRapK0Short;
1010 fTreeVariableRapLambda = lRapLambda;
1011 fTreeVariableAlphaV0 = lAlphaV0;
1012 fTreeVariablePtArmV0 = lPtArmV0;
1014 //Official means of acquiring N-sigmas
1015 fTreeVariableNSigmasPosProton = fPIDResponse->NumberOfSigmasTPC( pTrack, AliPID::kProton );
1016 fTreeVariableNSigmasPosPion = fPIDResponse->NumberOfSigmasTPC( pTrack, AliPID::kPion );
1017 fTreeVariableNSigmasNegProton = fPIDResponse->NumberOfSigmasTPC( nTrack, AliPID::kProton );
1018 fTreeVariableNSigmasNegPion = fPIDResponse->NumberOfSigmasTPC( nTrack, AliPID::kPion );
1020 //This requires an Invariant Mass Hypothesis afterwards
1021 fTreeVariableDistOverTotMom = TMath::Sqrt(
1022 TMath::Power( tDecayVertexV0[0] - lBestPrimaryVtxPos[0] , 2) +
1023 TMath::Power( tDecayVertexV0[1] - lBestPrimaryVtxPos[1] , 2) +
1024 TMath::Power( tDecayVertexV0[2] - lBestPrimaryVtxPos[2] , 2)
1026 fTreeVariableDistOverTotMom /= (lV0TotalMomentum+1e-10); //avoid division by zero, to be sure
1028 //Copy Multiplicity information
1029 fTreeVariableCentV0A = fCentrality_V0A;
1030 fTreeVariableCentV0C = fCentrality_V0C;
1031 fTreeVariableCentV0M = fCentrality_V0M;
1032 fTreeVariableCentV0AEq = fCentrality_V0AEq;
1033 fTreeVariableCentV0CEq = fCentrality_V0CEq;
1034 fTreeVariableCentV0MEq = fCentrality_V0MEq;
1035 fTreeVariableAmpV0A = fAmplitude_V0A;
1036 fTreeVariableAmpV0C = fAmplitude_V0C;
1037 fTreeVariableAmpV0AEq = fAmplitude_V0AEq;
1038 fTreeVariableAmpV0CEq = fAmplitude_V0CEq;
1039 fTreeVariableRefMultEta8 = fRefMultEta8;
1040 fTreeVariableRefMultEta5 = fRefMultEta5;
1041 fTreeVariableRunNumber = fRunNumber;
1043 //------------------------------------------------
1045 //------------------------------------------------
1047 // The conditionals are meant to decrease excessive
1050 //First Selection: Reject OnFly
1051 if( lOnFlyStatus == 0 ){
1052 //Second Selection: rough 20-sigma band, parametric.
1053 //K0Short: Enough to parametrize peak broadening with linear function.
1054 Double_t lUpperLimitK0Short = (5.63707e-01) + (1.14979e-02)*fTreeVariablePt;
1055 Double_t lLowerLimitK0Short = (4.30006e-01) - (1.10029e-02)*fTreeVariablePt;
1056 //Lambda: Linear (for higher pt) plus exponential (for low-pt broadening)
1057 //[0]+[1]*x+[2]*TMath::Exp(-[3]*x)
1058 Double_t lUpperLimitLambda = (1.13688e+00) + (5.27838e-03)*fTreeVariablePt + (8.42220e-02)*TMath::Exp(-(3.80595e+00)*fTreeVariablePt);
1059 Double_t lLowerLimitLambda = (1.09501e+00) - (5.23272e-03)*fTreeVariablePt - (7.52690e-02)*TMath::Exp(-(3.46339e+00)*fTreeVariablePt);
1061 if( (fTreeVariableInvMassLambda < lUpperLimitLambda && fTreeVariableInvMassLambda > lLowerLimitLambda ) ||
1062 (fTreeVariableInvMassAntiLambda < lUpperLimitLambda && fTreeVariableInvMassAntiLambda > lLowerLimitLambda ) ||
1063 (fTreeVariableInvMassK0s < lUpperLimitK0Short && fTreeVariableInvMassK0s > lLowerLimitK0Short ) ){
1064 //Pre-selection in case this is AA...
1065 if ( TMath::Abs(fTreeVariableNegEta)<0.8 && TMath::Abs(fTreeVariablePosEta)<0.8 && fkSaveV0Tree ) fTreeV0->Fill();
1068 }// This is the end of the V0 loop
1070 //------------------------------------------------
1071 // Fill V0 tree over.
1072 //------------------------------------------------
1076 //------------------------------------------------
1077 // Rerun cascade vertexer!
1078 //------------------------------------------------
1080 if( fkRunVertexers ){
1081 lESDevent->ResetCascades();
1082 lESDevent->ResetV0s();
1084 AliV0vertexer lV0vtxer;
1085 AliCascadeVertexer lCascVtxer;
1087 lV0vtxer.SetDefaultCuts(fV0VertexerSels);
1088 lCascVtxer.SetDefaultCuts(fCascadeVertexerSels);
1090 lV0vtxer.Tracks2V0vertices(lESDevent);
1091 lCascVtxer.V0sTracks2CascadeVertices(lESDevent);
1094 //------------------------------------------------
1095 // MAIN CASCADE LOOP STARTS HERE
1096 //------------------------------------------------
1097 // Code Credit: Antonin Maire (thanks^100)
1098 // ---> This is an adaptation
1100 Long_t ncascades = 0;
1101 ncascades = lESDevent->GetNumberOfCascades();
1103 for (Int_t iXi = 0; iXi < ncascades; iXi++){
1104 //------------------------------------------------
1106 //------------------------------------------------
1107 //Double_t lTrkgPrimaryVtxRadius3D = -500.0;
1108 //Double_t lBestPrimaryVtxRadius3D = -500.0;
1110 // - 1st part of initialisation : variables needed to store AliESDCascade data members
1111 Double_t lEffMassXi = 0. ;
1112 //Double_t lChi2Xi = -1. ;
1113 Double_t lDcaXiDaughters = -1. ;
1114 Double_t lXiCosineOfPointingAngle = -1. ;
1115 Double_t lPosXi[3] = { -1000.0, -1000.0, -1000.0 };
1116 Double_t lXiRadius = -1000. ;
1118 // - 2nd part of initialisation : Nbr of clusters within TPC for the 3 daughter cascade tracks
1119 Int_t lPosTPCClusters = -1; // For ESD only ...//FIXME : wait for availability in AOD
1120 Int_t lNegTPCClusters = -1; // For ESD only ...
1121 Int_t lBachTPCClusters = -1; // For ESD only ...
1123 // - 3rd part of initialisation : about V0 part in cascades
1124 Double_t lInvMassLambdaAsCascDghter = 0.;
1125 //Double_t lV0Chi2Xi = -1. ;
1126 Double_t lDcaV0DaughtersXi = -1.;
1128 Double_t lDcaBachToPrimVertexXi = -1., lDcaV0ToPrimVertexXi = -1.;
1129 Double_t lDcaPosToPrimVertexXi = -1.;
1130 Double_t lDcaNegToPrimVertexXi = -1.;
1131 Double_t lV0CosineOfPointingAngleXi = -1. ;
1132 Double_t lV0CosineOfPointingAngleXiSpecial = -1. ;
1133 Double_t lPosV0Xi[3] = { -1000. , -1000., -1000. }; // Position of VO coming from cascade
1134 Double_t lV0RadiusXi = -1000.0;
1135 Double_t lV0quality = 0.;
1137 // - 4th part of initialisation : Effective masses
1138 Double_t lInvMassXiMinus = 0.;
1139 Double_t lInvMassXiPlus = 0.;
1140 Double_t lInvMassOmegaMinus = 0.;
1141 Double_t lInvMassOmegaPlus = 0.;
1143 // - 6th part of initialisation : extra info for QA
1144 Double_t lXiMomX = 0. , lXiMomY = 0., lXiMomZ = 0.;
1145 Double_t lXiTransvMom = 0. ;
1146 //Double_t lXiTransvMomMC= 0. ;
1147 Double_t lXiTotMom = 0. ;
1149 Double_t lBachMomX = 0., lBachMomY = 0., lBachMomZ = 0.;
1150 //Double_t lBachTransvMom = 0.;
1151 //Double_t lBachTotMom = 0.;
1153 fTreeCascVarNegNSigmaPion = -100;
1154 fTreeCascVarNegNSigmaProton = -100;
1155 fTreeCascVarPosNSigmaPion = -100;
1156 fTreeCascVarPosNSigmaProton = -100;
1157 fTreeCascVarBachNSigmaPion = -100;
1158 fTreeCascVarBachNSigmaKaon = -100;
1160 Short_t lChargeXi = -2;
1161 //Double_t lV0toXiCosineOfPointingAngle = 0. ;
1163 Double_t lRapXi = -20.0, lRapOmega = -20.0; // lEta = -20.0, lTheta = 360., lPhi = 720. ;
1164 //Double_t lAlphaXi = -200., lPtArmXi = -200.0;
1166 // -------------------------------------
1167 // II.ESD - Calculation Part dedicated to Xi vertices (ESD)
1169 AliESDcascade *xi = lESDevent->GetCascade(iXi);
1172 // - II.Step 2 : Assigning the necessary variables for specific AliESDcascade data members (ESD)
1175 xi->ChangeMassHypothesis(lV0quality , 3312); // default working hypothesis : cascade = Xi- decay
1177 lEffMassXi = xi->GetEffMassXi();
1178 //lChi2Xi = xi->GetChi2Xi();
1179 lDcaXiDaughters = xi->GetDcaXiDaughters();
1180 lXiCosineOfPointingAngle = xi->GetCascadeCosineOfPointingAngle( lBestPrimaryVtxPos[0],
1181 lBestPrimaryVtxPos[1],
1182 lBestPrimaryVtxPos[2] );
1183 // Take care : the best available vertex should be used (like in AliCascadeVertexer)
1185 xi->GetXYZcascade( lPosXi[0], lPosXi[1], lPosXi[2] );
1186 lXiRadius = TMath::Sqrt( lPosXi[0]*lPosXi[0] + lPosXi[1]*lPosXi[1] );
1188 // - II.Step 3 : around the tracks : Bach + V0 (ESD)
1189 // ~ Necessary variables for ESDcascade data members coming from the ESDv0 part (inheritance)
1192 UInt_t lIdxPosXi = (UInt_t) TMath::Abs( xi->GetPindex() );
1193 UInt_t lIdxNegXi = (UInt_t) TMath::Abs( xi->GetNindex() );
1194 UInt_t lBachIdx = (UInt_t) TMath::Abs( xi->GetBindex() );
1195 // Care track label can be negative in MC production (linked with the track quality)
1196 // However = normally, not the case for track index ...
1198 // FIXME : rejection of a double use of a daughter track (nothing but just a crosscheck of what is done in the cascade vertexer)
1199 if(lBachIdx == lIdxNegXi) {
1200 AliWarning("Pb / Idx(Bach. track) = Idx(Neg. track) ... continue!"); continue;
1202 if(lBachIdx == lIdxPosXi) {
1203 AliWarning("Pb / Idx(Bach. track) = Idx(Pos. track) ... continue!"); continue;
1206 AliESDtrack *pTrackXi = lESDevent->GetTrack( lIdxPosXi );
1207 AliESDtrack *nTrackXi = lESDevent->GetTrack( lIdxNegXi );
1208 AliESDtrack *bachTrackXi = lESDevent->GetTrack( lBachIdx );
1210 if (!pTrackXi || !nTrackXi || !bachTrackXi ) {
1211 AliWarning("ERROR: Could not retrieve one of the 3 ESD daughter tracks of the cascade ...");
1215 fTreeCascVarPosEta = pTrackXi->Eta();
1216 fTreeCascVarNegEta = nTrackXi->Eta();
1217 fTreeCascVarBachEta = bachTrackXi->Eta();
1219 Double_t lBMom[3], lNMom[3], lPMom[3];
1220 xi->GetBPxPyPz( lBMom[0], lBMom[1], lBMom[2] );
1221 xi->GetPPxPyPz( lPMom[0], lPMom[1], lPMom[2] );
1222 xi->GetNPxPyPz( lNMom[0], lNMom[1], lNMom[2] );
1224 //fTreeCascVarBachTransMom = TMath::Sqrt( lBMom[0]*lBMom[0] + lBMom[1]*lBMom[1] );
1225 //fTreeCascVarPosTransMom = TMath::Sqrt( lPMom[0]*lPMom[0] + lPMom[1]*lPMom[1] );
1226 //fTreeCascVarNegTransMom = TMath::Sqrt( lNMom[0]*lNMom[0] + lNMom[1]*lNMom[1] );
1228 //------------------------------------------------
1229 // TPC dEdx information
1230 //------------------------------------------------
1231 fTreeCascVarNegNSigmaPion = fPIDResponse->NumberOfSigmasTPC( nTrackXi, AliPID::kPion );
1232 fTreeCascVarNegNSigmaProton = fPIDResponse->NumberOfSigmasTPC( nTrackXi, AliPID::kProton );
1233 fTreeCascVarPosNSigmaPion = fPIDResponse->NumberOfSigmasTPC( pTrackXi, AliPID::kPion );
1234 fTreeCascVarPosNSigmaProton = fPIDResponse->NumberOfSigmasTPC( pTrackXi, AliPID::kProton );
1235 fTreeCascVarBachNSigmaPion = fPIDResponse->NumberOfSigmasTPC( bachTrackXi, AliPID::kPion );
1236 fTreeCascVarBachNSigmaKaon = fPIDResponse->NumberOfSigmasTPC( bachTrackXi, AliPID::kKaon );
1238 //------------------------------------------------
1239 // TPC Number of clusters info
1240 // --- modified to save the smallest number
1241 // --- of TPC clusters for the 3 tracks
1242 //------------------------------------------------
1244 lPosTPCClusters = pTrackXi->GetTPCNcls();
1245 lNegTPCClusters = nTrackXi->GetTPCNcls();
1246 lBachTPCClusters = bachTrackXi->GetTPCNcls();
1248 // 1 - Poor quality related to TPCrefit
1249 ULong_t pStatus = pTrackXi->GetStatus();
1250 ULong_t nStatus = nTrackXi->GetStatus();
1251 ULong_t bachStatus = bachTrackXi->GetStatus();
1253 //fTreeCascVarkITSRefitBachelor = kTRUE;
1254 //fTreeCascVarkITSRefitNegative = kTRUE;
1255 //fTreeCascVarkITSRefitPositive = kTRUE;
1257 if ((pStatus&AliESDtrack::kTPCrefit) == 0) { AliWarning("Pb / V0 Pos. track has no TPCrefit ... continue!"); continue; }
1258 if ((nStatus&AliESDtrack::kTPCrefit) == 0) { AliWarning("Pb / V0 Neg. track has no TPCrefit ... continue!"); continue; }
1259 if ((bachStatus&AliESDtrack::kTPCrefit) == 0) { AliWarning("Pb / Bach. track has no TPCrefit ... continue!"); continue; }
1261 // 2 - Poor quality related to TPC clusters: lowest cut of 70 clusters
1262 if(lPosTPCClusters < 70) { AliWarning("Pb / V0 Pos. track has less than 70 TPC clusters ... continue!"); continue; }
1263 if(lNegTPCClusters < 70) { AliWarning("Pb / V0 Neg. track has less than 70 TPC clusters ... continue!"); continue; }
1264 if(lBachTPCClusters < 70) { AliWarning("Pb / Bach. track has less than 70 TPC clusters ... continue!"); continue; }
1265 Int_t leastnumberofclusters = 1000;
1266 if( lPosTPCClusters < leastnumberofclusters ) leastnumberofclusters = lPosTPCClusters;
1267 if( lNegTPCClusters < leastnumberofclusters ) leastnumberofclusters = lNegTPCClusters;
1268 if( lBachTPCClusters < leastnumberofclusters ) leastnumberofclusters = lBachTPCClusters;
1270 lInvMassLambdaAsCascDghter = xi->GetEffMass();
1271 // This value shouldn't change, whatever the working hyp. is : Xi-, Xi+, Omega-, Omega+
1272 lDcaV0DaughtersXi = xi->GetDcaV0Daughters();
1273 //lV0Chi2Xi = xi->GetChi2V0();
1275 lV0CosineOfPointingAngleXi = xi->GetV0CosineOfPointingAngle( lBestPrimaryVtxPos[0],
1276 lBestPrimaryVtxPos[1],
1277 lBestPrimaryVtxPos[2] );
1278 //Modification: V0 CosPA wrt to Cascade decay vertex
1279 lV0CosineOfPointingAngleXiSpecial = xi->GetV0CosineOfPointingAngle( lPosXi[0],
1283 lDcaV0ToPrimVertexXi = xi->GetD( lBestPrimaryVtxPos[0],
1284 lBestPrimaryVtxPos[1],
1285 lBestPrimaryVtxPos[2] );
1287 lDcaBachToPrimVertexXi = TMath::Abs( bachTrackXi->GetD( lBestPrimaryVtxPos[0],
1288 lBestPrimaryVtxPos[1],
1290 // Note : AliExternalTrackParam::GetD returns an algebraic value ...
1292 xi->GetXYZ( lPosV0Xi[0], lPosV0Xi[1], lPosV0Xi[2] );
1293 lV0RadiusXi = TMath::Sqrt( lPosV0Xi[0]*lPosV0Xi[0] + lPosV0Xi[1]*lPosV0Xi[1] );
1295 lDcaPosToPrimVertexXi = TMath::Abs( pTrackXi ->GetD( lBestPrimaryVtxPos[0],
1296 lBestPrimaryVtxPos[1],
1299 lDcaNegToPrimVertexXi = TMath::Abs( nTrackXi ->GetD( lBestPrimaryVtxPos[0],
1300 lBestPrimaryVtxPos[1],
1303 // - II.Step 4 : around effective masses (ESD)
1304 // ~ change mass hypotheses to cover all the possibilities : Xi-/+, Omega -/+
1306 if( bachTrackXi->Charge() < 0 ) {
1308 xi->ChangeMassHypothesis(lV0quality , 3312);
1309 // Calculate the effective mass of the Xi- candidate.
1310 // pdg code 3312 = Xi-
1311 lInvMassXiMinus = xi->GetEffMassXi();
1314 xi->ChangeMassHypothesis(lV0quality , 3334);
1315 // Calculate the effective mass of the Xi- candidate.
1316 // pdg code 3334 = Omega-
1317 lInvMassOmegaMinus = xi->GetEffMassXi();
1320 xi->ChangeMassHypothesis(lV0quality , 3312); // Back to default hyp.
1321 }// end if negative bachelor
1324 if( bachTrackXi->Charge() > 0 ){
1326 xi->ChangeMassHypothesis(lV0quality , -3312);
1327 // Calculate the effective mass of the Xi+ candidate.
1328 // pdg code -3312 = Xi+
1329 lInvMassXiPlus = xi->GetEffMassXi();
1332 xi->ChangeMassHypothesis(lV0quality , -3334);
1333 // Calculate the effective mass of the Xi+ candidate.
1334 // pdg code -3334 = Omega+
1335 lInvMassOmegaPlus = xi->GetEffMassXi();
1338 xi->ChangeMassHypothesis(lV0quality , -3312); // Back to "default" hyp.
1339 }// end if positive bachelor
1340 // - II.Step 6 : extra info for QA (ESD)
1341 // miscellaneous pieces of info that may help regarding data quality assessment.
1344 xi->GetPxPyPz( lXiMomX, lXiMomY, lXiMomZ );
1345 lXiTransvMom = TMath::Sqrt( lXiMomX*lXiMomX + lXiMomY*lXiMomY );
1346 lXiTotMom = TMath::Sqrt( lXiMomX*lXiMomX + lXiMomY*lXiMomY + lXiMomZ*lXiMomZ );
1348 xi->GetBPxPyPz( lBachMomX, lBachMomY, lBachMomZ );
1349 //lBachTransvMom = TMath::Sqrt( lBachMomX*lBachMomX + lBachMomY*lBachMomY );
1350 //lBachTotMom = TMath::Sqrt( lBachMomX*lBachMomX + lBachMomY*lBachMomY + lBachMomZ*lBachMomZ );
1352 lChargeXi = xi->Charge();
1354 //lV0toXiCosineOfPointingAngle = xi->GetV0CosineOfPointingAngle( lPosXi[0], lPosXi[1], lPosXi[2] );
1356 lRapXi = xi->RapXi();
1357 lRapOmega = xi->RapOmega();
1359 //lTheta = xi->Theta() *180.0/TMath::Pi();
1360 //lPhi = xi->Phi() *180.0/TMath::Pi();
1361 //lAlphaXi = xi->AlphaXi();
1362 //lPtArmXi = xi->PtArmXi();
1364 //------------------------------------------------
1365 // Set Variables for adding to tree
1366 //------------------------------------------------
1368 fTreeCascVarCharge = lChargeXi;
1369 if(lInvMassXiMinus!=0) fTreeCascVarMassAsXi = lInvMassXiMinus;
1370 if(lInvMassXiPlus!=0) fTreeCascVarMassAsXi = lInvMassXiPlus;
1371 if(lInvMassOmegaMinus!=0) fTreeCascVarMassAsOmega = lInvMassOmegaMinus;
1372 if(lInvMassOmegaPlus!=0) fTreeCascVarMassAsOmega = lInvMassOmegaPlus;
1373 fTreeCascVarPt = lXiTransvMom;
1374 fTreeCascVarRapXi = lRapXi ;
1375 fTreeCascVarRapOmega = lRapOmega ;
1376 fTreeCascVarDCACascDaughters = lDcaXiDaughters;
1377 fTreeCascVarDCABachToPrimVtx = lDcaBachToPrimVertexXi;
1378 fTreeCascVarDCAV0Daughters = lDcaV0DaughtersXi;
1379 fTreeCascVarDCAV0ToPrimVtx = lDcaV0ToPrimVertexXi;
1380 fTreeCascVarDCAPosToPrimVtx = lDcaPosToPrimVertexXi;
1381 fTreeCascVarDCANegToPrimVtx = lDcaNegToPrimVertexXi;
1382 fTreeCascVarCascCosPointingAngle = lXiCosineOfPointingAngle;
1383 fTreeCascVarCascRadius = lXiRadius;
1384 fTreeCascVarV0Mass = lInvMassLambdaAsCascDghter;
1385 fTreeCascVarV0CosPointingAngle = lV0CosineOfPointingAngleXi;
1386 fTreeCascVarV0CosPointingAngleSpecial = lV0CosineOfPointingAngleXiSpecial;
1387 fTreeCascVarV0Radius = lV0RadiusXi;
1388 fTreeCascVarLeastNbrClusters = leastnumberofclusters;
1390 //Copy Multiplicity information
1391 fTreeCascVarCentV0A = fCentrality_V0A;
1392 fTreeCascVarCentV0C = fCentrality_V0C;
1393 fTreeCascVarCentV0M = fCentrality_V0M;
1394 fTreeCascVarCentV0AEq = fCentrality_V0AEq;
1395 fTreeCascVarCentV0CEq = fCentrality_V0CEq;
1396 fTreeCascVarCentV0MEq = fCentrality_V0MEq;
1397 fTreeCascVarAmpV0A = fAmplitude_V0A;
1398 fTreeCascVarAmpV0C = fAmplitude_V0C;
1399 fTreeCascVarAmpV0AEq = fAmplitude_V0AEq;
1400 fTreeCascVarAmpV0CEq = fAmplitude_V0CEq;
1401 fTreeCascVarRefMultEta8 = fRefMultEta8;
1402 fTreeCascVarRefMultEta5 = fRefMultEta5;
1403 fTreeCascVarRunNumber = fRunNumber;
1405 fTreeCascVarDistOverTotMom = TMath::Sqrt(
1406 TMath::Power( lPosXi[0] - lBestPrimaryVtxPos[0] , 2) +
1407 TMath::Power( lPosXi[1] - lBestPrimaryVtxPos[1] , 2) +
1408 TMath::Power( lPosXi[2] - lBestPrimaryVtxPos[2] , 2)
1410 fTreeCascVarDistOverTotMom /= (lXiTotMom+1e-13);
1412 //All vars not specified here: specified elsewhere!
1414 //------------------------------------------------
1416 //------------------------------------------------
1418 // The conditional is meant to decrease excessive
1419 // memory usage! Be careful when loosening the
1422 //Xi Mass window: 150MeV wide
1423 //Omega mass window: 150MeV wide
1425 if( fkSaveCascadeTree && ( (fTreeCascVarMassAsXi<1.32+0.075&&fTreeCascVarMassAsXi>1.32-0.075) ||
1426 (fTreeCascVarMassAsOmega<1.68+0.075&&fTreeCascVarMassAsOmega>1.68-0.075) ) ){
1427 fTreeCascade->Fill();
1430 //------------------------------------------------
1432 //------------------------------------------------
1434 }// end of the Cascade loop (ESD or AOD)
1436 // Post output data.
1437 PostData(1, fListHist);
1438 PostData(2, fTreeEvent);
1439 PostData(3, fTreeV0);
1440 PostData(4, fTreeCascade);
1443 //________________________________________________________________________
1444 void AliAnalysisTaskStrangenessVsMultiplicity::Terminate(Option_t *)
1446 // Draw result to the screen
1447 // Called once at the end of the query
1449 TList *cRetrievedList = 0x0;
1450 cRetrievedList = (TList*)GetOutputData(1);
1451 if(!cRetrievedList){
1452 Printf("ERROR - AliAnalysisTaskStrangenessVsMultiplicity : ouput data container list not available\n");
1456 fHistEventCounter = dynamic_cast<TH1D*> ( cRetrievedList->FindObject("fHistEventCounter") );
1457 if (!fHistEventCounter) {
1458 Printf("ERROR - AliAnalysisTaskStrangenessVsMultiplicity : fHistEventCounter not available");
1462 TCanvas *canCheck = new TCanvas("AliAnalysisTaskStrangenessVsMultiplicity","V0 Multiplicity",10,10,510,510);
1463 canCheck->cd(1)->SetLogy();
1465 fHistEventCounter->SetMarkerStyle(22);
1466 fHistEventCounter->DrawCopy("E");
1469 //----------------------------------------------------------------------------
1471 Double_t AliAnalysisTaskStrangenessVsMultiplicity::MyRapidity(Double_t rE, Double_t rPz) const
1473 // Local calculation for rapidity
1474 Double_t ReturnValue = -100;
1475 if( (rE-rPz+1.e-13) != 0 && (rE+rPz) != 0 ){
1476 ReturnValue = 0.5*TMath::Log((rE+rPz)/(rE-rPz+1.e-13));