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