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