]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/SPECTRA/AliAnalysisTaskCheckCascade.cxx
update the define too...
[u/mrichter/AliRoot.git] / PWG2 / SPECTRA / AliAnalysisTaskCheckCascade.cxx
1
2 /**************************************************************************
3  *  Authors : Antonin Maire, Boris Hippolyte                              *
4  * Contributors are mentioned in the code where appropriate.              *
5  *                                                                        *
6  * Permission to use, copy, modify and distribute this software and its   *
7  * documentation strictly for non-commercial purposes is hereby granted   *
8  * without fee, provided that the above copyright notice appears in all   *
9  * copies and that both the copyright notice and this permission notice   *
10  * appear in the supporting documentation. The authors make no claims     *
11  * about the suitability of this software for any purpose. It is          *
12  * provided "as is" without express or implied warranty.                  *
13  **************************************************************************/
14
15 //-----------------------------------------------------------------
16 //                 AliAnalysisTaskCheckCascade class
17 //            This task is for QAing the Cascades from ESD and AOD 
18 //            Origin : Antonin Maire Fev2008, antonin.maire@ires.in2p3.fr
19 //            Modified : A.M Juillet 2008
20 //-----------------------------------------------------------------
21
22
23
24
25 class TTree;
26 class TParticle;
27 class TVector3;
28
29 //class AliMCEventHandler;
30 //class AliMCEvent;
31 //class AliStack;
32
33 class AliESDVertex;
34 class AliAODVertex;
35 class AliESDv0;
36 class AliAODv0;
37
38 #include <iostream>
39
40 #include "TList.h"
41 #include "TH1F.h"
42 #include "TH2F.h"
43 #include "TCanvas.h"
44 #include "TMath.h"
45
46 #include "AliLog.h"
47
48 #include "AliESDEvent.h"
49 #include "AliAODEvent.h"
50 //#include "AliCascadeVertexer.h"
51
52 #include "AliESDcascade.h"
53 #include "AliAODcascade.h"
54
55 #include "AliAnalysisTaskCheckCascade.h"
56
57 ClassImp(AliAnalysisTaskCheckCascade)
58
59
60
61 //________________________________________________________________________
62 AliAnalysisTaskCheckCascade::AliAnalysisTaskCheckCascade() 
63   : AliAnalysisTaskSE(), fAnalysisType("ESD"), fCollidingSystems(0),
64
65         // - Cascade part initialisation
66     fListHistCascade(0),
67     fHistTrackMultiplicity(0), fHistCascadeMultiplicity(0),
68     fHistVtxStatus(0),
69
70     fHistPosTrkgPrimaryVtxX(0), fHistPosTrkgPrimaryVtxY(0), fHistPosTrkgPrimaryVtxZ(0), fHistTrkgPrimaryVtxRadius(0),
71     fHistPosBestPrimaryVtxX(0), fHistPosBestPrimaryVtxY(0), fHistPosBestPrimaryVtxZ(0), fHistBestPrimaryVtxRadius(0),
72     f2dHistTrkgPrimVtxVsBestPrimVtx(0),
73
74     fHistEffMassXi(0),  fHistChi2Xi(0),  
75     fHistDcaXiDaughters(0), fHistDcaBachToPrimVertex(0), fHistXiCosineOfPointingAngle(0), fHistXiRadius(0),
76
77     fHistMassLambdaAsCascDghter(0),
78     fHistV0Chi2Xi(0),
79     fHistDcaV0DaughtersXi(0),
80     fHistDcaV0ToPrimVertexXi(0), 
81     fHistV0CosineOfPointingAngleXi(0),
82     fHistV0RadiusXi(0),
83     fHistDcaPosToPrimVertexXi(0), fHistDcaNegToPrimVertexXi(0), 
84
85     fHistMassXiMinus(0), fHistMassXiPlus(0),
86     fHistMassOmegaMinus(0), fHistMassOmegaPlus(0),
87
88     fHistXiTransvMom(0),    fHistXiTotMom(0),
89     fHistBachTransvMom(0),   fHistBachTotMom(0),
90
91     fHistChargeXi(0),
92     fHistV0toXiCosineOfPointingAngle(0),
93
94     fHistRapXi(0), fHistRapOmega(0), fHistEta(0),
95     fHistTheta(0), fHistPhi(0),
96
97     f2dHistArmenteros(0),                       
98     f2dHistEffMassLambdaVsEffMassXiMinus(0), f2dHistEffMassXiVsEffMassOmegaMinus(0),
99     f2dHistEffMassLambdaVsEffMassXiPlus(0), f2dHistEffMassXiVsEffMassOmegaPlus(0),
100     f2dHistXiRadiusVsEffMassXiMinus(0), f2dHistXiRadiusVsEffMassXiPlus(0)
101
102 {
103   // Dummy Constructor
104 }
105
106
107
108
109
110
111
112
113 //________________________________________________________________________
114 AliAnalysisTaskCheckCascade::AliAnalysisTaskCheckCascade(const char *name) 
115   : AliAnalysisTaskSE(name), fAnalysisType("ESD"), fCollidingSystems(0),
116      
117         // - Cascade part initialisation
118     fListHistCascade(0),
119     fHistTrackMultiplicity(0), fHistCascadeMultiplicity(0),
120     fHistVtxStatus(0),
121
122     fHistPosTrkgPrimaryVtxX(0), fHistPosTrkgPrimaryVtxY(0), fHistPosTrkgPrimaryVtxZ(0), fHistTrkgPrimaryVtxRadius(0),
123     fHistPosBestPrimaryVtxX(0), fHistPosBestPrimaryVtxY(0), fHistPosBestPrimaryVtxZ(0), fHistBestPrimaryVtxRadius(0),
124     f2dHistTrkgPrimVtxVsBestPrimVtx(0),
125
126     fHistEffMassXi(0),  fHistChi2Xi(0),  
127     fHistDcaXiDaughters(0), fHistDcaBachToPrimVertex(0), fHistXiCosineOfPointingAngle(0), fHistXiRadius(0),
128
129     fHistMassLambdaAsCascDghter(0),
130     fHistV0Chi2Xi(0),
131     fHistDcaV0DaughtersXi(0),
132     fHistDcaV0ToPrimVertexXi(0), 
133     fHistV0CosineOfPointingAngleXi(0),
134     fHistV0RadiusXi(0),
135     fHistDcaPosToPrimVertexXi(0), fHistDcaNegToPrimVertexXi(0), 
136
137     fHistMassXiMinus(0), fHistMassXiPlus(0),
138     fHistMassOmegaMinus(0), fHistMassOmegaPlus(0),
139
140     fHistXiTransvMom(0),    fHistXiTotMom(0),
141     fHistBachTransvMom(0),   fHistBachTotMom(0),
142
143     fHistChargeXi(0),
144     fHistV0toXiCosineOfPointingAngle(0),
145
146     fHistRapXi(0), fHistRapOmega(0), fHistEta(0),
147     fHistTheta(0), fHistPhi(0),
148
149     f2dHistArmenteros(0),                       
150     f2dHistEffMassLambdaVsEffMassXiMinus(0), f2dHistEffMassXiVsEffMassOmegaMinus(0),
151     f2dHistEffMassLambdaVsEffMassXiPlus(0), f2dHistEffMassXiVsEffMassOmegaPlus(0),
152     f2dHistXiRadiusVsEffMassXiMinus(0), f2dHistXiRadiusVsEffMassXiPlus(0)
153
154
155
156 {
157   // Constructor
158
159   // Define input and output slots here
160   // Input slot #0 works with a TChain
161   
162   // Output slot #0 writes into a TList container (Cascade)
163   DefineOutput(1, TList::Class());
164 }
165
166
167
168
169
170
171 //________________________________________________________________________
172 void AliAnalysisTaskCheckCascade::UserCreateOutputObjects()
173 {
174   // Create histograms
175   // Called once
176
177
178
179  fListHistCascade = new TList();
180
181   
182         // - General histos
183   
184 if(! fHistTrackMultiplicity) {  
185         if(fCollidingSystems)// AA collisions   
186                 fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", 
187                         "Multiplicity distribution;Number of tracks;Events", 
188                         200, 0, 40000);                 
189         else // pp collisions
190                 fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", 
191                         "Track Multiplicity;Nbr of tracks/Evt;Events", 
192                         200, 0, 200);
193         fListHistCascade->Add(fHistTrackMultiplicity);
194 }
195
196 if(! fHistCascadeMultiplicity) {
197         if(fCollidingSystems)// AA collisions
198                 fHistCascadeMultiplicity = new TH1F("fHistCascadeMultiplicity", 
199                         "Multiplicity distribution;Number of Cascades;Events", 
200                         25, 0, 25);             
201         else // pp collisions
202                 fHistCascadeMultiplicity = new TH1F("fHistCascadeMultiplicity", 
203                         "Cascades per event;Nbr of Cascades/Evt;Events", 
204                         10, 0, 10);
205         fListHistCascade->Add(fHistCascadeMultiplicity);
206 }
207
208
209
210 if(! fHistVtxStatus ){
211         fHistVtxStatus   = new TH1F( "fHistVtxStatus" , "Does a Trckg Prim.vtx exist ?; true=1 or false=0; Nb of Events" , 4, -1.0, 3.0 );
212         fListHistCascade->Add(fHistVtxStatus);
213 }
214
215
216
217
218
219
220         // - Vertex Positions
221   
222 if(! fHistPosTrkgPrimaryVtxX ){
223         fHistPosTrkgPrimaryVtxX   = new TH1F( "fHistPosTrkgPrimaryVtxX" , "Trkg Prim. Vertex Position in x; x (cm); Events" , 200, -0.5, 0.5 );
224         fListHistCascade->Add(fHistPosTrkgPrimaryVtxX);
225 }
226
227
228 if(! fHistPosTrkgPrimaryVtxY){
229         fHistPosTrkgPrimaryVtxY   = new TH1F( "fHistPosTrkgPrimaryVtxY" , "Trkg Prim. Vertex Position in y; y (cm); Events" , 200, -0.5, 0.5 );
230         fListHistCascade->Add(fHistPosTrkgPrimaryVtxY);
231 }
232
233 if(! fHistPosTrkgPrimaryVtxZ ){
234         fHistPosTrkgPrimaryVtxZ   = new TH1F( "fHistPosTrkgPrimaryVtxZ" , "Trkg Prim. Vertex Position in z; z (cm); Events" , 100, -15.0, 15.0 );
235         fListHistCascade->Add(fHistPosTrkgPrimaryVtxZ);
236 }
237
238 if(! fHistTrkgPrimaryVtxRadius ){
239         fHistTrkgPrimaryVtxRadius  = new TH1F( "fHistTrkgPrimaryVtxRadius",  "Trkg Prim. Vertex radius; r (cm); Events" , 150, 0., 15.0 );
240         fListHistCascade->Add(fHistTrkgPrimaryVtxRadius);
241 }
242
243
244
245
246 if(! fHistPosBestPrimaryVtxX ){
247         fHistPosBestPrimaryVtxX   = new TH1F( "fHistPosBestPrimaryVtxX" , "Best Prim. Vertex Position in x; x (cm); Events" , 200, -0.5, 0.5 );
248         fListHistCascade->Add(fHistPosBestPrimaryVtxX);
249 }
250
251 if(! fHistPosBestPrimaryVtxY){
252         fHistPosBestPrimaryVtxY   = new TH1F( "fHistPosBestPrimaryVtxY" , "Best Prim. Vertex Position in y; y (cm); Events" , 200, -0.5, 0.5 );
253         fListHistCascade->Add(fHistPosBestPrimaryVtxY);
254 }
255
256 if(! fHistPosBestPrimaryVtxZ ){
257         fHistPosBestPrimaryVtxZ   = new TH1F( "fHistPosBestPrimaryVtxZ" , "Best Prim. Vertex Position in z; z (cm); Events" , 100, -15.0, 15.0 );
258         fListHistCascade->Add(fHistPosBestPrimaryVtxZ);
259 }
260
261 if(! fHistBestPrimaryVtxRadius ){
262         fHistBestPrimaryVtxRadius  = new TH1F( "fHistBestPrimaryVtxRadius",  "Best Prim.  vertex radius; r (cm); Events" , 150, 0., 15.0 );
263         fListHistCascade->Add(fHistBestPrimaryVtxRadius);
264 }
265
266 if(! f2dHistTrkgPrimVtxVsBestPrimVtx) {
267         f2dHistTrkgPrimVtxVsBestPrimVtx = new TH2F( "f2dHistTrkgPrimVtxVsBestPrimVtx", "r_{Trck Prim. Vtx} Vs r_{Best Prim. Vtx}; r_{Track Vtx} (cm); r_{Best Vtx} (cm)", 300, 0., 15.0, 300, 0., 15.);
268         fListHistCascade->Add(f2dHistTrkgPrimVtxVsBestPrimVtx);
269 }
270
271
272
273
274 // - Typical histos for cascades
275
276
277 if(! fHistEffMassXi) {
278      fHistEffMassXi = new TH1F("fHistEffMassXi", "Cascade candidates ; Invariant Mass (GeV/c^{2}) ; Counts", 200, 1.2, 2.0);
279      fListHistCascade->Add(fHistEffMassXi);
280 }
281    
282 if(! fHistChi2Xi ){
283         fHistChi2Xi = new TH1F("fHistChi2Xi", "Cascade #chi^{2}; #chi^{2}; Number of Cascades", 160, 0, 40);
284         fListHistCascade->Add(fHistChi2Xi);
285 }
286   
287 if(! fHistDcaXiDaughters ){
288         fHistDcaXiDaughters = new TH1F( "fHistDcaXiDaughters",  "DCA between Xi Daughters; DCA (cm) ; Number of Cascades", 100, 0., 0.5);
289         fListHistCascade->Add(fHistDcaXiDaughters);
290 }
291
292 if(! fHistDcaBachToPrimVertex) {
293         fHistDcaBachToPrimVertex = new TH1F("fHistDcaBachToPrimVertex", "DCA of Bach. to Prim. Vertex;DCA (cm);Number of Cascades", 250, 0., 0.25);
294         fListHistCascade->Add(fHistDcaBachToPrimVertex);
295 }
296
297 if(! fHistXiCosineOfPointingAngle) {
298         fHistXiCosineOfPointingAngle = new TH1F("fHistXiCosineOfPointingAngle", "Cosine of Xi Pointing Angle; Cos (Xi Point.Angl);Number of Xis", 200, 0.98, 1.0);
299         fListHistCascade->Add(fHistXiCosineOfPointingAngle);
300 }
301
302 if(! fHistXiRadius ){
303         fHistXiRadius  = new TH1F( "fHistXiRadius",  "Casc. decay transv. radius; r (cm); Counts" , 200, 0., 20.0 );
304         fListHistCascade->Add(fHistXiRadius);
305 }
306
307
308 // - Histos about ~ the "V0 part" of the cascade,  coming by inheritance from AliESDv0
309
310
311
312 if (! fHistMassLambdaAsCascDghter) {
313      fHistMassLambdaAsCascDghter = new TH1F("fHistMassLambdaAsCascDghter","#Lambda associated to Casc. candidates;Eff. Mass (GeV/c^{2});Counts", 160,1.00,1.8);
314     fListHistCascade->Add(fHistMassLambdaAsCascDghter);
315 }
316
317 if (! fHistV0Chi2Xi) {
318         fHistV0Chi2Xi = new TH1F("fHistV0Chi2Xi", "V0 #chi^{2}, in cascade; #chi^{2};Counts", 160, 0, 40);
319         fListHistCascade->Add(fHistV0Chi2Xi);
320 }
321
322 if (! fHistDcaV0DaughtersXi) {
323         fHistDcaV0DaughtersXi = new TH1F("fHistDcaV0DaughtersXi", "DCA between V0 daughters, in cascade;DCA (cm);Number of V0s", 120, 0., 0.6);
324         fListHistCascade->Add(fHistDcaV0DaughtersXi);
325 }
326
327 if (! fHistDcaV0ToPrimVertexXi) {
328         fHistDcaV0ToPrimVertexXi = new TH1F("fHistDcaV0ToPrimVertexXi", "DCA of V0 to Prim. Vertex, in cascade;DCA (cm);Number of Cascades", 200, 0., 1.);
329         fListHistCascade->Add(fHistDcaV0ToPrimVertexXi);
330 }
331
332 if (! fHistV0CosineOfPointingAngleXi) {
333         fHistV0CosineOfPointingAngleXi = new TH1F("fHistV0CosineOfPointingAngleXi", "Cosine of V0 Pointing Angle, in cascade;Cos(V0 Point. Angl); Counts", 200, 0.98, 1.0);
334         fListHistCascade->Add(fHistV0CosineOfPointingAngleXi);
335 }
336
337 if (! fHistV0RadiusXi) {
338         fHistV0RadiusXi = new TH1F("fHistV0RadiusXi", "V0 decay radius, in cascade; radius (cm); Counts", 200, 0, 20);
339         fListHistCascade->Add(fHistV0RadiusXi);
340 }
341
342 if (! fHistDcaPosToPrimVertexXi) {
343         fHistDcaPosToPrimVertexXi = new TH1F("fHistDcaPosToPrimVertexXi", "DCA of V0 pos daughter to Prim. Vertex;DCA (cm);Counts", 300, 0, 3);
344         fListHistCascade->Add(fHistDcaPosToPrimVertexXi);
345 }
346
347 if (! fHistDcaNegToPrimVertexXi) {
348         fHistDcaNegToPrimVertexXi = new TH1F("fHistDcaNegToPrimVertexXi", "DCA of V0 neg daughter to Prim. Vertex;DCA (cm);Counts", 300, 0, 3);
349         fListHistCascade->Add(fHistDcaNegToPrimVertexXi);
350 }
351
352
353
354
355         // - Effective mass histos for cascades.
356   
357 if (! fHistMassXiMinus) {
358     fHistMassXiMinus = new TH1F("fHistMassXiMinus","#Xi^{-} candidates;M( #Lambda , #pi^{-} ) (GeV/c^{2});Counts", 200,1.2,2.0);
359     fListHistCascade->Add(fHistMassXiMinus);
360 }
361   
362 if (! fHistMassXiPlus) {
363     fHistMassXiPlus = new TH1F("fHistMassXiPlus","#Xi^{+} candidates;M( #bar{#Lambda}^{0} , #pi^{+} ) (GeV/c^{2});Counts",200,1.2,2.0);
364     fListHistCascade->Add(fHistMassXiPlus);
365 }
366
367 if (! fHistMassOmegaMinus) {
368         fHistMassOmegaMinus = new TH1F("fHistMassOmegaMinus","#Omega^{-} candidates;M( #Lambda , K^{-} ) (GeV/c^{2});Counts", 250,1.5,2.5);
369     fListHistCascade->Add(fHistMassOmegaMinus);
370 }
371  
372 if (! fHistMassOmegaPlus) {
373         fHistMassOmegaPlus = new TH1F("fHistMassOmegaPlus","#Omega^{+} candidates;M( #bar{#Lambda}^{0} , K^{+} ) (GeV/c^{2});Counts", 250,1.5,2.5);
374     fListHistCascade->Add(fHistMassOmegaPlus);
375 }
376
377
378
379         // - Complements for QA
380
381 if(! fHistXiTransvMom ){
382         fHistXiTransvMom  = new TH1F( "fHistXiTransvMom" , "Xi transverse momentum ; p_{t}(#Xi) (GeV/c); Counts", 100, 0.0, 10.0);
383         fListHistCascade->Add(fHistXiTransvMom);
384 }
385
386 if(! fHistXiTotMom ){
387         fHistXiTotMom  = new TH1F( "fHistXiTotMom" , "Xi momentum norm; p_{tot}(#Xi) (GeV/c); Counts", 150, 0.0, 15.0);
388         fListHistCascade->Add(fHistXiTotMom);
389 }
390
391
392 if(! fHistBachTransvMom ){
393         fHistBachTransvMom  = new TH1F( "fHistBachTransvMom" , "Bach. transverse momentum ; p_{t}(Bach.) (GeV/c); Counts", 100, 0.0, 5.0);
394         fListHistCascade->Add(fHistBachTransvMom);
395 }
396
397 if(! fHistBachTotMom ){
398         fHistBachTotMom  = new TH1F( "fHistBachTotMom" , "Bach. momentum norm; p_{tot}(Bach.) (GeV/c); Counts", 100, 0.0, 5.0);
399         fListHistCascade->Add(fHistBachTotMom);
400 }
401
402
403 if(! fHistChargeXi ){
404         fHistChargeXi  = new TH1F( "fHistChargeXi" , "Charge of casc. candidates ; Sign ; Counts", 5, -2.0, 3.0);
405         fListHistCascade->Add(fHistChargeXi);
406 }
407
408
409 if (! fHistV0toXiCosineOfPointingAngle) {
410         fHistV0toXiCosineOfPointingAngle = new TH1F("fHistV0toXiCosineOfPointingAngle", "Cos. of V0 Ptng Angl / Xi vtx ;Cos(V0 Point. Angl / Xi vtx); Counts", 100, 0.99, 1.0);
411         fListHistCascade->Add(fHistV0toXiCosineOfPointingAngle);
412 }
413
414
415 if(! fHistRapXi ){
416         fHistRapXi  = new TH1F( "fHistRapXi" , "Rapidity of Xi candidates ; y ; Counts", 200, -5.0, 5.0);
417         fListHistCascade->Add(fHistRapXi);
418 }
419
420 if(! fHistRapOmega ){
421         fHistRapOmega  = new TH1F( "fHistRapOmega" , "Rapidity of Omega candidates ; y ; Counts", 200, -5.0, 5.0);
422         fListHistCascade->Add(fHistRapOmega);
423 }
424
425 if(! fHistEta ){
426         fHistEta  = new TH1F( "fHistEta" , "Pseudo-rap. of casc. candidates ; #eta ; Counts", 120, -3.0, 3.0);
427         fListHistCascade->Add(fHistEta);
428 }
429
430 if(! fHistTheta ){
431         fHistTheta  = new TH1F( "fHistTheta" , "#theta of casc. candidates ; #theta (deg) ; Counts", 180, 0., 180.0);
432         fListHistCascade->Add(fHistTheta);
433 }
434
435 if(! fHistPhi ){
436         fHistPhi  = new TH1F( "fHistPhi" , "#phi of casc. candidates ; #phi (deg) ; Counts", 360, 0., 360.);
437         fListHistCascade->Add(fHistPhi);
438 }
439
440
441 if(! f2dHistArmenteros) {
442         f2dHistArmenteros = new TH2F( "f2dHistArmenteros", "#alpha_{Arm}(casc. cand.) Vs Pt_{Arm}(casc. cand.); #alpha_{Arm} ; Pt_{Arm} (GeV/c)", 140, -1.2, 1.2, 300, 0., 0.3);
443         fListHistCascade->Add(f2dHistArmenteros);
444 }
445
446 if(! f2dHistEffMassLambdaVsEffMassXiMinus) {
447         f2dHistEffMassLambdaVsEffMassXiMinus = new TH2F( "f2dHistEffMassLambdaVsEffMassXiMinus", "M_{#Lambda} Vs M_{#Xi^{-} candidates} ; Inv. M_{#Lambda^{0}} (GeV/c^{2}) ; M( #Lambda , #pi^{-} ) (GeV/c^{2})", 300, 1.1,1.13, 200, 1.2, 2.0);
448         fListHistCascade->Add(f2dHistEffMassLambdaVsEffMassXiMinus);
449 }
450
451 if(! f2dHistEffMassXiVsEffMassOmegaMinus) {
452         f2dHistEffMassXiVsEffMassOmegaMinus = new TH2F( "f2dHistEffMassXiVsEffMassOmegaMinus", "M_{#Xi^{-} candidates} Vs M_{#Omega^{-} candidates} ; M( #Lambda , #pi^{-} ) (GeV/c^{2}) ; M( #Lambda , K^{-} ) (GeV/c^{2})", 200, 1.2, 2.0, 250, 1.5, 2.5);
453         fListHistCascade->Add(f2dHistEffMassXiVsEffMassOmegaMinus);
454 }
455
456 if(! f2dHistEffMassLambdaVsEffMassXiPlus) {
457         f2dHistEffMassLambdaVsEffMassXiPlus = new TH2F( "f2dHistEffMassLambdaVsEffMassXiPlus", "M_{#Lambda} Vs M_{#Xi^{+} candidates} ; Inv. M_{#Lambda^{0}} (GeV/c^{2}) ; M( #Lambda , #pi^{+} ) (GeV/c^{2})", 300, 1.1,1.13, 200, 1.2, 2.0);
458         fListHistCascade->Add(f2dHistEffMassLambdaVsEffMassXiPlus);
459 }
460
461 if(! f2dHistEffMassXiVsEffMassOmegaPlus) {
462         f2dHistEffMassXiVsEffMassOmegaPlus = new TH2F( "f2dHistEffMassXiVsEffMassOmegaPlus", "M_{#Xi^{+} candidates} Vs M_{#Omega^{+} candidates} ; M( #Lambda , #pi^{+} ) (GeV/c^{2}) ; M( #Lambda , K^{+} ) (GeV/c^{2})", 200, 1.2, 2.0, 250, 1.5, 2.5);
463         fListHistCascade->Add(f2dHistEffMassXiVsEffMassOmegaPlus);
464 }
465
466 if(! f2dHistXiRadiusVsEffMassXiMinus) {
467         f2dHistXiRadiusVsEffMassXiMinus = new TH2F( "f2dHistXiRadiusVsEffMassXiMinus", "Transv. R_{Xi Decay} Vs M_{#Xi^{-} candidates}; r_{Xi} (cm); M( #Lambda , #pi^{-} ) (GeV/c^{2}) ", 450, 0., 45.0, 200, 1.2, 2.0);
468         fListHistCascade->Add(f2dHistXiRadiusVsEffMassXiMinus);
469 }
470
471 if(! f2dHistXiRadiusVsEffMassXiPlus) {
472         f2dHistXiRadiusVsEffMassXiPlus = new TH2F( "f2dHistXiRadiusVsEffMassXiPlus", "Transv. R_{Xi Decay} Vs M_{#Xi^{+} candidates}; r_{Xi} (cm); M( #Lambda , #pi^{+} ) (GeV/c^{2}) ", 450, 0., 45.0, 200, 1.2, 2.0);
473         fListHistCascade->Add(f2dHistXiRadiusVsEffMassXiPlus);
474 }
475
476
477
478
479 }// end UserCreateOutputObjects
480
481
482
483
484
485
486 //________________________________________________________________________
487 void AliAnalysisTaskCheckCascade::UserExec(Option_t *) 
488 {
489   // Main loop
490   // Called for each event
491         
492         AliESDEvent *lESDevent = 0x0;
493         AliAODEvent *lAODevent = 0x0;
494         Int_t ncascades = -1;
495         
496         
497   // Connect to the InputEvent  
498   // After these lines, we should have an ESD/AOD event + the number of cascades in it.
499                 
500   if(fAnalysisType == "ESD"){
501         lESDevent = dynamic_cast<AliESDEvent*>( InputEvent() );
502         if (!lESDevent) {
503                 Printf("ERROR: lESDevent not available \n");
504                 return;
505         }
506         ncascades = lESDevent->GetNumberOfCascades();
507   }
508   
509   if(fAnalysisType == "AOD"){  
510           lAODevent = dynamic_cast<AliAODEvent*>( InputEvent() ); 
511         if (!lAODevent) {
512                 Printf("ERROR: lAODevent not available \n");
513                 return;
514         }
515         ncascades = lAODevent->GetNumberOfCascades();
516   }
517         
518   // ---------------------------------------------------------------
519   // I - Initialisation of the local variables that will be needed
520
521   
522         // - 0th part of initialisation : around primary vertex ...
523   Short_t  lStatusTrackingPrimVtx  = -2;
524   Double_t lTrkgPrimaryVtxPos[3]   = {-100.0, -100.0, -100.0};
525   Double_t lTrkgPrimaryVtxRadius3D = -500.0;
526   
527   Double_t lBestPrimaryVtxPos[3]   = {-100.0, -100.0, -100.0};
528   Double_t lBestPrimaryVtxRadius3D = -500.0;
529   
530         // - 1st part of initialisation : variables needed to store AliESDCascade data members
531   Double_t lEffMassXi      = 0. ;
532   Double_t lChi2Xi         = 0. ;
533   Double_t lDcaXiDaughters = 0. ;
534   Double_t lXiCosineOfPointingAngle = 0. ;
535   Double_t lPosXi[3] = { -1000.0, -1000.0, -1000.0 };
536   Double_t lXiRadius = 0. ;
537         
538                 
539         // - 2nd part of initialisation : about V0 part in cascades
540
541   Double_t lInvMassLambdaAsCascDghter = 0.;
542   Double_t lV0Chi2Xi         = 0. ;
543   Double_t lDcaV0DaughtersXi = 0.;
544         
545   Double_t lDcaBachToPrimVertexXi = 0., lDcaV0ToPrimVertexXi = 0.;
546   Double_t lDcaPosToPrimVertexXi  = 0.;
547   Double_t lDcaNegToPrimVertexXi  = 0.;
548   Double_t lV0CosineOfPointingAngleXi = 0. ;
549   Double_t lPosV0Xi[3] = { -1000. , -1000., -1000. }; // Position of VO coming from cascade
550   Double_t lV0RadiusXi = -1000.0;
551   Double_t lV0quality  = 0.;
552
553         
554         // - 3rd part of initialisation : Effective masses
555   
556   Double_t lInvMassXiMinus    = 0.;
557   Double_t lInvMassXiPlus     = 0.;
558   Double_t lInvMassOmegaMinus = 0.;
559   Double_t lInvMassOmegaPlus  = 0.;
560
561
562         // - 4th part of initialisation : extra info for QA
563   
564   Double_t lXiMomX       = 0., lXiMomY = 0., lXiMomZ = 0.;
565   Double_t lXiTransvMom  = 0. ;
566   Double_t lXiTotMom     = 0. ;
567         
568   Double_t lBachMomX       = 0., lBachMomY  = 0., lBachMomZ   = 0.;
569   Double_t lBachTransvMom  = 0.;
570   Double_t lBachTotMom     = 0.;
571
572   Short_t  lChargeXi = -2;
573   Double_t lV0toXiCosineOfPointingAngle = 0. ;
574
575   Double_t lRapXi   = -20.0, lRapOmega = -20.0,  lEta = -20.0, lTheta = 360., lPhi = 720. ;
576   Double_t lAlphaXi = -200., lPtArmXi  = -200.0;
577
578
579   
580   //-------------------------------------------------
581   // O - Cascade vertexer (ESD)
582
583         //   if(fAnalysisType == "ESD" ){
584         //      lESDevent->ResetCascades();
585         //      AliCascadeVertexer CascVtxer;
586         //      CascVtxer.V0sTracks2CascadeVertices(lESDevent);
587         //     }
588   
589   
590   // ---------------------------------------------------------------
591   // I - General histos (filled for any event) - (ESD)
592
593   fHistTrackMultiplicity  ->Fill( (InputEvent())->GetNumberOfTracks() );
594   fHistCascadeMultiplicity->Fill( ncascades );
595   
596   
597   
598   
599   for (Int_t iXi = 0; iXi < ncascades; iXi++)
600   {// This is the begining of the Cascade loop (ESD or AOD)
601            
602           
603   if(fAnalysisType == "ESD"){ 
604   
605   // -------------------------------------
606   // II - Calcultaion Part dedicated to Xi vertices (ESD)
607   
608         AliESDcascade *xi = lESDevent->GetCascade(iXi);
609         if (!xi) continue;
610
611         // Just to know which file is currently open : locate the file containing Xi
612         // cout << "Name of the file containing Xi candidate(s) :" 
613         //      << fInputHandler->GetTree()->GetCurrentFile()->GetName() 
614         //      << endl;
615         
616
617                 // - II.Step 1 : Characteristics of the event : prim. Vtx + magnetic field (ESD)
618                 //-------------
619         
620         // For new code (v4-16-Release-Rev06 or trunk)
621         
622         const AliESDVertex *lPrimaryTrackingVtx = lESDevent->GetPrimaryVertexTracks();  
623                 // get the vtx stored in ESD found with tracks
624         
625                 lPrimaryTrackingVtx->GetXYZ( lTrkgPrimaryVtxPos );
626                 lTrkgPrimaryVtxRadius3D = TMath::Sqrt(  lTrkgPrimaryVtxPos[0] * lTrkgPrimaryVtxPos[0] +
627                                                         lTrkgPrimaryVtxPos[1] * lTrkgPrimaryVtxPos[1] +
628                                                         lTrkgPrimaryVtxPos[2] * lTrkgPrimaryVtxPos[2] );
629
630                 lStatusTrackingPrimVtx = lPrimaryTrackingVtx->GetStatus();
631
632
633         const AliESDVertex *lPrimaryBestVtx = lESDevent->GetPrimaryVertex();    
634                 // get the best primary vertex available for the event
635                 // As done in AliCascadeVertexer, we keep the one which is the best one available.
636                 // between : Tracking vertex > SPD vertex > TPC vertex > default SPD vertex
637                 // This one will be used for next calculations (DCA essentially)
638         
639                 lPrimaryBestVtx->GetXYZ( lBestPrimaryVtxPos );
640                 lBestPrimaryVtxRadius3D = TMath::Sqrt(  lBestPrimaryVtxPos[0] * lBestPrimaryVtxPos[0] +
641                                                         lBestPrimaryVtxPos[1] * lBestPrimaryVtxPos[1] +
642                                                         lBestPrimaryVtxPos[2] * lBestPrimaryVtxPos[2] );
643                 
644         
645         // For older evts
646         
647         //      const AliESDVertex *lPrimaryTrackingVtx = lESDevent->GetPrimaryVertexTracks();  
648         //              get the vtx stored in ESD found with tracks
649         //      Double_t lTrkgPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
650         //              lPrimaryTrackingVtx->GetXYZ( lTrkgPrimaryVtxPos );
651         //              Double_t lTrkgPrimaryVtxRadius3D = TMath::Sqrt( lTrkgPrimaryVtxPos[0]*lTrkgPrimaryVtxPos[0] +
652         //                                              lTrkgPrimaryVtxPos[1] * lTrkgPrimaryVtxPos[1] +
653         //                                              lTrkgPrimaryVtxPos[2] * lTrkgPrimaryVtxPos[2] );
654         // 
655         //      const AliESDVertex *lPrimarySPDVtx = lESDevent->GetPrimaryVertexSPD();  
656         //              get the vtx found by exclusive use of SPD
657         //      Double_t lSPDPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
658         //              lPrimarySPDVtx->GetXYZ( lSPDPrimaryVtxPos );
659         //              
660         //      // As done in AliCascadeVertexer, we keep, between both retrieved vertices, 
661         //      // the one which is the best one available.
662         //      // This one will be used for next calculations (DCA essentially)
663         //              
664         //              Double_t lBestPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
665         //              if( lPrimaryTrackingVtx->GetStatus() ) { // if tracking vtx = ok
666         //                      lBestPrimaryVtxPos[0] = lTrkgPrimaryVtxPos[0];
667         //                      lBestPrimaryVtxPos[1] = lTrkgPrimaryVtxPos[1];
668         //                      lBestPrimaryVtxPos[2] = lTrkgPrimaryVtxPos[2];
669         //              }
670         //              else{
671         //                      lBestPrimaryVtxPos[0] = lSPDPrimaryVtxPos[0];
672         //                      lBestPrimaryVtxPos[1] = lSPDPrimaryVtxPos[1];
673         //                      lBestPrimaryVtxPos[2] = lSPDPrimaryVtxPos[2];
674         //              }
675         //
676         //      Double_t lBestPrimaryVtxRadius3D = TMath::Sqrt( lBestPrimaryVtxPos[0]*lBestPrimaryVtxPos[0] +
677         //                                              lBestPrimaryVtxPos[1] * lBestPrimaryVtxPos[1] +
678         //                                              lBestPrimaryVtxPos[2] * lBestPrimaryVtxPos[2] );
679         
680         // Back to normal
681                         
682         Double_t lMagneticField = lESDevent->GetMagneticField( );
683         
684         
685         
686                 // - II.Step 2 : Assigning the necessary variables for specific AliESDcascade data members (ESD)        
687                 //-------------
688         lV0quality = 0.;
689         xi->ChangeMassHypothesis(lV0quality , 3312); // default working hypothesis : cascade = Xi- decay
690
691         lEffMassXi                      = xi->GetEffMassXi();
692         lChi2Xi                         = xi->GetChi2Xi();
693         lDcaXiDaughters                 = xi->GetDcaXiDaughters();
694         lXiCosineOfPointingAngle        = xi->GetCascadeCosineOfPointingAngle( lBestPrimaryVtxPos[0],
695                                                                                 lBestPrimaryVtxPos[1],
696                                                                                 lBestPrimaryVtxPos[2] );
697                 // Take care : the best available vertex should be used (like in AliCascadeVertexer)
698         
699         xi->GetXYZcascade( lPosXi[0],  lPosXi[1], lPosXi[2] ); 
700         lXiRadius                       = TMath::Sqrt( lPosXi[0]*lPosXi[0]  +  lPosXi[1]*lPosXi[1] );
701                 
702                 
703
704                 // - II.Step 3 : around the tracks : Bach + V0 (ESD)
705                 // ~ Necessary variables for ESDcascade data members coming from the ESDv0 part (inheritance)
706                 //-------------
707                 
708                 UInt_t lIdxPosXi        = (UInt_t) TMath::Abs( xi->GetPindex() );
709                 UInt_t lIdxNegXi        = (UInt_t) TMath::Abs( xi->GetNindex() );
710                 UInt_t lBachIdx         = (UInt_t) TMath::Abs( xi->GetBindex() );
711                         // Care track label can be negative in MC production (linked with the track quality)
712                         // However = normally, not the case for track index ...
713
714         AliESDtrack *pTrackXi           = lESDevent->GetTrack( lIdxPosXi );
715         AliESDtrack *nTrackXi           = lESDevent->GetTrack( lIdxNegXi );
716         AliESDtrack *bachTrackXi        = lESDevent->GetTrack( lBachIdx );
717         if (!pTrackXi || !nTrackXi || !bachTrackXi ) {
718                 Printf("ERROR: Could not retrieve one of the 3 daughter tracks of the cascade ...");
719                 continue;
720         }
721         
722         lInvMassLambdaAsCascDghter      = xi->GetEffMass();
723                 // This value shouldn't change, whatever the working hyp. is : Xi-, Xi+, Omega-, Omega+
724         lDcaV0DaughtersXi               = xi->GetDcaV0Daughters(); 
725         lV0Chi2Xi                       = xi->GetChi2V0();
726         
727         lV0CosineOfPointingAngleXi      = xi->GetV0CosineOfPointingAngle( lBestPrimaryVtxPos[0],
728                                                                           lBestPrimaryVtxPos[1],
729                                                                           lBestPrimaryVtxPos[2] );
730
731         lDcaV0ToPrimVertexXi            = xi->GetD( lBestPrimaryVtxPos[0], 
732                                                     lBestPrimaryVtxPos[1], 
733                                                     lBestPrimaryVtxPos[2] );
734                 
735         lDcaBachToPrimVertexXi = TMath::Abs( bachTrackXi->GetD( lBestPrimaryVtxPos[0], 
736                                                                 lBestPrimaryVtxPos[1], 
737                                                                 lMagneticField  ) ); 
738                                         // Note : AliExternalTrackParam::GetD returns an algebraic value ...
739                 
740                 xi->GetXYZ( lPosV0Xi[0],  lPosV0Xi[1], lPosV0Xi[2] ); 
741         lV0RadiusXi             = TMath::Sqrt( lPosV0Xi[0]*lPosV0Xi[0]  +  lPosV0Xi[1]*lPosV0Xi[1] );
742         
743         lDcaPosToPrimVertexXi   = TMath::Abs( pTrackXi  ->GetD( lBestPrimaryVtxPos[0], 
744                                                                 lBestPrimaryVtxPos[1], 
745                                                                 lMagneticField  )     ); 
746         
747         lDcaNegToPrimVertexXi   = TMath::Abs( nTrackXi  ->GetD( lBestPrimaryVtxPos[0], 
748                                                                 lBestPrimaryVtxPos[1], 
749                                                                 lMagneticField  )     ); 
750         
751         
752                 // - II.Step 4 : around effective masses (ESD)
753                 // ~ change mass hypotheses to cover all the possibilities :  Xi-/+, Omega -/+
754                 //-------------
755
756         lV0quality = 0.;
757         xi->ChangeMassHypothesis(lV0quality , 3312);    // Calculate the effective mass of the Xi- candidate. 
758                                                         // pdg code 3312 = Xi-
759                 if( bachTrackXi->Charge() < 0 ) 
760                         lInvMassXiMinus = xi->GetEffMassXi();
761                         
762         lV0quality = 0.;
763         xi->ChangeMassHypothesis(lV0quality , -3312);   // Calculate the effective mass of the Xi+ candidate. 
764                                                         // pdg code -3312 = Xi+
765                 if( bachTrackXi->Charge() >  0 )
766                         lInvMassXiPlus = xi->GetEffMassXi();
767                         
768         lV0quality = 0.;
769         xi->ChangeMassHypothesis(lV0quality , 3334);    // Calculate the effective mass of the Xi- candidate. 
770                                                         // pdg code 3334 = Omega-
771                 if( bachTrackXi->Charge() < 0 )
772                         lInvMassOmegaMinus = xi->GetEffMassXi();
773                 
774         lV0quality = 0.;
775         xi->ChangeMassHypothesis(lV0quality , -3334);   // Calculate the effective mass of the Xi+ candidate. 
776                                                         // pdg code -3334  = Omega+
777                 if( bachTrackXi->Charge() >  0 )
778                         lInvMassOmegaPlus = xi->GetEffMassXi();
779
780         lV0quality = 0.;
781         xi->ChangeMassHypothesis(lV0quality , 3312);    // Back to default hyp.
782                                                         
783
784                 // - II.Step 5 : extra info for QA (ESD)
785                 // miscellaneous pieces onf info that may help regarding data quality assessment.
786                 //-------------
787
788         xi->GetPxPyPz( lXiMomX, lXiMomY, lXiMomZ );
789                 lXiTransvMom    = TMath::Sqrt( lXiMomX*lXiMomX   + lXiMomY*lXiMomY );
790                 lXiTotMom       = TMath::Sqrt( lXiMomX*lXiMomX   + lXiMomY*lXiMomY   + lXiMomZ*lXiMomZ );
791                 
792         xi->GetBPxPyPz(  lBachMomX,  lBachMomY,  lBachMomZ );
793                 lBachTransvMom  = TMath::Sqrt( lBachMomX*lBachMomX   + lBachMomY*lBachMomY );
794                 lBachTotMom     = TMath::Sqrt( lBachMomX*lBachMomX   + lBachMomY*lBachMomY  +  lBachMomZ*lBachMomZ  );
795
796         lChargeXi = xi->Charge();
797
798         lV0toXiCosineOfPointingAngle = xi->GetV0CosineOfPointingAngle( lPosXi[0], lPosXi[1], lPosXi[2] );
799         
800         lRapXi    = xi->RapXi();
801         lRapOmega = xi->RapOmega();
802         lEta      = xi->Eta();
803         lTheta    = xi->Theta() *180.0/TMath::Pi();
804         lPhi      = xi->Phi()   *180.0/TMath::Pi();
805         lAlphaXi  = xi->AlphaXi();
806         lPtArmXi  = xi->PtArmXi();
807         
808         
809   }// end of ESD treatment
810   
811  
812   if(fAnalysisType == "AOD"){
813         
814         // -------------------------------------
815         // II - Calcultaion Part dedicated to Xi vertices (ESD)
816         
817         const AliAODcascade *xi = lAODevent->GetCascade(iXi);
818         if (!xi) continue;
819                 
820         // Just to know which file is currently open : locate the file containing Xi
821         // cout << "Name of the file containing Xi candidate(s) :" <<  fesdH->GetTree()->GetCurrentFile()->GetName() << endl;
822         
823         
824                 // - II.Step 1 : Characteristics of the event : prim. Vtx + magnetic field (AOD)
825                 //-------------
826         
827         lTrkgPrimaryVtxPos[0]   = -100.0;
828         lTrkgPrimaryVtxPos[1]   = -100.0;
829         lTrkgPrimaryVtxPos[2]   = -100.0;
830         lTrkgPrimaryVtxRadius3D = -500. ;
831         // We don't have the different prim. vertex at the AOD level -> nothing to do.
832
833         const AliAODVertex *lPrimaryBestVtx = lAODevent->GetPrimaryVertex();    
834         // get the best primary vertex available for the event
835         // We may keep the one which is the best one available = GetVertex(0)
836         // Pb with pile-up to expect
837         // This one will be used for next calculations (DCA essentially)
838
839         lPrimaryBestVtx->GetXYZ( lBestPrimaryVtxPos );
840         lBestPrimaryVtxRadius3D = TMath::Sqrt(  lBestPrimaryVtxPos[0] * lBestPrimaryVtxPos[0] +
841                                                 lBestPrimaryVtxPos[1] * lBestPrimaryVtxPos[1] +
842                                                 lBestPrimaryVtxPos[2] * lBestPrimaryVtxPos[2] );
843                 
844         
845                 // - II.Step 2 : Assigning the necessary variables for specific AliAODcascade data members (AOD)        
846                 //-------------
847         
848         lEffMassXi                      = xi->MassXi(); // default working hypothesis : cascade = Xi- decay
849         lChi2Xi                         = xi->Chi2Xi();
850         lDcaXiDaughters                 = xi->DcaXiDaughters();
851         lXiCosineOfPointingAngle        = xi->CosPointingAngleXi( lBestPrimaryVtxPos[0], 
852                                                                   lBestPrimaryVtxPos[1], 
853                                                                   lBestPrimaryVtxPos[2] );
854                                         // Take care : 
855                                         // the best available vertex should be used (like in AliCascadeVertexer)
856
857                 lPosXi[0] = xi->DecayVertexXiX();
858                 lPosXi[1] = xi->DecayVertexXiY();
859                 lPosXi[2] = xi->DecayVertexXiZ();
860         lXiRadius = TMath::Sqrt( lPosXi[0]*lPosXi[0]  +  lPosXi[1]*lPosXi[1] );
861                 
862
863                 // - II.Step 3 : around the tracks : Bach + V0 (AOD)
864                 // ~ Necessary variables for AODcascade data members coming from the AODv0 part (inheritance)
865                 //-------------
866                 
867         lChargeXi                       = xi->ChargeXi();
868         
869         if( lChargeXi < 0)      
870           lInvMassLambdaAsCascDghter    = xi->MassLambda();
871         else                    
872           lInvMassLambdaAsCascDghter    = xi->MassAntiLambda();
873
874         lDcaV0DaughtersXi               = xi->DcaV0Daughters(); 
875         lV0Chi2Xi                       = xi->Chi2V0();
876
877         lV0CosineOfPointingAngleXi      = xi->CosPointingAngle( lBestPrimaryVtxPos );
878         lDcaV0ToPrimVertexXi            = xi->DcaV0ToPrimVertex();
879         
880         lDcaBachToPrimVertexXi          = xi->DcaBachToPrimVertex(); 
881         
882         
883                 lPosV0Xi[0] = xi->DecayVertexV0X();
884                 lPosV0Xi[1] = xi->DecayVertexV0Y();
885                 lPosV0Xi[2] = xi->DecayVertexV0Z(); 
886         lV0RadiusXi     = TMath::Sqrt( lPosV0Xi[0]*lPosV0Xi[0]  +  lPosV0Xi[1]*lPosV0Xi[1] );
887
888         lDcaPosToPrimVertexXi           = xi->DcaPosToPrimVertex(); 
889         lDcaNegToPrimVertexXi           = xi->DcaNegToPrimVertex(); 
890
891
892                 // - II.Step 4 : around effective masses (AOD)
893                 // ~ change mass hypotheses to cover all the possibilities :  Xi-/+, Omega -/+
894                 //-------------
895
896         if( lChargeXi < 0 )             lInvMassXiMinus         = xi->MassXi();
897         if( lChargeXi > 0 )             lInvMassXiPlus          = xi->MassXi();
898         if( lChargeXi < 0 )             lInvMassOmegaMinus      = xi->MassOmega();
899         if( lChargeXi > 0 )             lInvMassOmegaPlus       = xi->MassOmega();
900
901         
902                 // - II.Step 5 : extra info for QA (AOD)
903                 // miscellaneous pieces onf info that may help regarding data quality assessment.
904                 //-------------
905
906                 lXiMomX = xi->MomXiX();
907                 lXiMomY = xi->MomXiY();
908                 lXiMomZ = xi->MomXiZ();
909         lXiTransvMom    = TMath::Sqrt( lXiMomX*lXiMomX   + lXiMomY*lXiMomY );
910         lXiTotMom       = TMath::Sqrt( lXiMomX*lXiMomX   + lXiMomY*lXiMomY   + lXiMomZ*lXiMomZ );
911         
912                 lBachMomX = xi->MomBachX();
913                 lBachMomY = xi->MomBachY();
914                 lBachMomZ = xi->MomBachZ();             
915         lBachTransvMom  = TMath::Sqrt( lBachMomX*lBachMomX   + lBachMomY*lBachMomY );
916         lBachTotMom     = TMath::Sqrt( lBachMomX*lBachMomX   + lBachMomY*lBachMomY  +  lBachMomZ*lBachMomZ  );
917
918         
919         lV0toXiCosineOfPointingAngle = xi->CosPointingAngle( xi->GetDecayVertexXi() );
920         
921         lRapXi    = xi->RapXi();
922         lRapOmega = xi->RapOmega();
923         lEta      = xi->Eta();                          // Will not work ! need a method Pz(), Py() Px() 
924         lTheta    = xi->Theta() *180.0/TMath::Pi();     // in AODcascade.
925         lPhi      = xi->Phi()   *180.0/TMath::Pi();     // Here, we will get eta, theta, phi for the V0 ...
926         lAlphaXi  = xi->AlphaXi();
927         lPtArmXi  = xi->PtArmXi();
928
929          
930   }// end of AOD treatment
931
932
933   // -------------------------------------
934   // III - Filling the TH1Fs
935   
936
937         // - III.Step 1 
938
939         fHistVtxStatus                  ->Fill( lStatusTrackingPrimVtx   );  // 1 if tracking vtx = ok
940
941         if( lStatusTrackingPrimVtx ){
942                 fHistPosTrkgPrimaryVtxX ->Fill( lTrkgPrimaryVtxPos[0]    );
943                 fHistPosTrkgPrimaryVtxY ->Fill( lTrkgPrimaryVtxPos[1]    );     
944                 fHistPosTrkgPrimaryVtxZ ->Fill( lTrkgPrimaryVtxPos[2]    ); 
945                 fHistTrkgPrimaryVtxRadius->Fill( lTrkgPrimaryVtxRadius3D );
946         }
947
948         fHistPosBestPrimaryVtxX         ->Fill( lBestPrimaryVtxPos[0]    );
949         fHistPosBestPrimaryVtxY         ->Fill( lBestPrimaryVtxPos[1]    );
950         fHistPosBestPrimaryVtxZ         ->Fill( lBestPrimaryVtxPos[2]    );
951         fHistBestPrimaryVtxRadius       ->Fill( lBestPrimaryVtxRadius3D  );
952         
953         f2dHistTrkgPrimVtxVsBestPrimVtx->Fill( lTrkgPrimaryVtxRadius3D, lBestPrimaryVtxRadius3D );
954
955                 
956         // - III.Step 2
957         fHistEffMassXi                  ->Fill( lEffMassXi               );
958         fHistChi2Xi                     ->Fill( lChi2Xi                  );     // Flag CascadeVtxer: Cut Variable a
959         fHistDcaXiDaughters             ->Fill( lDcaXiDaughters          );     // Flag CascadeVtxer: Cut Variable e 
960         fHistDcaBachToPrimVertex        ->Fill( lDcaBachToPrimVertexXi   );     // Flag CascadeVtxer: Cut Variable d
961         fHistXiCosineOfPointingAngle    ->Fill( lXiCosineOfPointingAngle );     // Flag CascadeVtxer: Cut Variable f
962         fHistXiRadius                   ->Fill( lXiRadius                );     // Flag CascadeVtxer: Cut Variable g+h
963         
964         
965         // - III.Step 3
966         fHistMassLambdaAsCascDghter     ->Fill( lInvMassLambdaAsCascDghter );   // Flag CascadeVtxer: Cut Variable c
967         fHistV0Chi2Xi                   ->Fill( lV0Chi2Xi                  );   
968         fHistDcaV0DaughtersXi           ->Fill( lDcaV0DaughtersXi          );
969         fHistV0CosineOfPointingAngleXi  ->Fill( lV0CosineOfPointingAngleXi ); 
970         fHistV0RadiusXi                 ->Fill( lV0RadiusXi                );
971         
972         fHistDcaV0ToPrimVertexXi        ->Fill( lDcaV0ToPrimVertexXi       );   // Flag CascadeVtxer: Cut Variable b
973         fHistDcaPosToPrimVertexXi       ->Fill( lDcaPosToPrimVertexXi      );
974         fHistDcaNegToPrimVertexXi       ->Fill( lDcaNegToPrimVertexXi      );
975                 
976         
977         // - III.Step 4
978         if( lChargeXi < 0 )     fHistMassXiMinus        ->Fill( lInvMassXiMinus );
979         if( lChargeXi > 0 )     fHistMassXiPlus         ->Fill( lInvMassXiPlus );
980         if( lChargeXi < 0 )     fHistMassOmegaMinus     ->Fill( lInvMassOmegaMinus );
981         if( lChargeXi > 0 )     fHistMassOmegaPlus      ->Fill( lInvMassOmegaPlus );
982
983 //      if( bachTrackXi->Charge() < 0 ) fHistMassXiMinus        ->Fill( lInvMassXiMinus    );
984 //      if( bachTrackXi->Charge() > 0 ) fHistMassXiPlus         ->Fill( lInvMassXiPlus     );
985 //      if( bachTrackXi->Charge() < 0 ) fHistMassOmegaMinus     ->Fill( lInvMassOmegaMinus );
986 //      if( bachTrackXi->Charge() > 0 ) fHistMassOmegaPlus      ->Fill( lInvMassOmegaPlus  );
987
988
989         // - III.Step 5
990         fHistXiTransvMom        ->Fill( lXiTransvMom   );
991         fHistXiTotMom           ->Fill( lXiTotMom      );
992                 
993         fHistBachTransvMom      ->Fill( lBachTransvMom );
994         fHistBachTotMom         ->Fill( lBachTotMom    );
995
996         fHistChargeXi           ->Fill( lChargeXi      );
997         fHistV0toXiCosineOfPointingAngle->Fill( lV0toXiCosineOfPointingAngle );
998
999         fHistRapXi              ->Fill( lRapXi               );
1000         fHistRapOmega           ->Fill( lRapOmega            );
1001         fHistEta                ->Fill( lEta                 );
1002         fHistTheta              ->Fill( lTheta               );
1003         fHistPhi                ->Fill( lPhi                 );
1004
1005         f2dHistArmenteros       ->Fill( lAlphaXi, lPtArmXi   );
1006                 
1007         if( lChargeXi < 0 ) {
1008                 f2dHistEffMassLambdaVsEffMassXiMinus->Fill( lInvMassLambdaAsCascDghter, lInvMassXiMinus ); 
1009                 f2dHistEffMassXiVsEffMassOmegaMinus ->Fill( lInvMassXiMinus, lInvMassOmegaMinus );
1010                 f2dHistXiRadiusVsEffMassXiMinus     ->Fill( lXiRadius, lInvMassXiMinus );
1011         }
1012         else{
1013                 f2dHistEffMassLambdaVsEffMassXiPlus ->Fill( lInvMassLambdaAsCascDghter, lInvMassXiPlus );
1014                 f2dHistEffMassXiVsEffMassOmegaPlus  ->Fill( lInvMassXiPlus, lInvMassOmegaPlus );
1015                 f2dHistXiRadiusVsEffMassXiPlus->Fill( lXiRadius, lInvMassXiPlus);
1016         }
1017                 
1018    
1019     }// end of the Cascade loop (ESD or AOD)
1020     
1021   
1022   // Post output data.
1023  PostData(1, fListHistCascade);
1024 }
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035 //________________________________________________________________________
1036 void AliAnalysisTaskCheckCascade::Terminate(Option_t *) 
1037 {
1038   // Draw result to the screen
1039   // Called once at the end of the query
1040
1041   fHistTrackMultiplicity = dynamic_cast<TH1F*> (   ((TList*)GetOutputData(1))->FindObject("fHistTrackMultiplicity") );
1042   if (!fHistTrackMultiplicity) {
1043     Printf("ERROR: fHistTrackMultiplicity not available");
1044     return;
1045   }
1046  
1047   fHistCascadeMultiplicity = dynamic_cast<TH1F*> ( ((TList*)GetOutputData(1))->FindObject("fHistCascadeMultiplicity"));
1048   if (!fHistCascadeMultiplicity) {
1049     Printf("ERROR: fHistCascadeMultiplicity not available");
1050     return;
1051   }
1052    
1053   TCanvas *c2 = new TCanvas("AliAnalysisTaskCheckCascade","Multiplicity",10,10,510,510);
1054   c2->cd(1)->SetLogy();
1055
1056   fHistTrackMultiplicity->SetMarkerStyle(22);
1057   fHistTrackMultiplicity->DrawCopy("E");
1058   
1059   fHistCascadeMultiplicity->SetMarkerStyle(26);
1060   fHistCascadeMultiplicity->DrawCopy("ESAME");
1061   
1062 }