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