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