]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/STRANGENESS/Cascades/AliAnalysisTaskStrangenessVsMultiplicity.cxx
New MFT analysis tools
[u/mrichter/AliRoot.git] / PWGLF / STRANGENESS / Cascades / AliAnalysisTaskStrangenessVsMultiplicity.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 // +-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
17 //
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:
21 //
22 //  1) (under development) it should act as an auxiliary task and provide a
23 //     calibrated estimator
24 //
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:
28 //
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)
37 //      --- Run Number
38 //
39 // +-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
40
41 class TTree;
42 class TParticle;
43 class TVector3;
44
45 //class AliMCEventHandler;
46 //class AliMCEvent;
47 //class AliStack;
48
49 class AliESDVertex;
50 class AliAODVertex;
51 class AliESDv0;
52 class AliAODv0;
53
54 #include <Riostream.h>
55 #include "TList.h"
56 #include "TH1.h"
57 #include "TH2.h"
58 #include "TH3.h"
59 #include "TFile.h"
60 #include "THnSparse.h"
61 #include "TVector3.h"
62 #include "TCanvas.h"
63 #include "TMath.h"
64 #include "TLegend.h"
65 //#include "AliLog.h"
66
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"
78 #include "AliStack.h"
79 #include "AliCentrality.h"
80 #include "AliPPVsMultUtils.h"
81
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"
92
93 using std::cout;
94 using std::endl;
95
96 ClassImp(AliAnalysisTaskStrangenessVsMultiplicity)
97
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 fkApplyTrackletsVsClustersCut(kTRUE),
105 //---> Variables for fTreeEvent
106 fAmplitude_V0A   (0),
107 fAmplitude_V0C   (0),
108 fAmplitude_V0AEq (0),
109 fAmplitude_V0CEq (0),
110 fCentrality_V0A(0),
111 fCentrality_V0C(0),
112 fCentrality_V0M(0),
113 fCentrality_V0AEq(0),
114 fCentrality_V0CEq(0),
115 fCentrality_V0MEq(0),
116 fCustomCentrality_V0M(0),
117 fCustomCentrality_V0MEq(0),
118 fRefMultEta5(0),
119 fRefMultEta8(0),
120 fRunNumber(0),
121 fEvSel_HasAtLeastSPDVertex(0),
122 fEvSel_VtxZCut(0),
123 fEvSel_IsNotPileup(0),
124 fEvSel_IsNotPileupMV(0),
125 fEvSel_IsNotPileupInMultBins(0),
126 fEvSel_Triggered(0),
127 fEvSel_nTracklets(0),
128 fEvSel_nSPDClusters(0),
129 fEvSel_VtxZ(0),
130 fEvSel_nSPDPrimVertices(0),
131 fEvSel_distZ(0),
132 fEvSel_nContributors(0),
133 fEvSel_nContributorsPileup(0),
134
135 //---> Variables for fTreeV0
136 fTreeVariableChi2V0(0),
137 fTreeVariableDcaV0Daughters(0),
138 fTreeVariableDcaV0ToPrimVertex(0),
139 fTreeVariableDcaPosToPrimVertex(0),
140 fTreeVariableDcaNegToPrimVertex(0),
141 fTreeVariableV0CosineOfPointingAngle(0),
142 fTreeVariableV0Radius(0),
143 fTreeVariablePt(0),
144 fTreeVariableRapK0Short(0),
145 fTreeVariableRapLambda(0),
146 fTreeVariableInvMassK0s(0),
147 fTreeVariableInvMassLambda(0),
148 fTreeVariableInvMassAntiLambda(0),
149 fTreeVariableAlphaV0(0),
150 fTreeVariablePtArmV0(0),
151 fTreeVariableNegEta(0),
152 fTreeVariablePosEta(0),
153
154 fTreeVariableNSigmasPosProton(0),
155 fTreeVariableNSigmasPosPion(0),
156 fTreeVariableNSigmasNegProton(0),
157 fTreeVariableNSigmasNegPion(0),
158
159 fTreeVariableDistOverTotMom(0),
160 fTreeVariableLeastNbrCrossedRows(0),
161 fTreeVariableLeastRatioCrossedRowsOverFindable(0),
162
163 fTreeVariableCentV0A(0),
164 fTreeVariableCentV0C(0),
165 fTreeVariableCentV0M(0),
166 fTreeVariableCentV0AEq(0),
167 fTreeVariableCentV0CEq(0),
168 fTreeVariableCentV0MEq(0),
169 fTreeVariableCustomCentV0M(0),
170 fTreeVariableAmpV0A(0),
171 fTreeVariableAmpV0C(0),
172 fTreeVariableAmpV0AEq(0),
173 fTreeVariableAmpV0CEq(0),
174 fTreeVariableRefMultEta8(0),
175 fTreeVariableRefMultEta5(0),
176 fTreeVariableRunNumber(0),
177 //---> Variables for fTreeCascade
178 fTreeCascVarCharge(0),
179 fTreeCascVarMassAsXi(0),
180 fTreeCascVarMassAsOmega(0),
181 fTreeCascVarPt(0),
182 fTreeCascVarRapXi(0),
183 fTreeCascVarRapOmega(0),
184 fTreeCascVarNegEta(0),
185 fTreeCascVarPosEta(0),
186 fTreeCascVarBachEta(0),
187 fTreeCascVarDCACascDaughters(0),
188 fTreeCascVarDCABachToPrimVtx(0),
189 fTreeCascVarDCAV0Daughters(0),
190 fTreeCascVarDCAV0ToPrimVtx(0),
191 fTreeCascVarDCAPosToPrimVtx(0),
192 fTreeCascVarDCANegToPrimVtx(0),
193 fTreeCascVarCascCosPointingAngle(0),
194 fTreeCascVarCascRadius(0),
195 fTreeCascVarV0Mass(0),
196 fTreeCascVarV0CosPointingAngle(0),
197 fTreeCascVarV0CosPointingAngleSpecial(0),
198 fTreeCascVarV0Radius(0),
199 fTreeCascVarLeastNbrClusters(0),
200 fTreeCascVarDistOverTotMom(0),
201 fTreeCascVarNegNSigmaPion(0),
202 fTreeCascVarNegNSigmaProton(0),
203 fTreeCascVarPosNSigmaPion(0),
204 fTreeCascVarPosNSigmaProton(0),
205 fTreeCascVarBachNSigmaPion(0),
206 fTreeCascVarBachNSigmaKaon(0),
207 fTreeCascVarCentV0A(0),
208 fTreeCascVarCentV0C(0),
209 fTreeCascVarCentV0M(0),
210 fTreeCascVarCentV0AEq(0),
211 fTreeCascVarCentV0CEq(0),
212 fTreeCascVarCentV0MEq(0),
213 fTreeCascVarCustomCentV0M(0),
214 fTreeCascVarAmpV0A(0),
215 fTreeCascVarAmpV0C(0),
216 fTreeCascVarAmpV0AEq(0),
217 fTreeCascVarAmpV0CEq(0),
218 fTreeCascVarRefMultEta8(0),
219 fTreeCascVarRefMultEta5(0),
220 fTreeCascVarRunNumber(0),
221 //Histos
222 fHistEventCounter(0)
223 //------------------------------------------------
224 // Tree Variables
225 {
226     
227 }
228
229 AliAnalysisTaskStrangenessVsMultiplicity::AliAnalysisTaskStrangenessVsMultiplicity(const char *name)
230 : AliAnalysisTaskSE(name), fListHist(0), fTreeEvent(0), fTreeV0(0), fTreeCascade(0), fPIDResponse(0), fESDtrackCuts(0), fPPVsMultUtils(0), fUtils(0),
231 fkSaveV0Tree      ( kFALSE ),
232 fkSaveCascadeTree ( kTRUE  ),
233 fkRunVertexers    ( kTRUE  ),
234 fkSkipEventSelection( kFALSE ),
235 fkApplyTrackletsVsClustersCut(kTRUE),
236 //---> Variables for fTreeEvent
237 fAmplitude_V0A (0),
238 fAmplitude_V0C (0),
239 fAmplitude_V0AEq (0),
240 fAmplitude_V0CEq (0),
241 fCentrality_V0A(0),
242 fCentrality_V0C(0),
243 fCentrality_V0M(0),
244 fCentrality_V0AEq(0),
245 fCentrality_V0CEq(0),
246 fCentrality_V0MEq(0),
247 fCustomCentrality_V0M(0),
248 fCustomCentrality_V0MEq(0),
249 fRefMultEta5(0),
250 fRefMultEta8(0),
251 fRunNumber(0),
252 fEvSel_HasAtLeastSPDVertex(0),
253 fEvSel_VtxZCut(0),
254 fEvSel_IsNotPileup(0),
255 fEvSel_IsNotPileupMV(0),
256 fEvSel_IsNotPileupInMultBins(0),
257 fEvSel_Triggered(0),
258 fEvSel_nTracklets(0),
259 fEvSel_nSPDClusters(0),
260 fEvSel_VtxZ(0),
261 fEvSel_nSPDPrimVertices(0),
262 fEvSel_distZ(0),
263 fEvSel_nContributors(0),
264 fEvSel_nContributorsPileup(0),
265
266 //---> Variables for fTreeV0
267 fTreeVariableChi2V0(0),
268 fTreeVariableDcaV0Daughters(0),
269 fTreeVariableDcaV0ToPrimVertex(0),
270 fTreeVariableDcaPosToPrimVertex(0),
271 fTreeVariableDcaNegToPrimVertex(0),
272 fTreeVariableV0CosineOfPointingAngle(0),
273 fTreeVariableV0Radius(0),
274 fTreeVariablePt(0),
275 fTreeVariableRapK0Short(0),
276 fTreeVariableRapLambda(0),
277 fTreeVariableInvMassK0s(0),
278 fTreeVariableInvMassLambda(0),
279 fTreeVariableInvMassAntiLambda(0),
280 fTreeVariableAlphaV0(0),
281 fTreeVariablePtArmV0(0),
282 fTreeVariableNegEta(0),
283 fTreeVariablePosEta(0),
284
285 fTreeVariableNSigmasPosProton(0),
286 fTreeVariableNSigmasPosPion(0),
287 fTreeVariableNSigmasNegProton(0),
288 fTreeVariableNSigmasNegPion(0),
289
290 fTreeVariableDistOverTotMom(0),
291 fTreeVariableLeastNbrCrossedRows(0),
292 fTreeVariableLeastRatioCrossedRowsOverFindable(0),
293
294 fTreeVariableCentV0A(0),
295 fTreeVariableCentV0C(0),
296 fTreeVariableCentV0M(0),
297 fTreeVariableCentV0AEq(0),
298 fTreeVariableCentV0CEq(0),
299 fTreeVariableCentV0MEq(0),
300 fTreeVariableCustomCentV0M(0),
301 fTreeVariableAmpV0A(0),
302 fTreeVariableAmpV0C(0),
303 fTreeVariableAmpV0AEq(0),
304 fTreeVariableAmpV0CEq(0),
305 fTreeVariableRefMultEta8(0),
306 fTreeVariableRefMultEta5(0),
307 fTreeVariableRunNumber(0),
308 //---> Variables for fTreeCascade
309 fTreeCascVarCharge(0),
310 fTreeCascVarMassAsXi(0),
311 fTreeCascVarMassAsOmega(0),
312 fTreeCascVarPt(0),
313 fTreeCascVarRapXi(0),
314 fTreeCascVarRapOmega(0),
315 fTreeCascVarNegEta(0),
316 fTreeCascVarPosEta(0),
317 fTreeCascVarBachEta(0),
318 fTreeCascVarDCACascDaughters(0),
319 fTreeCascVarDCABachToPrimVtx(0),
320 fTreeCascVarDCAV0Daughters(0),
321 fTreeCascVarDCAV0ToPrimVtx(0),
322 fTreeCascVarDCAPosToPrimVtx(0),
323 fTreeCascVarDCANegToPrimVtx(0),
324 fTreeCascVarCascCosPointingAngle(0),
325 fTreeCascVarCascRadius(0),
326 fTreeCascVarV0Mass(0),
327 fTreeCascVarV0CosPointingAngle(0),
328 fTreeCascVarV0CosPointingAngleSpecial(0),
329 fTreeCascVarV0Radius(0),
330 fTreeCascVarLeastNbrClusters(0),
331 fTreeCascVarDistOverTotMom(0),
332 fTreeCascVarNegNSigmaPion(0),
333 fTreeCascVarNegNSigmaProton(0),
334 fTreeCascVarPosNSigmaPion(0),
335 fTreeCascVarPosNSigmaProton(0),
336 fTreeCascVarBachNSigmaPion(0),
337 fTreeCascVarBachNSigmaKaon(0),
338 fTreeCascVarCentV0A(0),
339 fTreeCascVarCentV0C(0),
340 fTreeCascVarCentV0M(0),
341 fTreeCascVarCentV0AEq(0),
342 fTreeCascVarCentV0CEq(0),
343 fTreeCascVarCentV0MEq(0),
344 fTreeCascVarCustomCentV0M(0),
345 fTreeCascVarAmpV0A(0),
346 fTreeCascVarAmpV0C(0),
347 fTreeCascVarAmpV0AEq(0),
348 fTreeCascVarAmpV0CEq(0),
349 fTreeCascVarRefMultEta8(0),
350 fTreeCascVarRefMultEta5(0),
351 fTreeCascVarRunNumber(0),
352 //Histos
353 fHistEventCounter(0)
354 {
355     
356     //Re-vertex: Will only apply for cascade candidates
357     
358     fV0VertexerSels[0] =  33.  ;  // max allowed chi2
359     fV0VertexerSels[1] =   0.02;  // min allowed impact parameter for the 1st daughter (LHC09a4 : 0.05)
360     fV0VertexerSels[2] =   0.02;  // min allowed impact parameter for the 2nd daughter (LHC09a4 : 0.05)
361     fV0VertexerSels[3] =   2.0 ;  // max allowed DCA between the daughter tracks       (LHC09a4 : 0.5)
362     fV0VertexerSels[4] =   0.95;  // min allowed cosine of V0's pointing angle         (LHC09a4 : 0.99)
363     fV0VertexerSels[5] =   1.0 ;  // min radius of the fiducial volume                 (LHC09a4 : 0.2)
364     fV0VertexerSels[6] = 200.  ;  // max radius of the fiducial volume                 (LHC09a4 : 100.0)
365     
366     fCascadeVertexerSels[0] =  33.   ;  // max allowed chi2 (same as PDC07)
367     fCascadeVertexerSels[1] =   0.05 ;  // min allowed V0 impact parameter                    (PDC07 : 0.05   / LHC09a4 : 0.025 )
368     fCascadeVertexerSels[2] =   0.010;  // "window" around the Lambda mass                    (PDC07 : 0.008  / LHC09a4 : 0.010 )
369     fCascadeVertexerSels[3] =   0.03 ;  // min allowed bachelor's impact parameter            (PDC07 : 0.035  / LHC09a4 : 0.025 )
370     fCascadeVertexerSels[4] =   2.0  ;  // max allowed DCA between the V0 and the bachelor    (PDC07 : 0.1    / LHC09a4 : 0.2   )
371     fCascadeVertexerSels[5] =   0.95 ;  // min allowed cosine of the cascade pointing angle   (PDC07 : 0.9985 / LHC09a4 : 0.998 )
372     fCascadeVertexerSels[6] =   0.4  ;  // min radius of the fiducial volume                  (PDC07 : 0.9    / LHC09a4 : 0.2   )
373     fCascadeVertexerSels[7] = 100.   ;  // max radius of the fiducial volume                  (PDC07 : 100    / LHC09a4 : 100   )
374     
375     
376     DefineOutput(1, TList::Class()); // Event Counter Histo
377     DefineOutput(2, TTree::Class()); // Event Tree
378     DefineOutput(3, TTree::Class()); // V0 Tree
379     DefineOutput(4, TTree::Class()); // Cascade Tree
380 }
381
382
383 AliAnalysisTaskStrangenessVsMultiplicity::~AliAnalysisTaskStrangenessVsMultiplicity()
384 {
385     //------------------------------------------------
386     // DESTRUCTOR
387     //------------------------------------------------
388     
389     if (fListHist){
390         delete fListHist;
391         fListHist = 0x0;
392     }
393     if (fTreeEvent){
394         delete fTreeEvent;
395         fTreeEvent = 0x0;
396     }
397     if (fTreeV0){
398         delete fTreeV0;
399         fTreeV0 = 0x0;
400     }
401     if (fTreeCascade){
402         delete fTreeCascade;
403         fTreeCascade = 0x0;
404     }
405     if (fPPVsMultUtils){
406         delete fPPVsMultUtils;
407         fPPVsMultUtils = 0x0;
408     }
409     if (fUtils){
410         delete fUtils;
411         fUtils = 0x0;
412     }
413 }
414
415 //________________________________________________________________________
416 void AliAnalysisTaskStrangenessVsMultiplicity::UserCreateOutputObjects()
417 {
418     
419     OpenFile(2);
420     // Called once
421     
422     //------------------------------------------------
423     
424     fTreeEvent = new TTree("fTreeEvent","Event");
425     
426     //------------------------------------------------
427     // fTree Branch definitions - Event by Event info
428     //------------------------------------------------
429     
430     //-----------BASIC-INFO---------------------------
431     
432     //--- VZERO Data (Integrated)
433     fTreeEvent->Branch("fAmplitude_V0A",&fAmplitude_V0A,"fAmplitude_V0A/F");
434     fTreeEvent->Branch("fAmplitude_V0C",&fAmplitude_V0C,"fAmplitude_V0C/F");
435     fTreeEvent->Branch("fAmplitude_V0AEq",&fAmplitude_V0AEq,"fAmplitude_V0AEq/F");
436     fTreeEvent->Branch("fAmplitude_V0CEq",&fAmplitude_V0CEq,"fAmplitude_V0CEq/F");
437     
438     //Info from AliCentrality (not necessarily 'centrality' per se)
439     fTreeEvent->Branch("fCentrality_V0A",&fCentrality_V0A,"fCentrality_V0A/F");
440     fTreeEvent->Branch("fCentrality_V0C",&fCentrality_V0C,"fCentrality_V0C/F");
441     fTreeEvent->Branch("fCentrality_V0M",&fCentrality_V0M,"fCentrality_V0M/F");
442     fTreeEvent->Branch("fCentrality_V0AEq",&fCentrality_V0AEq,"fCentrality_V0AEq/F");
443     fTreeEvent->Branch("fCentrality_V0CEq",&fCentrality_V0CEq,"fCentrality_V0CEq/F");
444     fTreeEvent->Branch("fCentrality_V0MEq",&fCentrality_V0MEq,"fCentrality_V0MEq/F");
445     fTreeEvent->Branch("fCustomCentrality_V0M",&fCustomCentrality_V0M,"fCustomCentrality_V0M/F");
446     fTreeEvent->Branch("fCustomCentrality_V0MEq",&fCustomCentrality_V0MEq,"fCustomCentrality_V0MEq/F");
447     
448     //Official GetReferenceMultiplicity
449     fTreeEvent->Branch("fRefMultEta5",&fRefMultEta5,"fRefMultEta5/I");
450     fTreeEvent->Branch("fRefMultEta8",&fRefMultEta8,"fRefMultEta8/I");
451     
452     //Run Number
453     fTreeEvent->Branch("fRunNumber", &fRunNumber, "fRunNumber/I");
454     
455     //Booleans for Event Selection
456     fTreeEvent->Branch("fEvSel_HasAtLeastSPDVertex", &fEvSel_HasAtLeastSPDVertex, "fEvSel_HasAtLeastSPDVertex/O");
457     fTreeEvent->Branch("fEvSel_VtxZCut", &fEvSel_VtxZCut, "fEvSel_VtxZCut/O");
458     fTreeEvent->Branch("fEvSel_IsNotPileup", &fEvSel_IsNotPileup, "fEvSel_IsNotPileup/O");
459     fTreeEvent->Branch("fEvSel_IsNotPileupMV", &fEvSel_IsNotPileupMV, "fEvSel_IsNotPileupMV/O");
460     fTreeEvent->Branch("fEvSel_IsNotPileupInMultBins", &fEvSel_IsNotPileupInMultBins, "fEvSel_IsNotPileupInMultBins/O");
461     fTreeEvent->Branch("fEvSel_Triggered", &fEvSel_Triggered, "fEvSel_Triggered/O");
462     
463     //Tracklets vs clusters test
464     fTreeEvent->Branch("fEvSel_nTracklets", &fEvSel_nTracklets, "fEvSel_nTracklets/I");
465     fTreeEvent->Branch("fEvSel_nSPDClusters", &fEvSel_nSPDClusters, "fEvSel_nSPDClusters/I");
466     
467     fTreeEvent->Branch("fEvSel_VtxZ", &fEvSel_VtxZ, "fEvSel_VtxZ/F");
468     fTreeEvent->Branch("fEvSel_nSPDPrimVertices", &fEvSel_nSPDPrimVertices, "fEvSel_nSPDPrimVertices/I");
469     fTreeEvent->Branch("fEvSel_distZ", &fEvSel_distZ, "fEvSel_distZ/F");
470     fTreeEvent->Branch("fEvSel_nContributors", &fEvSel_nContributors, "fEvSel_nContributors/I");
471     fTreeEvent->Branch("fEvSel_nContributorsPileup", &fEvSel_nContributorsPileup, "fEvSel_nContributorsPileup/I");
472     
473     //Create Basic V0 Output Tree
474     fTreeV0 = new TTree( "fTreeV0", "V0 Candidates");
475     
476     //------------------------------------------------
477     // fTreeV0 Branch definitions
478     //------------------------------------------------
479     
480     //-----------BASIC-INFO---------------------------
481     fTreeV0->Branch("fTreeVariableChi2V0",&fTreeVariableChi2V0,"fTreeVariableChi2V0/F");
482     fTreeV0->Branch("fTreeVariableDcaV0Daughters",&fTreeVariableDcaV0Daughters,"fTreeVariableDcaV0Daughters/F");
483     fTreeV0->Branch("fTreeVariableDcaV0ToPrimVertex",&fTreeVariableDcaV0ToPrimVertex,"fTreeVariableDcaV0ToPrimVertex/F");
484     fTreeV0->Branch("fTreeVariableDcaPosToPrimVertex",&fTreeVariableDcaPosToPrimVertex,"fTreeVariableDcaPosToPrimVertex/F");
485     fTreeV0->Branch("fTreeVariableDcaNegToPrimVertex",&fTreeVariableDcaNegToPrimVertex,"fTreeVariableDcaNegToPrimVertex/F");
486     fTreeV0->Branch("fTreeVariableV0Radius",&fTreeVariableV0Radius,"fTreeVariableV0Radius/F");
487     fTreeV0->Branch("fTreeVariablePt",&fTreeVariablePt,"fTreeVariablePt/F");
488     fTreeV0->Branch("fTreeVariableRapK0Short",&fTreeVariableRapK0Short,"fTreeVariableRapK0Short/F");
489     fTreeV0->Branch("fTreeVariableRapLambda",&fTreeVariableRapLambda,"fTreeVariableRapLambda/F");
490     fTreeV0->Branch("fTreeVariableInvMassK0s",&fTreeVariableInvMassK0s,"fTreeVariableInvMassK0s/F");
491     fTreeV0->Branch("fTreeVariableInvMassLambda",&fTreeVariableInvMassLambda,"fTreeVariableInvMassLambda/F");
492     fTreeV0->Branch("fTreeVariableInvMassAntiLambda",&fTreeVariableInvMassAntiLambda,"fTreeVariableInvMassAntiLambda/F");
493     fTreeV0->Branch("fTreeVariableV0CosineOfPointingAngle",&fTreeVariableV0CosineOfPointingAngle,"fTreeVariableV0CosineOfPointingAngle/F");
494     fTreeV0->Branch("fTreeVariableAlphaV0",&fTreeVariableAlphaV0,"fTreeVariableAlphaV0/F");
495     fTreeV0->Branch("fTreeVariablePtArmV0",&fTreeVariablePtArmV0,"fTreeVariablePtArmV0/F");
496     fTreeV0->Branch("fTreeVariableLeastNbrCrossedRows",&fTreeVariableLeastNbrCrossedRows,"fTreeVariableLeastNbrCrossedRows/I");
497     fTreeV0->Branch("fTreeVariableLeastRatioCrossedRowsOverFindable",&fTreeVariableLeastRatioCrossedRowsOverFindable,"fTreeVariableLeastRatioCrossedRowsOverFindable/F");
498     fTreeV0->Branch("fTreeVariableDistOverTotMom",&fTreeVariableDistOverTotMom,"fTreeVariableDistOverTotMom/F");
499     fTreeV0->Branch("fTreeVariableNSigmasPosProton",&fTreeVariableNSigmasPosProton,"fTreeVariableNSigmasPosProton/F");
500     fTreeV0->Branch("fTreeVariableNSigmasPosPion",&fTreeVariableNSigmasPosPion,"fTreeVariableNSigmasPosPion/F");
501     fTreeV0->Branch("fTreeVariableNSigmasNegProton",&fTreeVariableNSigmasNegProton,"fTreeVariableNSigmasNegProton/F");
502     fTreeV0->Branch("fTreeVariableNSigmasNegPion",&fTreeVariableNSigmasNegPion,"fTreeVariableNSigmasNegPion/F");
503     fTreeV0->Branch("fTreeVariableNegEta",&fTreeVariableNegEta,"fTreeVariableNegEta/F");
504     fTreeV0->Branch("fTreeVariablePosEta",&fTreeVariablePosEta,"fTreeVariablePosEta/F");
505     //-----------MULTIPLICITY-INFO--------------------
506     fTreeV0->Branch("fTreeVariableCentV0A",&fTreeVariableCentV0A,"fTreeVariableCentV0A/F");
507     fTreeV0->Branch("fTreeVariableCentV0C",&fTreeVariableCentV0C,"fTreeVariableCentV0C/F");
508     fTreeV0->Branch("fTreeVariableCentV0M",&fTreeVariableCentV0M,"fTreeVariableCentV0M/F");
509     fTreeV0->Branch("fTreeVariableCentV0AEq",&fTreeVariableCentV0AEq,"fTreeVariableCentV0AEq/F");
510     fTreeV0->Branch("fTreeVariableCentV0CEq",&fTreeVariableCentV0CEq,"fTreeVariableCentV0CEq/F");
511     fTreeV0->Branch("fTreeVariableCentV0MEq",&fTreeVariableCentV0MEq,"fTreeVariableCentV0MEq/F");
512     fTreeV0->Branch("fTreeVariableCustomCentV0M",&fTreeVariableCustomCentV0M,"fTreeVariableCustomCentV0M/F");
513     fTreeV0->Branch("fTreeVariableAmpV0A",&fTreeVariableAmpV0A,"fTreeVariableAmpV0A/F");
514     fTreeV0->Branch("fTreeVariableAmpV0C",&fTreeVariableAmpV0C,"fTreeVariableAmpV0C/F");
515     fTreeV0->Branch("fTreeVariableAmpV0AEq",&fTreeVariableAmpV0AEq,"fTreeVariableAmpV0AEq/F");
516     fTreeV0->Branch("fTreeVariableAmpV0CEq",&fTreeVariableAmpV0CEq,"fTreeVariableAmpV0CEq/F");
517     fTreeV0->Branch("fTreeVariableRefMultEta8",&fTreeVariableRefMultEta8,"fTreeVariableRefMultEta8/I");
518     fTreeV0->Branch("fTreeVariableRefMultEta5",&fTreeVariableRefMultEta5,"fTreeVariableRefMultEta5/I");
519     fTreeV0->Branch("fTreeVariableRunNumber",&fTreeVariableRunNumber,"fTreeVariableRunNumber/I");
520     //------------------------------------------------
521     
522     //Create Cascade output tree
523     fTreeCascade = new TTree("fTreeCascade","CascadeCandidates");
524     
525     //------------------------------------------------
526     // fTreeCascade Branch definitions - Cascade Tree
527     //------------------------------------------------
528     
529     //-----------BASIC-INFO---------------------------
530         fTreeCascade->Branch("fTreeCascVarCharge",&fTreeCascVarCharge,"fTreeCascVarCharge/I");
531     fTreeCascade->Branch("fTreeCascVarMassAsXi",&fTreeCascVarMassAsXi,"fTreeCascVarMassAsXi/F");
532     fTreeCascade->Branch("fTreeCascVarMassAsOmega",&fTreeCascVarMassAsOmega,"fTreeCascVarMassAsOmega/F");
533     fTreeCascade->Branch("fTreeCascVarPt",&fTreeCascVarPt,"fTreeCascVarPt/F");
534     fTreeCascade->Branch("fTreeCascVarRapXi",&fTreeCascVarRapXi,"fTreeCascVarRapXi/F");
535     fTreeCascade->Branch("fTreeCascVarRapOmega",&fTreeCascVarRapOmega,"fTreeCascVarRapOmega/F");
536     fTreeCascade->Branch("fTreeCascVarNegEta",&fTreeCascVarNegEta,"fTreeCascVarNegEta/F");
537     fTreeCascade->Branch("fTreeCascVarPosEta",&fTreeCascVarPosEta,"fTreeCascVarPosEta/F");
538     fTreeCascade->Branch("fTreeCascVarBachEta",&fTreeCascVarBachEta,"fTreeCascVarBachEta/F");
539     //-----------INFO-FOR-CUTS------------------------
540     fTreeCascade->Branch("fTreeCascVarDCACascDaughters",&fTreeCascVarDCACascDaughters,"fTreeCascVarDCACascDaughters/F");
541     fTreeCascade->Branch("fTreeCascVarDCABachToPrimVtx",&fTreeCascVarDCABachToPrimVtx,"fTreeCascVarDCABachToPrimVtx/F");
542     fTreeCascade->Branch("fTreeCascVarDCAV0Daughters",&fTreeCascVarDCAV0Daughters,"fTreeCascVarDCAV0Daughters/F");
543     fTreeCascade->Branch("fTreeCascVarDCAV0ToPrimVtx",&fTreeCascVarDCAV0ToPrimVtx,"fTreeCascVarDCAV0ToPrimVtx/F");
544     fTreeCascade->Branch("fTreeCascVarDCAPosToPrimVtx",&fTreeCascVarDCAPosToPrimVtx,"fTreeCascVarDCAPosToPrimVtx/F");
545     fTreeCascade->Branch("fTreeCascVarDCANegToPrimVtx",&fTreeCascVarDCANegToPrimVtx,"fTreeCascVarDCANegToPrimVtx/F");
546     fTreeCascade->Branch("fTreeCascVarCascCosPointingAngle",&fTreeCascVarCascCosPointingAngle,"fTreeCascVarCascCosPointingAngle/F");
547     fTreeCascade->Branch("fTreeCascVarCascRadius",&fTreeCascVarCascRadius,"fTreeCascVarCascRadius/F");
548     fTreeCascade->Branch("fTreeCascVarV0Mass",&fTreeCascVarV0Mass,"fTreeCascVarV0Mass/F");
549     fTreeCascade->Branch("fTreeCascVarV0CosPointingAngle",&fTreeCascVarV0CosPointingAngle,"fTreeCascVarV0CosPointingAngle/F");
550     fTreeCascade->Branch("fTreeCascVarV0CosPointingAngleSpecial",&fTreeCascVarV0CosPointingAngleSpecial,"fTreeCascVarV0CosPointingAngleSpecial/F");
551     fTreeCascade->Branch("fTreeCascVarV0Radius",&fTreeCascVarV0Radius,"fTreeCascVarV0Radius/F");
552     fTreeCascade->Branch("fTreeCascVarLeastNbrClusters",&fTreeCascVarLeastNbrClusters,"fTreeCascVarLeastNbrClusters/I");
553     //-----------MULTIPLICITY-INFO--------------------
554     fTreeCascade->Branch("fTreeCascVarCentV0A",&fTreeCascVarCentV0A,"fTreeCascVarCentV0A/F");
555     fTreeCascade->Branch("fTreeCascVarCentV0C",&fTreeCascVarCentV0C,"fTreeCascVarCentV0C/F");
556     fTreeCascade->Branch("fTreeCascVarCentV0M",&fTreeCascVarCentV0M,"fTreeCascVarCentV0M/F");
557     fTreeCascade->Branch("fTreeCascVarCentV0AEq",&fTreeCascVarCentV0AEq,"fTreeCascVarCentV0AEq/F");
558     fTreeCascade->Branch("fTreeCascVarCentV0CEq",&fTreeCascVarCentV0CEq,"fTreeCascVarCentV0CEq/F");
559     fTreeCascade->Branch("fTreeCascVarCentV0MEq",&fTreeCascVarCentV0MEq,"fTreeCascVarCentV0MEq/F");
560     fTreeCascade->Branch("fTreeCascVarCustomCentV0M",&fTreeCascVarCustomCentV0M,"fTreeCascVarCustomCentV0M/F");
561     fTreeCascade->Branch("fTreeCascVarAmpV0A",&fTreeCascVarAmpV0A,"fTreeCascVarAmpV0A/F");
562     fTreeCascade->Branch("fTreeCascVarAmpV0C",&fTreeCascVarAmpV0C,"fTreeCascVarAmpV0C/F");
563     fTreeCascade->Branch("fTreeCascVarAmpV0AEq",&fTreeCascVarAmpV0AEq,"fTreeCascVarAmpV0AEq/F");
564     fTreeCascade->Branch("fTreeCascVarAmpV0CEq",&fTreeCascVarAmpV0CEq,"fTreeCascVarAmpV0CEq/F");
565     fTreeCascade->Branch("fTreeCascVarRefMultEta8",&fTreeCascVarRefMultEta8,"fTreeCascVarRefMultEta8/I");
566     fTreeCascade->Branch("fTreeCascVarRefMultEta5",&fTreeCascVarRefMultEta5,"fTreeCascVarRefMultEta5/I");
567     fTreeCascade->Branch("fTreeCascVarRunNumber",&fTreeCascVarRunNumber,"fTreeCascVarRunNumber/I");
568     //-----------DECAY-LENGTH-INFO--------------------
569     fTreeCascade->Branch("fTreeCascVarDistOverTotMom",&fTreeCascVarDistOverTotMom,"fTreeCascVarDistOverTotMom/F");
570     //------------------------------------------------
571     fTreeCascade->Branch("fTreeCascVarNegNSigmaPion",&fTreeCascVarNegNSigmaPion,"fTreeCascVarNegNSigmaPion/F");
572     fTreeCascade->Branch("fTreeCascVarNegNSigmaProton",&fTreeCascVarNegNSigmaProton,"fTreeCascVarNegNSigmaProton/F");
573     fTreeCascade->Branch("fTreeCascVarPosNSigmaPion",&fTreeCascVarPosNSigmaPion,"fTreeCascVarPosNSigmaPion/F");
574     fTreeCascade->Branch("fTreeCascVarPosNSigmaProton",&fTreeCascVarPosNSigmaProton,"fTreeCascVarPosNSigmaProton/F");
575     fTreeCascade->Branch("fTreeCascVarBachNSigmaPion",&fTreeCascVarBachNSigmaPion,"fTreeCascVarBachNSigmaPion/F");
576     fTreeCascade->Branch("fTreeCascVarBachNSigmaKaon",&fTreeCascVarBachNSigmaKaon,"fTreeCascVarBachNSigmaKaon/F");
577     
578     //------------------------------------------------
579     // Particle Identification Setup
580     //------------------------------------------------
581     
582     AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
583     AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
584     fPIDResponse = inputHandler->GetPIDResponse();
585     
586     // Multiplicity
587     if(! fESDtrackCuts ){
588         fESDtrackCuts = new AliESDtrackCuts();
589     }
590     //Helper
591     if(! fPPVsMultUtils ){
592         fPPVsMultUtils = new AliPPVsMultUtils();
593     }
594     //Analysis Utils
595     if(! fUtils ){
596         fUtils = new AliAnalysisUtils();
597     }
598     
599     //------------------------------------------------
600     // V0 Multiplicity Histograms
601     //------------------------------------------------
602     
603     // Create histograms
604     OpenFile(1);
605     fListHist = new TList();
606     fListHist->SetOwner();  // See http://root.cern.ch/root/html/TCollection.html#TCollection:SetOwner
607     
608     if(! fHistEventCounter ) {
609         //Histogram Output: Event-by-Event
610         fHistEventCounter = new TH1D( "fHistEventCounter", ";Evt. Sel. Step;Count",5,0,5);
611         fHistEventCounter->GetXaxis()->SetBinLabel(1, "Processed");
612         fHistEventCounter->GetXaxis()->SetBinLabel(2, "Phys-Sel");
613         fHistEventCounter->GetXaxis()->SetBinLabel(3, "Has Vtx");
614         fHistEventCounter->GetXaxis()->SetBinLabel(4, "Vtx |z|<10cm");
615         fHistEventCounter->GetXaxis()->SetBinLabel(5, "Isn't Pileup");
616         fListHist->Add(fHistEventCounter);
617     }
618     
619     //List of Histograms: Normal
620     PostData(1, fListHist);
621     
622     //TTree Object: Saved to base directory. Should cache to disk while saving.
623     //(Important to avoid excessive memory usage, particularly when merging)
624     PostData(2, fTreeEvent);
625     PostData(3, fTreeV0);
626     PostData(4, fTreeCascade);
627     
628 }// end UserCreateOutputObjects
629
630
631 //________________________________________________________________________
632 void AliAnalysisTaskStrangenessVsMultiplicity::UserExec(Option_t *)
633 {
634     // Main loop
635     // Called for each event
636     
637     AliESDEvent *lESDevent = 0x0;
638     
639     //Zero all booleans, etc: safe initialization per event
640     fEvSel_HasAtLeastSPDVertex    = kFALSE;
641     fEvSel_VtxZCut                = kFALSE;
642     fEvSel_IsNotPileup            = kFALSE;
643     fEvSel_IsNotPileupMV          = kFALSE;
644     fEvSel_IsNotPileupInMultBins  = kFALSE;
645     fEvSel_Triggered              = kFALSE;
646     
647     fEvSel_nTracklets   = -1;
648     fEvSel_nSPDClusters = -1;
649     fEvSel_nContributors = -1;
650     fEvSel_nContributorsPileup = -1;
651     fEvSel_nSPDPrimVertices = -1;
652     fEvSel_distZ = -100;
653     fEvSel_VtxZ = -100;
654     
655     // Connect to the InputEvent
656     // After these lines, we should have an ESD/AOD event + the number of V0s in it.
657     
658     // Appropriate for ESD analysis!
659     
660     lESDevent = dynamic_cast<AliESDEvent*>( InputEvent() );
661     if (!lESDevent) {
662         AliWarning("ERROR: lESDevent not available \n");
663         return;
664     }
665     
666     //Get VZERO Information for multiplicity later
667     AliVVZERO* esdV0 = lESDevent->GetVZEROData();
668     if (!esdV0) {
669         AliError("AliVVZERO not available");
670         return;
671     }
672     
673     fRunNumber = lESDevent->GetRunNumber();
674     
675     Double_t lMagneticField = -10;
676     lMagneticField = lESDevent->GetMagneticField( );
677     
678     //------------------------------------------------
679     // Variable Definition
680     //------------------------------------------------
681     
682     //------------------------------------------------
683     // Physics Selection
684     //------------------------------------------------
685     
686     fHistEventCounter->Fill(0.5);
687     
688     UInt_t maskIsSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
689     Bool_t isSelected = 0;
690     isSelected = (maskIsSelected & AliVEvent::kMB) == AliVEvent::kMB;
691     fEvSel_Triggered = isSelected;
692     
693     //Standard Min-Bias Selection - always do this!
694     if ( (! isSelected) && (! fkSkipEventSelection ) ) {
695         PostData(1, fListHist);
696         PostData(2, fTreeEvent);
697         PostData(3, fTreeV0);
698         PostData(4, fTreeCascade);
699         return;
700     }
701     
702     //Tracklets vs Clusters cut via AliAnalysisUtils
703     if ( fkApplyTrackletsVsClustersCut && (! fkSkipEventSelection ) ) {
704         if( fUtils->IsSPDClusterVsTrackletBG( lESDevent ) ){
705             PostData(1, fListHist);
706             PostData(2, fTreeEvent);
707             PostData(3, fTreeV0);
708             PostData(4, fTreeCascade);
709             return;
710         }
711     }
712     
713     
714     fHistEventCounter->Fill(1.5);
715     
716     //------------------------------------------------
717     // Primary Vertex Requirements Section:
718     //  ---> pp: has vertex, |z|<10cm
719     //------------------------------------------------
720     
721     //classical Proton-proton like selection
722     const AliESDVertex *lPrimaryBestESDVtx     = lESDevent->GetPrimaryVertex();
723     const AliESDVertex *lPrimaryTrackingESDVtx = lESDevent->GetPrimaryVertexTracks();
724     const AliESDVertex *lPrimarySPDVtx         = lESDevent->GetPrimaryVertexSPD();
725     
726     Double_t lBestPrimaryVtxPos[3]          = {-100.0, -100.0, -100.0};
727     lPrimaryBestESDVtx->GetXYZ( lBestPrimaryVtxPos );
728     
729     //Only accept if Tracking or SPD vertex is fine
730     if (!lPrimarySPDVtx->GetStatus() && !lPrimaryTrackingESDVtx->GetStatus() && !fkSkipEventSelection ){
731         AliWarning("Pb / No SPD prim. vertex nor prim. Tracking vertex ... return !");
732         PostData(1, fListHist);
733         PostData(2, fTreeEvent);
734         PostData(3, fTreeV0);
735         PostData(4, fTreeCascade);
736         return;
737     }
738     
739     if(! (!lPrimarySPDVtx->GetStatus() && !lPrimaryTrackingESDVtx->GetStatus()) ){
740         //Passed selection!
741         fEvSel_HasAtLeastSPDVertex = kTRUE;
742     }
743     
744     //Has SPD or Tracking Vertex
745     fHistEventCounter -> Fill(2.5);
746     
747     //Always do Primary Vertex Selection
748     if(TMath::Abs(lBestPrimaryVtxPos[2]) > 10.0 && !fkSkipEventSelection ) {
749         AliWarning("Pb / | Z position of Best Prim Vtx | > 10.0 cm ... return !");
750         PostData(1, fListHist);
751         PostData(2, fTreeEvent);
752         PostData(3, fTreeV0);
753         PostData(4, fTreeCascade);
754         return;
755     }
756     
757     if(TMath::Abs(lBestPrimaryVtxPos[2]) <= 10.0 ){
758         //Passed selection!
759         fEvSel_VtxZCut = kTRUE;
760     }
761     fEvSel_VtxZ = lBestPrimaryVtxPos[2] ; //Set
762     
763     //Fill Event selected counter
764     fHistEventCounter -> Fill(3.5);
765     
766     //------------------------------------------------
767     // Check if this isn't pileup
768     //------------------------------------------------
769     
770     if(lESDevent->IsPileupFromSPD() && !fkSkipEventSelection ){
771         // minContributors=3, minZdist=0.8, nSigmaZdist=3., nSigmaDiamXY=2., nSigmaDiamZ=5.
772         //-> see http://alisoft.cern.ch/viewvc/trunk/STEER/AliESDEvent.h?root=AliRoot&r1=41914&r2=42199&pathrev=42199
773         AliWarning("Pb / Event tagged as pile-up by SPD... return !");
774         PostData(1, fListHist);
775         PostData(2, fTreeEvent);
776         PostData(3, fTreeV0);
777         PostData(4, fTreeCascade);
778         return;
779     }
780     
781     if( !lESDevent->IsPileupFromSPD()           ) fEvSel_IsNotPileup           = kTRUE;
782     if( !lESDevent->IsPileupFromSPDInMultBins() ) fEvSel_IsNotPileupInMultBins = kTRUE;
783     
784     //Acquire information to compute residual pileup
785     fEvSel_nSPDPrimVertices = lESDevent->GetNumberOfPileupVerticesSPD();
786     
787     //Long_t lNcontributorsSPDvtx = lPrimarySPDVtx -> GetNContributors();
788     Long_t lNcontributorsSecondLargest = -1;
789     Long_t lIndexSecondLargest = -1;
790     //Look for the two events with the largest numbers of contributors...
791     for(Int_t i=0; i<fEvSel_nSPDPrimVertices;i++){
792         const AliESDVertex* pv=lESDevent -> GetPileupVertexSPD(i);
793         if( pv->GetNContributors() > lNcontributorsSecondLargest ){
794             lNcontributorsSecondLargest = pv->GetNContributors();
795             lIndexSecondLargest = i;
796         }
797     }
798     fEvSel_nContributors = lPrimaryBestESDVtx -> GetNContributors();
799     if( fEvSel_nSPDPrimVertices > 0 && lIndexSecondLargest > -1){
800         const AliESDVertex* largestpv=lESDevent ->GetPileupVertexSPD(lIndexSecondLargest);
801         fEvSel_distZ = lPrimarySPDVtx->GetZ() - largestpv->GetZ();
802         fEvSel_nContributorsPileup = largestpv -> GetNContributors();
803     }
804     
805     //First implementation of pileup from multi-vertexer (simple use of analysis utils)
806     if ( !fUtils->IsPileUpMV( lESDevent ) ) fEvSel_IsNotPileupMV = kTRUE;
807     
808     //Fill Event isn't pileup counter
809     fHistEventCounter -> Fill(4.5);
810     
811     //------------------------------------------------
812     // Multiplicity Information Acquistion
813     //------------------------------------------------
814     
815     //Standard GetReferenceMultiplicity Estimator (0.5 and 0.8)
816     fRefMultEta5 = fESDtrackCuts->GetReferenceMultiplicity(lESDevent, AliESDtrackCuts::kTrackletsITSTPC,0.5);
817     fRefMultEta8 = fESDtrackCuts->GetReferenceMultiplicity(lESDevent, AliESDtrackCuts::kTrackletsITSTPC,0.8);
818     
819     // VZERO PART
820     Float_t  multV0A  = 0;            //  multiplicity from V0 reco side A
821     Float_t  multV0C  = 0;            //  multiplicity from V0 reco side C
822     Float_t  multV0AEq  = 0;          //  multiplicity from V0 reco side A
823     Float_t  multV0CEq  = 0;          //  multiplicity from V0 reco side C
824     Float_t  multV0ACorr  = 0;            //  multiplicity from V0 reco side A
825     Float_t  multV0CCorr  = 0;            //  multiplicity from V0 reco side C
826     
827     //Non-Equalized Signal: copy of multV0ACorr and multV0CCorr from AliCentralitySelectionTask
828     //Getters for uncorrected multiplicity
829     multV0A=esdV0->GetMTotV0A();
830     multV0C=esdV0->GetMTotV0C();
831     
832     //Get Z vertex position of SPD vertex (why not Tracking if available?)
833     Float_t zvtx = lPrimarySPDVtx->GetZ();
834     
835     //Acquire Corrected multV0A
836     multV0ACorr = AliESDUtils::GetCorrV0A(multV0A,zvtx);
837     multV0CCorr = AliESDUtils::GetCorrV0C(multV0C,zvtx);
838     
839     //Copy to Event Tree for extra information
840     fAmplitude_V0A = multV0ACorr;
841     fAmplitude_V0C = multV0CCorr;
842     
843     // Equalized signals // From AliCentralitySelectionTask // Updated
844     for(Int_t iCh = 32; iCh < 64; ++iCh) {
845         Double_t mult = lESDevent->GetVZEROEqMultiplicity(iCh);
846         multV0AEq += mult;
847     }
848     for(Int_t iCh = 0; iCh < 32; ++iCh) {
849         Double_t mult = lESDevent->GetVZEROEqMultiplicity(iCh);
850         multV0CEq += mult;
851     }
852     fAmplitude_V0AEq = multV0AEq;
853     fAmplitude_V0CEq = multV0CEq;
854     
855     fCentrality_V0A   = -100;
856     fCentrality_V0C   = -100;
857     fCentrality_V0M   = -100;
858     fCentrality_V0AEq = -100;
859     fCentrality_V0CEq = -100;
860     fCentrality_V0MEq = -100;
861     
862     //AliCentrality... Check if working?
863     AliCentrality* centrality;
864     centrality = lESDevent->GetCentrality();
865     if ( !(centrality->GetQuality()>1) ){
866         fCentrality_V0A   = centrality->GetCentralityPercentile( "V0A"   );
867         fCentrality_V0C   = centrality->GetCentralityPercentile( "V0C"   );
868         fCentrality_V0M   = centrality->GetCentralityPercentile( "V0M"   );
869         fCentrality_V0AEq = centrality->GetCentralityPercentile( "V0AEq" );
870         fCentrality_V0CEq = centrality->GetCentralityPercentile( "V0CEq" );
871         fCentrality_V0MEq = centrality->GetCentralityPercentile( "V0MEq" );
872     }
873     
874     //Cross-check/debug
875     fCustomCentrality_V0M   = fPPVsMultUtils->GetMultiplicityPercentile(lESDevent,"V0M");
876     fCustomCentrality_V0MEq = fPPVsMultUtils->GetMultiplicityPercentile(lESDevent,"V0MEq");
877     
878     //Attribution for analysis entries
879     fTreeVariableCustomCentV0M = fCustomCentrality_V0M; 
880     fTreeCascVarCustomCentV0M  = fCustomCentrality_V0M; 
881     
882     //Tracklets vs Clusters Exploratory data
883     fEvSel_nTracklets     = lESDevent->GetMultiplicity()->GetNumberOfTracklets();
884     fEvSel_nSPDClusters   = lESDevent->GetNumberOfITSClusters(0) + lESDevent->GetNumberOfITSClusters(1);
885     
886     //Event-level fill
887     fTreeEvent->Fill() ;
888     
889     //STOP HERE if skipping event selections (no point in doing the rest...)
890     if( fkSkipEventSelection ){
891         PostData(1, fListHist);
892         PostData(2, fTreeEvent);
893         PostData(3, fTreeV0);
894         PostData(4, fTreeCascade);
895         return;
896     }
897     
898     //------------------------------------------------
899     
900     //------------------------------------------------
901     // Fill V0 Tree as needed
902     //------------------------------------------------
903     
904     //Variable definition
905     Int_t    lOnFlyStatus = 0;// nv0sOn = 0, nv0sOff = 0;
906     Double_t lChi2V0 = 0;
907     Double_t lDcaV0Daughters = 0, lDcaV0ToPrimVertex = 0;
908     Double_t lDcaPosToPrimVertex = 0, lDcaNegToPrimVertex = 0;
909     Double_t lV0CosineOfPointingAngle = 0;
910     Double_t lV0Radius = 0, lPt = 0;
911     Double_t lRapK0Short = 0, lRapLambda = 0;
912     Double_t lInvMassK0s = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
913     Double_t lAlphaV0 = 0, lPtArmV0 = 0;
914     
915     Double_t fMinV0Pt = 0;
916     Double_t fMaxV0Pt = 100;
917     
918     Int_t nv0s = 0;
919     nv0s = lESDevent->GetNumberOfV0s();
920     
921     for (Int_t iV0 = 0; iV0 < nv0s; iV0++) //extra-crazy test
922     {// This is the begining of the V0 loop
923         AliESDv0 *v0 = ((AliESDEvent*)lESDevent)->GetV0(iV0);
924         if (!v0) continue;
925         
926         Double_t tDecayVertexV0[3]; v0->GetXYZ(tDecayVertexV0[0],tDecayVertexV0[1],tDecayVertexV0[2]);
927         
928         Double_t tV0mom[3];
929         v0->GetPxPyPz( tV0mom[0],tV0mom[1],tV0mom[2] );
930         Double_t lV0TotalMomentum = TMath::Sqrt(
931                                                 tV0mom[0]*tV0mom[0]+tV0mom[1]*tV0mom[1]+tV0mom[2]*tV0mom[2] );
932         
933         lV0Radius = TMath::Sqrt(tDecayVertexV0[0]*tDecayVertexV0[0]+tDecayVertexV0[1]*tDecayVertexV0[1]);
934         
935         lPt = v0->Pt();
936         lRapK0Short = v0->RapK0Short();
937         lRapLambda  = v0->RapLambda();
938         if ((lPt<fMinV0Pt)||(fMaxV0Pt<lPt)) continue;
939         
940         UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
941         UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
942         
943         Double_t lMomPos[3]; v0->GetPPxPyPz(lMomPos[0],lMomPos[1],lMomPos[2]);
944         Double_t lMomNeg[3]; v0->GetNPxPyPz(lMomNeg[0],lMomNeg[1],lMomNeg[2]);
945         
946         AliESDtrack *pTrack=((AliESDEvent*)lESDevent)->GetTrack(lKeyPos);
947         AliESDtrack *nTrack=((AliESDEvent*)lESDevent)->GetTrack(lKeyNeg);
948         if (!pTrack || !nTrack) {
949             Printf("ERROR: Could not retreive one of the daughter track");
950             continue;
951         }
952         
953         //Daughter Eta for Eta selection, afterwards
954         fTreeVariableNegEta = nTrack->Eta();
955         fTreeVariablePosEta = pTrack->Eta();
956         
957         // Filter like-sign V0 (next: add counter and distribution)
958         if ( pTrack->GetSign() == nTrack->GetSign()){
959             continue;
960         }
961         
962         //________________________________________________________________________
963         // Track quality cuts
964         Float_t lPosTrackCrossedRows = pTrack->GetTPCClusterInfo(2,1);
965         Float_t lNegTrackCrossedRows = nTrack->GetTPCClusterInfo(2,1);
966         fTreeVariableLeastNbrCrossedRows = (Int_t) lPosTrackCrossedRows;
967         if( lNegTrackCrossedRows < fTreeVariableLeastNbrCrossedRows )
968             fTreeVariableLeastNbrCrossedRows = (Int_t) lNegTrackCrossedRows;
969         
970         // TPC refit condition (done during reconstruction for Offline but not for On-the-fly)
971         if( !(pTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
972         if( !(nTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
973         
974         
975         if ( ( ( pTrack->GetTPCClusterInfo(2,1) ) < 70 ) || ( ( nTrack->GetTPCClusterInfo(2,1) ) < 70 ) ) continue;
976         
977         //GetKinkIndex condition
978         if( pTrack->GetKinkIndex(0)>0 || nTrack->GetKinkIndex(0)>0 ) continue;
979         
980         //Findable clusters > 0 condition
981         if( pTrack->GetTPCNclsF()<=0 || nTrack->GetTPCNclsF()<=0 ) continue;
982         
983         //Compute ratio Crossed Rows / Findable clusters
984         //Note: above test avoids division by zero!
985         Float_t lPosTrackCrossedRowsOverFindable = lPosTrackCrossedRows / ((double)(pTrack->GetTPCNclsF()));
986         Float_t lNegTrackCrossedRowsOverFindable = lNegTrackCrossedRows / ((double)(nTrack->GetTPCNclsF()));
987         
988         fTreeVariableLeastRatioCrossedRowsOverFindable = lPosTrackCrossedRowsOverFindable;
989         if( lNegTrackCrossedRowsOverFindable < fTreeVariableLeastRatioCrossedRowsOverFindable )
990             fTreeVariableLeastRatioCrossedRowsOverFindable = lNegTrackCrossedRowsOverFindable;
991         
992         //Lowest Cut Level for Ratio Crossed Rows / Findable = 0.8, set here
993         if ( fTreeVariableLeastRatioCrossedRowsOverFindable < 0.8 ) continue;
994         
995         //End track Quality Cuts
996         //________________________________________________________________________
997         
998         lDcaPosToPrimVertex = TMath::Abs(pTrack->GetD(lBestPrimaryVtxPos[0],
999                                                       lBestPrimaryVtxPos[1],
1000                                                       lMagneticField) );
1001         
1002         lDcaNegToPrimVertex = TMath::Abs(nTrack->GetD(lBestPrimaryVtxPos[0],
1003                                                       lBestPrimaryVtxPos[1],
1004                                                       lMagneticField) );
1005         
1006         lOnFlyStatus = v0->GetOnFlyStatus();
1007         lChi2V0 = v0->GetChi2V0();
1008         lDcaV0Daughters = v0->GetDcaV0Daughters();
1009         lDcaV0ToPrimVertex = v0->GetD(lBestPrimaryVtxPos[0],lBestPrimaryVtxPos[1],lBestPrimaryVtxPos[2]);
1010         lV0CosineOfPointingAngle = v0->GetV0CosineOfPointingAngle(lBestPrimaryVtxPos[0],lBestPrimaryVtxPos[1],lBestPrimaryVtxPos[2]);
1011         fTreeVariableV0CosineOfPointingAngle=lV0CosineOfPointingAngle;
1012         
1013         // Getting invariant mass infos directly from ESD
1014         v0->ChangeMassHypothesis(310);
1015         lInvMassK0s = v0->GetEffMass();
1016         v0->ChangeMassHypothesis(3122);
1017         lInvMassLambda = v0->GetEffMass();
1018         v0->ChangeMassHypothesis(-3122);
1019         lInvMassAntiLambda = v0->GetEffMass();
1020         lAlphaV0 = v0->AlphaV0();
1021         lPtArmV0 = v0->PtArmV0();
1022         
1023         fTreeVariablePt = v0->Pt();
1024         fTreeVariableChi2V0 = lChi2V0;
1025         fTreeVariableDcaV0ToPrimVertex = lDcaV0ToPrimVertex;
1026         fTreeVariableDcaV0Daughters = lDcaV0Daughters;
1027         fTreeVariableV0CosineOfPointingAngle = lV0CosineOfPointingAngle;
1028         fTreeVariableV0Radius = lV0Radius;
1029         fTreeVariableDcaPosToPrimVertex = lDcaPosToPrimVertex;
1030         fTreeVariableDcaNegToPrimVertex = lDcaNegToPrimVertex;
1031         fTreeVariableInvMassK0s = lInvMassK0s;
1032         fTreeVariableInvMassLambda = lInvMassLambda;
1033         fTreeVariableInvMassAntiLambda = lInvMassAntiLambda;
1034         fTreeVariableRapK0Short = lRapK0Short;
1035         fTreeVariableRapLambda = lRapLambda;
1036         fTreeVariableAlphaV0 = lAlphaV0;
1037         fTreeVariablePtArmV0 = lPtArmV0;
1038         
1039         //Official means of acquiring N-sigmas
1040         fTreeVariableNSigmasPosProton = fPIDResponse->NumberOfSigmasTPC( pTrack, AliPID::kProton );
1041         fTreeVariableNSigmasPosPion   = fPIDResponse->NumberOfSigmasTPC( pTrack, AliPID::kPion );
1042         fTreeVariableNSigmasNegProton = fPIDResponse->NumberOfSigmasTPC( nTrack, AliPID::kProton );
1043         fTreeVariableNSigmasNegPion   = fPIDResponse->NumberOfSigmasTPC( nTrack, AliPID::kPion );
1044         
1045         //This requires an Invariant Mass Hypothesis afterwards
1046         fTreeVariableDistOverTotMom = TMath::Sqrt(
1047                                                   TMath::Power( tDecayVertexV0[0] - lBestPrimaryVtxPos[0] , 2) +
1048                                                   TMath::Power( tDecayVertexV0[1] - lBestPrimaryVtxPos[1] , 2) +
1049                                                   TMath::Power( tDecayVertexV0[2] - lBestPrimaryVtxPos[2] , 2)
1050                                                   );
1051         fTreeVariableDistOverTotMom /= (lV0TotalMomentum+1e-10); //avoid division by zero, to be sure
1052         
1053         //Copy Multiplicity information
1054         fTreeVariableCentV0A = fCentrality_V0A;
1055         fTreeVariableCentV0C = fCentrality_V0C;
1056         fTreeVariableCentV0M = fCentrality_V0M;
1057         fTreeVariableCentV0AEq = fCentrality_V0AEq;
1058         fTreeVariableCentV0CEq = fCentrality_V0CEq;
1059         fTreeVariableCentV0MEq = fCentrality_V0MEq;
1060         fTreeVariableAmpV0A = fAmplitude_V0A;
1061         fTreeVariableAmpV0C = fAmplitude_V0C;
1062         fTreeVariableAmpV0AEq = fAmplitude_V0AEq;
1063         fTreeVariableAmpV0CEq = fAmplitude_V0CEq;
1064         fTreeVariableRefMultEta8 = fRefMultEta8;
1065         fTreeVariableRefMultEta5 = fRefMultEta5;
1066         fTreeVariableRunNumber = fRunNumber;
1067         
1068         //------------------------------------------------
1069         // Fill Tree!
1070         //------------------------------------------------
1071         
1072         // The conditionals are meant to decrease excessive
1073         // memory usage!
1074         
1075         //First Selection: Reject OnFly
1076         if( lOnFlyStatus == 0 ){
1077             //Second Selection: rough 20-sigma band, parametric.
1078             //K0Short: Enough to parametrize peak broadening with linear function.
1079             Double_t lUpperLimitK0Short = (5.63707e-01) + (1.14979e-02)*fTreeVariablePt;
1080             Double_t lLowerLimitK0Short = (4.30006e-01) - (1.10029e-02)*fTreeVariablePt;
1081             //Lambda: Linear (for higher pt) plus exponential (for low-pt broadening)
1082             //[0]+[1]*x+[2]*TMath::Exp(-[3]*x)
1083             Double_t lUpperLimitLambda = (1.13688e+00) + (5.27838e-03)*fTreeVariablePt + (8.42220e-02)*TMath::Exp(-(3.80595e+00)*fTreeVariablePt);
1084             Double_t lLowerLimitLambda = (1.09501e+00) - (5.23272e-03)*fTreeVariablePt - (7.52690e-02)*TMath::Exp(-(3.46339e+00)*fTreeVariablePt);
1085             //Do Selection
1086             if( (fTreeVariableInvMassLambda    < lUpperLimitLambda  && fTreeVariableInvMassLambda     > lLowerLimitLambda     ) ||
1087                (fTreeVariableInvMassAntiLambda < lUpperLimitLambda  && fTreeVariableInvMassAntiLambda > lLowerLimitLambda     ) ||
1088                (fTreeVariableInvMassK0s        < lUpperLimitK0Short && fTreeVariableInvMassK0s        > lLowerLimitK0Short    ) ){
1089                 //Pre-selection in case this is AA...
1090                 if ( TMath::Abs(fTreeVariableNegEta)<0.8 && TMath::Abs(fTreeVariablePosEta)<0.8 && fkSaveV0Tree ) fTreeV0->Fill();
1091             }
1092         }
1093     }// This is the end of the V0 loop
1094     
1095     //------------------------------------------------
1096     // Fill V0 tree over.
1097     //------------------------------------------------
1098     
1099     
1100     
1101     //------------------------------------------------
1102     // Rerun cascade vertexer!
1103     //------------------------------------------------
1104     
1105     if( fkRunVertexers ){
1106         lESDevent->ResetCascades();
1107         lESDevent->ResetV0s();
1108         
1109         AliV0vertexer lV0vtxer;
1110         AliCascadeVertexer lCascVtxer;
1111         
1112         lV0vtxer.SetDefaultCuts(fV0VertexerSels);
1113         lCascVtxer.SetDefaultCuts(fCascadeVertexerSels);
1114         
1115         lV0vtxer.Tracks2V0vertices(lESDevent);
1116         lCascVtxer.V0sTracks2CascadeVertices(lESDevent);
1117     }
1118     
1119     //------------------------------------------------
1120     // MAIN CASCADE LOOP STARTS HERE
1121     //------------------------------------------------
1122     // Code Credit: Antonin Maire (thanks^100)
1123     // ---> This is an adaptation
1124     
1125     Long_t ncascades = 0;
1126         ncascades = lESDevent->GetNumberOfCascades();
1127     
1128     for (Int_t iXi = 0; iXi < ncascades; iXi++){
1129         //------------------------------------------------
1130         // Initializations
1131         //------------------------------------------------
1132         //Double_t lTrkgPrimaryVtxRadius3D = -500.0;
1133         //Double_t lBestPrimaryVtxRadius3D = -500.0;
1134         
1135         // - 1st part of initialisation : variables needed to store AliESDCascade data members
1136         Double_t lEffMassXi      = 0. ;
1137         //Double_t lChi2Xi         = -1. ;
1138         Double_t lDcaXiDaughters = -1. ;
1139         Double_t lXiCosineOfPointingAngle = -1. ;
1140         Double_t lPosXi[3] = { -1000.0, -1000.0, -1000.0 };
1141         Double_t lXiRadius = -1000. ;
1142         
1143         // - 2nd part of initialisation : Nbr of clusters within TPC for the 3 daughter cascade tracks
1144         Int_t    lPosTPCClusters    = -1; // For ESD only ...//FIXME : wait for availability in AOD
1145         Int_t    lNegTPCClusters    = -1; // For ESD only ...
1146         Int_t    lBachTPCClusters   = -1; // For ESD only ...
1147         
1148         // - 3rd part of initialisation : about V0 part in cascades
1149         Double_t lInvMassLambdaAsCascDghter = 0.;
1150         //Double_t lV0Chi2Xi         = -1. ;
1151         Double_t lDcaV0DaughtersXi = -1.;
1152                 
1153         Double_t lDcaBachToPrimVertexXi = -1., lDcaV0ToPrimVertexXi = -1.;
1154         Double_t lDcaPosToPrimVertexXi  = -1.;
1155         Double_t lDcaNegToPrimVertexXi  = -1.;
1156         Double_t lV0CosineOfPointingAngleXi = -1. ;
1157         Double_t lV0CosineOfPointingAngleXiSpecial = -1. ;
1158         Double_t lPosV0Xi[3] = { -1000. , -1000., -1000. }; // Position of VO coming from cascade
1159         Double_t lV0RadiusXi = -1000.0;
1160         Double_t lV0quality  = 0.;
1161         
1162         // - 4th part of initialisation : Effective masses
1163         Double_t lInvMassXiMinus    = 0.;
1164         Double_t lInvMassXiPlus     = 0.;
1165         Double_t lInvMassOmegaMinus = 0.;
1166         Double_t lInvMassOmegaPlus  = 0.;
1167         
1168         // - 6th part of initialisation : extra info for QA
1169         Double_t lXiMomX       = 0. , lXiMomY = 0., lXiMomZ = 0.;
1170         Double_t lXiTransvMom  = 0. ;
1171         //Double_t lXiTransvMomMC= 0. ;
1172         Double_t lXiTotMom     = 0. ;
1173                 
1174         Double_t lBachMomX       = 0., lBachMomY  = 0., lBachMomZ   = 0.;
1175         //Double_t lBachTransvMom  = 0.;
1176         //Double_t lBachTotMom     = 0.;
1177         
1178         fTreeCascVarNegNSigmaPion   = -100;
1179         fTreeCascVarNegNSigmaProton = -100;
1180         fTreeCascVarPosNSigmaPion   = -100;
1181         fTreeCascVarPosNSigmaProton = -100;
1182         fTreeCascVarBachNSigmaPion  = -100;
1183         fTreeCascVarBachNSigmaKaon  = -100;
1184         
1185         Short_t  lChargeXi = -2;
1186         //Double_t lV0toXiCosineOfPointingAngle = 0. ;
1187         
1188         Double_t lRapXi   = -20.0, lRapOmega = -20.0; //  lEta = -20.0, lTheta = 360., lPhi = 720. ;
1189         //Double_t lAlphaXi = -200., lPtArmXi  = -200.0;
1190             
1191         // -------------------------------------
1192         // II.ESD - Calculation Part dedicated to Xi vertices (ESD)
1193         
1194         AliESDcascade *xi = lESDevent->GetCascade(iXi);
1195         if (!xi) continue;
1196         
1197                 // - II.Step 2 : Assigning the necessary variables for specific AliESDcascade data members (ESD)
1198                 //-------------
1199         lV0quality = 0.;
1200         xi->ChangeMassHypothesis(lV0quality , 3312); // default working hypothesis : cascade = Xi- decay
1201         
1202         lEffMassXi                      = xi->GetEffMassXi();
1203         //lChi2Xi                           = xi->GetChi2Xi();
1204         lDcaXiDaughters         = xi->GetDcaXiDaughters();
1205         lXiCosineOfPointingAngle                    = xi->GetCascadeCosineOfPointingAngle( lBestPrimaryVtxPos[0],
1206                                                                                       lBestPrimaryVtxPos[1],
1207                                                                                       lBestPrimaryVtxPos[2] );
1208         // Take care : the best available vertex should be used (like in AliCascadeVertexer)
1209         
1210         xi->GetXYZcascade( lPosXi[0],  lPosXi[1], lPosXi[2] );
1211         lXiRadius                       = TMath::Sqrt( lPosXi[0]*lPosXi[0]  +  lPosXi[1]*lPosXi[1] );
1212         
1213                 // - II.Step 3 : around the tracks : Bach + V0 (ESD)
1214                 // ~ Necessary variables for ESDcascade data members coming from the ESDv0 part (inheritance)
1215                 //-------------
1216                 
1217         UInt_t lIdxPosXi        = (UInt_t) TMath::Abs( xi->GetPindex() );
1218         UInt_t lIdxNegXi        = (UInt_t) TMath::Abs( xi->GetNindex() );
1219         UInt_t lBachIdx         = (UInt_t) TMath::Abs( xi->GetBindex() );
1220         // Care track label can be negative in MC production (linked with the track quality)
1221         // However = normally, not the case for track index ...
1222         
1223         // FIXME : rejection of a double use of a daughter track (nothing but just a crosscheck of what is done in the cascade vertexer)
1224         if(lBachIdx == lIdxNegXi) {
1225             AliWarning("Pb / Idx(Bach. track) = Idx(Neg. track) ... continue!"); continue;
1226         }
1227         if(lBachIdx == lIdxPosXi) {
1228             AliWarning("Pb / Idx(Bach. track) = Idx(Pos. track) ... continue!"); continue;
1229         }
1230         
1231         AliESDtrack *pTrackXi           = lESDevent->GetTrack( lIdxPosXi );
1232         AliESDtrack *nTrackXi           = lESDevent->GetTrack( lIdxNegXi );
1233         AliESDtrack *bachTrackXi        = lESDevent->GetTrack( lBachIdx );
1234         
1235         if (!pTrackXi || !nTrackXi || !bachTrackXi ) {
1236             AliWarning("ERROR: Could not retrieve one of the 3 ESD daughter tracks of the cascade ...");
1237             continue;
1238         }
1239         
1240         fTreeCascVarPosEta = pTrackXi->Eta();
1241         fTreeCascVarNegEta = nTrackXi->Eta();
1242         fTreeCascVarBachEta = bachTrackXi->Eta();
1243         
1244         Double_t lBMom[3], lNMom[3], lPMom[3];
1245         xi->GetBPxPyPz( lBMom[0], lBMom[1], lBMom[2] );
1246         xi->GetPPxPyPz( lPMom[0], lPMom[1], lPMom[2] );
1247         xi->GetNPxPyPz( lNMom[0], lNMom[1], lNMom[2] );
1248         
1249         //fTreeCascVarBachTransMom = TMath::Sqrt( lBMom[0]*lBMom[0] + lBMom[1]*lBMom[1] );
1250         //fTreeCascVarPosTransMom  = TMath::Sqrt( lPMom[0]*lPMom[0] + lPMom[1]*lPMom[1] );
1251         //fTreeCascVarNegTransMom  = TMath::Sqrt( lNMom[0]*lNMom[0] + lNMom[1]*lNMom[1] );
1252         
1253         //------------------------------------------------
1254         // TPC dEdx information
1255         //------------------------------------------------
1256         fTreeCascVarNegNSigmaPion   = fPIDResponse->NumberOfSigmasTPC( nTrackXi, AliPID::kPion   );
1257         fTreeCascVarNegNSigmaProton = fPIDResponse->NumberOfSigmasTPC( nTrackXi, AliPID::kProton );
1258         fTreeCascVarPosNSigmaPion   = fPIDResponse->NumberOfSigmasTPC( pTrackXi, AliPID::kPion );
1259         fTreeCascVarPosNSigmaProton = fPIDResponse->NumberOfSigmasTPC( pTrackXi, AliPID::kProton );
1260         fTreeCascVarBachNSigmaPion  = fPIDResponse->NumberOfSigmasTPC( bachTrackXi, AliPID::kPion );
1261         fTreeCascVarBachNSigmaKaon  = fPIDResponse->NumberOfSigmasTPC( bachTrackXi, AliPID::kKaon );
1262         
1263         //------------------------------------------------
1264         // TPC Number of clusters info
1265         // --- modified to save the smallest number
1266         // --- of TPC clusters for the 3 tracks
1267         //------------------------------------------------
1268         
1269         lPosTPCClusters   = pTrackXi->GetTPCNcls();
1270         lNegTPCClusters   = nTrackXi->GetTPCNcls();
1271         lBachTPCClusters  = bachTrackXi->GetTPCNcls();
1272         
1273         // 1 - Poor quality related to TPCrefit
1274         ULong_t pStatus    = pTrackXi->GetStatus();
1275         ULong_t nStatus    = nTrackXi->GetStatus();
1276         ULong_t bachStatus = bachTrackXi->GetStatus();
1277         
1278         //fTreeCascVarkITSRefitBachelor = kTRUE;
1279         //fTreeCascVarkITSRefitNegative = kTRUE;
1280         //fTreeCascVarkITSRefitPositive = kTRUE;
1281         
1282         if ((pStatus&AliESDtrack::kTPCrefit)    == 0) { AliWarning("Pb / V0 Pos. track has no TPCrefit ... continue!"); continue; }
1283         if ((nStatus&AliESDtrack::kTPCrefit)    == 0) { AliWarning("Pb / V0 Neg. track has no TPCrefit ... continue!"); continue; }
1284         if ((bachStatus&AliESDtrack::kTPCrefit) == 0) { AliWarning("Pb / Bach.   track has no TPCrefit ... continue!"); continue; }
1285         
1286         // 2 - Poor quality related to TPC clusters: lowest cut of 70 clusters
1287         if(lPosTPCClusters  < 70) { AliWarning("Pb / V0 Pos. track has less than 70 TPC clusters ... continue!"); continue; }
1288         if(lNegTPCClusters  < 70) { AliWarning("Pb / V0 Neg. track has less than 70 TPC clusters ... continue!"); continue; }
1289         if(lBachTPCClusters < 70) { AliWarning("Pb / Bach.   track has less than 70 TPC clusters ... continue!"); continue; }
1290         Int_t leastnumberofclusters = 1000;
1291         if( lPosTPCClusters < leastnumberofclusters ) leastnumberofclusters = lPosTPCClusters;
1292         if( lNegTPCClusters < leastnumberofclusters ) leastnumberofclusters = lNegTPCClusters;
1293         if( lBachTPCClusters < leastnumberofclusters ) leastnumberofclusters = lBachTPCClusters;
1294         
1295         lInvMassLambdaAsCascDghter      = xi->GetEffMass();
1296         // This value shouldn't change, whatever the working hyp. is : Xi-, Xi+, Omega-, Omega+
1297         lDcaV0DaughtersXi               = xi->GetDcaV0Daughters();
1298         //lV0Chi2Xi                     = xi->GetChi2V0();
1299         
1300         lV0CosineOfPointingAngleXi      = xi->GetV0CosineOfPointingAngle( lBestPrimaryVtxPos[0],
1301                                                                      lBestPrimaryVtxPos[1],
1302                                                                      lBestPrimaryVtxPos[2] );
1303         //Modification: V0 CosPA wrt to Cascade decay vertex
1304         lV0CosineOfPointingAngleXiSpecial       = xi->GetV0CosineOfPointingAngle( lPosXi[0],
1305                                                                              lPosXi[1],
1306                                                                              lPosXi[2] );
1307         
1308         lDcaV0ToPrimVertexXi            = xi->GetD( lBestPrimaryVtxPos[0],
1309                                                lBestPrimaryVtxPos[1],
1310                                                lBestPrimaryVtxPos[2] );
1311                 
1312         lDcaBachToPrimVertexXi = TMath::Abs( bachTrackXi->GetD( lBestPrimaryVtxPos[0],
1313                                                                lBestPrimaryVtxPos[1],
1314                                                                lMagneticField  ) );
1315         // Note : AliExternalTrackParam::GetD returns an algebraic value ...
1316                 
1317         xi->GetXYZ( lPosV0Xi[0],  lPosV0Xi[1], lPosV0Xi[2] );
1318         lV0RadiusXi             = TMath::Sqrt( lPosV0Xi[0]*lPosV0Xi[0]  +  lPosV0Xi[1]*lPosV0Xi[1] );
1319         
1320         lDcaPosToPrimVertexXi   = TMath::Abs( pTrackXi  ->GetD( lBestPrimaryVtxPos[0],
1321                                                                lBestPrimaryVtxPos[1],
1322                                                                lMagneticField  )     );
1323         
1324         lDcaNegToPrimVertexXi   = TMath::Abs( nTrackXi  ->GetD( lBestPrimaryVtxPos[0],
1325                                                                lBestPrimaryVtxPos[1],
1326                                                                lMagneticField  )     );
1327                 
1328         // - II.Step 4 : around effective masses (ESD)
1329         // ~ change mass hypotheses to cover all the possibilities :  Xi-/+, Omega -/+
1330                 
1331         if( bachTrackXi->Charge() < 0 ) {
1332             lV0quality = 0.;
1333             xi->ChangeMassHypothesis(lV0quality , 3312);
1334             // Calculate the effective mass of the Xi- candidate.
1335             // pdg code 3312 = Xi-
1336             lInvMassXiMinus = xi->GetEffMassXi();
1337             
1338             lV0quality = 0.;
1339             xi->ChangeMassHypothesis(lV0quality , 3334);
1340             // Calculate the effective mass of the Xi- candidate.
1341             // pdg code 3334 = Omega-
1342             lInvMassOmegaMinus = xi->GetEffMassXi();
1343             
1344             lV0quality = 0.;
1345             xi->ChangeMassHypothesis(lV0quality , 3312);        // Back to default hyp.
1346         }// end if negative bachelor
1347         
1348         
1349         if( bachTrackXi->Charge() >  0 ){
1350             lV0quality = 0.;
1351             xi->ChangeMassHypothesis(lV0quality , -3312);
1352             // Calculate the effective mass of the Xi+ candidate.
1353             // pdg code -3312 = Xi+
1354             lInvMassXiPlus = xi->GetEffMassXi();
1355             
1356             lV0quality = 0.;
1357             xi->ChangeMassHypothesis(lV0quality , -3334);
1358             // Calculate the effective mass of the Xi+ candidate.
1359             // pdg code -3334  = Omega+
1360             lInvMassOmegaPlus = xi->GetEffMassXi();
1361             
1362             lV0quality = 0.;
1363             xi->ChangeMassHypothesis(lV0quality , -3312);       // Back to "default" hyp.
1364         }// end if positive bachelor
1365         // - II.Step 6 : extra info for QA (ESD)
1366         // miscellaneous pieces of info that may help regarding data quality assessment.
1367         //-------------
1368         
1369         xi->GetPxPyPz( lXiMomX, lXiMomY, lXiMomZ );
1370         lXiTransvMom    = TMath::Sqrt( lXiMomX*lXiMomX   + lXiMomY*lXiMomY );
1371         lXiTotMom       = TMath::Sqrt( lXiMomX*lXiMomX   + lXiMomY*lXiMomY   + lXiMomZ*lXiMomZ );
1372                 
1373         xi->GetBPxPyPz(  lBachMomX,  lBachMomY,  lBachMomZ );
1374         //lBachTransvMom  = TMath::Sqrt( lBachMomX*lBachMomX   + lBachMomY*lBachMomY );
1375         //lBachTotMom   = TMath::Sqrt( lBachMomX*lBachMomX   + lBachMomY*lBachMomY  +  lBachMomZ*lBachMomZ  );
1376         
1377         lChargeXi = xi->Charge();
1378         
1379         //lV0toXiCosineOfPointingAngle = xi->GetV0CosineOfPointingAngle( lPosXi[0], lPosXi[1], lPosXi[2] );
1380         
1381         lRapXi    = xi->RapXi();
1382         lRapOmega = xi->RapOmega();
1383         //lEta      = xi->Eta();
1384         //lTheta    = xi->Theta() *180.0/TMath::Pi();
1385         //lPhi      = xi->Phi()   *180.0/TMath::Pi();
1386         //lAlphaXi  = xi->AlphaXi();
1387         //lPtArmXi  = xi->PtArmXi();
1388         
1389         //------------------------------------------------
1390         // Set Variables for adding to tree
1391         //------------------------------------------------
1392         
1393         fTreeCascVarCharge      = lChargeXi;
1394         if(lInvMassXiMinus!=0)    fTreeCascVarMassAsXi = lInvMassXiMinus;
1395         if(lInvMassXiPlus!=0)     fTreeCascVarMassAsXi = lInvMassXiPlus;
1396         if(lInvMassOmegaMinus!=0) fTreeCascVarMassAsOmega = lInvMassOmegaMinus;
1397         if(lInvMassOmegaPlus!=0)  fTreeCascVarMassAsOmega = lInvMassOmegaPlus;
1398         fTreeCascVarPt = lXiTransvMom;
1399         fTreeCascVarRapXi = lRapXi ;
1400         fTreeCascVarRapOmega = lRapOmega ;
1401         fTreeCascVarDCACascDaughters = lDcaXiDaughters;
1402         fTreeCascVarDCABachToPrimVtx = lDcaBachToPrimVertexXi;
1403         fTreeCascVarDCAV0Daughters = lDcaV0DaughtersXi;
1404         fTreeCascVarDCAV0ToPrimVtx = lDcaV0ToPrimVertexXi;
1405         fTreeCascVarDCAPosToPrimVtx = lDcaPosToPrimVertexXi;
1406         fTreeCascVarDCANegToPrimVtx = lDcaNegToPrimVertexXi;
1407         fTreeCascVarCascCosPointingAngle = lXiCosineOfPointingAngle;
1408         fTreeCascVarCascRadius = lXiRadius;
1409         fTreeCascVarV0Mass = lInvMassLambdaAsCascDghter;
1410         fTreeCascVarV0CosPointingAngle = lV0CosineOfPointingAngleXi;
1411         fTreeCascVarV0CosPointingAngleSpecial = lV0CosineOfPointingAngleXiSpecial;
1412         fTreeCascVarV0Radius = lV0RadiusXi;
1413         fTreeCascVarLeastNbrClusters = leastnumberofclusters;
1414         
1415         //Copy Multiplicity information 
1416         fTreeCascVarCentV0A = fCentrality_V0A; 
1417         fTreeCascVarCentV0C = fCentrality_V0C; 
1418         fTreeCascVarCentV0M = fCentrality_V0M; 
1419         fTreeCascVarCentV0AEq = fCentrality_V0AEq; 
1420         fTreeCascVarCentV0CEq = fCentrality_V0CEq; 
1421         fTreeCascVarCentV0MEq = fCentrality_V0MEq; 
1422         fTreeCascVarAmpV0A = fAmplitude_V0A; 
1423         fTreeCascVarAmpV0C = fAmplitude_V0C; 
1424         fTreeCascVarAmpV0AEq = fAmplitude_V0AEq; 
1425         fTreeCascVarAmpV0CEq = fAmplitude_V0CEq; 
1426         fTreeCascVarRefMultEta8 = fRefMultEta8;
1427         fTreeCascVarRefMultEta5 = fRefMultEta5;
1428         fTreeCascVarRunNumber = fRunNumber; 
1429         
1430         fTreeCascVarDistOverTotMom = TMath::Sqrt(
1431                                                  TMath::Power( lPosXi[0] - lBestPrimaryVtxPos[0] , 2) +
1432                                                  TMath::Power( lPosXi[1] - lBestPrimaryVtxPos[1] , 2) +
1433                                                  TMath::Power( lPosXi[2] - lBestPrimaryVtxPos[2] , 2)
1434                                                  );
1435         fTreeCascVarDistOverTotMom /= (lXiTotMom+1e-13);
1436         
1437         //All vars not specified here: specified elsewhere!
1438         
1439         //------------------------------------------------
1440         // Fill Tree! 
1441         //------------------------------------------------
1442         
1443         // The conditional is meant to decrease excessive
1444         // memory usage! Be careful when loosening the 
1445         // cut!
1446         
1447         //Xi    Mass window: 150MeV wide
1448         //Omega mass window: 150MeV wide
1449         
1450         if( fkSaveCascadeTree && ( (fTreeCascVarMassAsXi<1.32+0.075&&fTreeCascVarMassAsXi>1.32-0.075) ||
1451                                   (fTreeCascVarMassAsOmega<1.68+0.075&&fTreeCascVarMassAsOmega>1.68-0.075) ) ){
1452             fTreeCascade->Fill();
1453         }
1454         
1455         //------------------------------------------------
1456         // Fill tree over.
1457         //------------------------------------------------
1458         
1459         }// end of the Cascade loop (ESD or AOD)
1460     
1461     // Post output data.
1462     PostData(1, fListHist); 
1463     PostData(2, fTreeEvent);
1464     PostData(3, fTreeV0);
1465     PostData(4, fTreeCascade);
1466 }
1467
1468 //________________________________________________________________________
1469 void AliAnalysisTaskStrangenessVsMultiplicity::Terminate(Option_t *)
1470 {
1471     // Draw result to the screen
1472     // Called once at the end of the query
1473     
1474     TList *cRetrievedList = 0x0;
1475     cRetrievedList = (TList*)GetOutputData(1);
1476     if(!cRetrievedList){
1477         Printf("ERROR - AliAnalysisTaskStrangenessVsMultiplicity : ouput data container list not available\n");
1478         return;
1479     }   
1480         
1481     fHistEventCounter = dynamic_cast<TH1D*> (  cRetrievedList->FindObject("fHistEventCounter")  );
1482     if (!fHistEventCounter) {
1483         Printf("ERROR - AliAnalysisTaskStrangenessVsMultiplicity : fHistEventCounter not available");
1484         return;
1485     }
1486     
1487     TCanvas *canCheck = new TCanvas("AliAnalysisTaskStrangenessVsMultiplicity","V0 Multiplicity",10,10,510,510);
1488     canCheck->cd(1)->SetLogy();
1489     
1490     fHistEventCounter->SetMarkerStyle(22);
1491     fHistEventCounter->DrawCopy("E");
1492 }
1493
1494 //----------------------------------------------------------------------------
1495
1496 Double_t AliAnalysisTaskStrangenessVsMultiplicity::MyRapidity(Double_t rE, Double_t rPz) const
1497 {
1498     // Local calculation for rapidity
1499     Double_t ReturnValue = -100;
1500     if( (rE-rPz+1.e-13) != 0 && (rE+rPz) != 0 ){ 
1501         ReturnValue =  0.5*TMath::Log((rE+rPz)/(rE-rPz+1.e-13));
1502     }
1503     return ReturnValue;
1504