]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/SPECTRA/AliAnalysisTaskESDCheckCascade.cxx
Compilation warning
[u/mrichter/AliRoot.git] / PWG2 / SPECTRA / AliAnalysisTaskESDCheckCascade.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 //                 AliAnalysisTaskESDCheckCascade class
17 //            This task is for QAing the Cascades from the ESD
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 AliESDv0;
35
36 #include <iostream>
37
38 #include "TChain.h"
39 #include "TFile.h"
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 #include "AliAnalysisManager.h"
48 #include "AliESDInputHandler.h"
49
50 #include "AliESDEvent.h"
51 //#include "AliCascadeVertexer.h"
52
53 #include "AliESDcascade.h"
54
55 #include "AliAnalysisTaskESDCheckCascade.h"
56
57 ClassImp(AliAnalysisTaskESDCheckCascade)
58
59
60
61 //________________________________________________________________________
62 AliAnalysisTaskESDCheckCascade::AliAnalysisTaskESDCheckCascade() 
63   : AliAnalysisTask(), 
64
65     fesdH(0), fMCevent(0), fESD(0),
66
67         // - Cascade part initialisation
68     fListHistCascade(0),
69     fHistTrackMultiplicity(0), fHistCascadeMultiplicity(0),
70     fHistVtxStatus(0),
71
72     fHistPosTrkgPrimaryVtxX(0), fHistPosTrkgPrimaryVtxY(0), fHistPosTrkgPrimaryVtxZ(0), fHistTrkgPrimaryVtxRadius(0),
73     fHistPosBestPrimaryVtxX(0), fHistPosBestPrimaryVtxY(0), fHistPosBestPrimaryVtxZ(0), fHistBestPrimaryVtxRadius(0),
74     f2dHistTrkgPrimVtxVsBestPrimVtx(0),
75
76     fHistEffMassXi(0),  fHistChi2Xi(0),  
77     fHistDcaXiDaughters(0), fHistDcaBachToPrimVertex(0), fHistXiCosineOfPointingAngle(0), fHistXiRadius(0),
78
79     fHistMassLambdaAsCascDghter(0),
80     fHistV0Chi2Xi(0),
81     fHistDcaV0DaughtersXi(0),
82     fHistDcaV0ToPrimVertexXi(0), 
83     fHistV0CosineOfPointingAngleXi(0),
84     fHistV0RadiusXi(0),
85     fHistDcaPosToPrimVertexXi(0), fHistDcaNegToPrimVertexXi(0), 
86
87     fHistMassXiMinus(0), fHistMassXiPlus(0),
88     fHistMassOmegaMinus(0), fHistMassOmegaPlus(0),
89
90     fHistXiTransvMom(0),    fHistXiTotMom(0),
91     fHistBachTransvMom(0),   fHistBachTotMom(0),
92
93     fHistChargeXi(0),
94     fHistV0toXiCosineOfPointingAngle(0),
95
96     fHistRapXi(0), fHistRapOmega(0), fHistEta(0),
97     fHistTheta(0), fHistPhi(0),
98
99     f2dHistArmenteros(0),                       
100     f2dHistEffMassXiVsEffMassLambda(0),
101     f2dHistXiRadiusVsEffMassXi(0)
102
103 {
104   // Dummy Constructor
105 }
106
107
108
109
110
111
112
113
114 //________________________________________________________________________
115 AliAnalysisTaskESDCheckCascade::AliAnalysisTaskESDCheckCascade(const char *name) 
116   : AliAnalysisTask(name, ""), 
117
118     fesdH(0), fMCevent(0), fESD(0),    
119     
120         // - Cascade part initialisation
121     fListHistCascade(0),
122     fHistTrackMultiplicity(0), fHistCascadeMultiplicity(0),
123     fHistVtxStatus(0),
124
125     fHistPosTrkgPrimaryVtxX(0), fHistPosTrkgPrimaryVtxY(0), fHistPosTrkgPrimaryVtxZ(0), fHistTrkgPrimaryVtxRadius(0),
126     fHistPosBestPrimaryVtxX(0), fHistPosBestPrimaryVtxY(0), fHistPosBestPrimaryVtxZ(0), fHistBestPrimaryVtxRadius(0),
127     f2dHistTrkgPrimVtxVsBestPrimVtx(0),
128
129     fHistEffMassXi(0),  fHistChi2Xi(0),  
130     fHistDcaXiDaughters(0), fHistDcaBachToPrimVertex(0), fHistXiCosineOfPointingAngle(0), fHistXiRadius(0),
131
132     fHistMassLambdaAsCascDghter(0),
133     fHistV0Chi2Xi(0),
134     fHistDcaV0DaughtersXi(0),
135     fHistDcaV0ToPrimVertexXi(0), 
136     fHistV0CosineOfPointingAngleXi(0),
137     fHistV0RadiusXi(0),
138     fHistDcaPosToPrimVertexXi(0), fHistDcaNegToPrimVertexXi(0), 
139
140     fHistMassXiMinus(0), fHistMassXiPlus(0),
141     fHistMassOmegaMinus(0), fHistMassOmegaPlus(0),
142
143     fHistXiTransvMom(0),    fHistXiTotMom(0),
144     fHistBachTransvMom(0),   fHistBachTotMom(0),
145
146     fHistChargeXi(0),
147     fHistV0toXiCosineOfPointingAngle(0),
148
149     fHistRapXi(0), fHistRapOmega(0), fHistEta(0),
150     fHistTheta(0), fHistPhi(0),
151
152     f2dHistArmenteros(0),                       
153     f2dHistEffMassXiVsEffMassLambda(0),
154     f2dHistXiRadiusVsEffMassXi(0)
155
156
157
158 {
159   // Constructor
160
161   // Define input and output slots here
162   // Input slot #0 works with a TChain
163   DefineInput(0, TChain::Class());
164   // Output slot #0 writes into a TList container (Cascade)
165   DefineOutput(0, TList::Class());
166 }
167
168
169
170
171
172
173
174 //________________________________________________________________________
175 void AliAnalysisTaskESDCheckCascade::ConnectInputData(Option_t *) 
176 {
177   // Connect ESD or AOD here
178   // Called once
179
180   AliLog::SetGlobalLogLevel(AliLog::kError); // to suppress the extensive info prompted by a run with MC
181
182   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
183   if (!tree) {
184     Printf("ERROR: Could not read chain from input slot 0");
185   } else {
186         
187         fesdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
188     //AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
189
190         if (!fesdH) {
191                 Printf("ERROR: Could not get ESDInputHandler");
192         } else
193                 fESD = fesdH->GetEvent();
194   }
195 }
196
197
198
199
200
201
202
203 //________________________________________________________________________
204 void AliAnalysisTaskESDCheckCascade::CreateOutputObjects()
205 {
206   // Create histograms
207   // Called once
208
209
210
211  fListHistCascade = new TList();
212
213   
214         // - General histos
215   
216 if(! fHistTrackMultiplicity) {
217         fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Track Multiplicity;Nbr of tracks/Evt;Events", 200, 0, 200);
218         //  fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 200, 0, 40000); //HERE for Pb-Pb
219         fListHistCascade->Add(fHistTrackMultiplicity);
220 }
221
222 if(! fHistCascadeMultiplicity) {
223         fHistCascadeMultiplicity = new TH1F("fHistCascadeMultiplicity", "Cascades per event;Nbr of Cascades/Evt;Events", 10, 0, 10);
224         //  fHistCascadeMultiplicity = new TH1F("fHistCascadeMultiplicity", "Multiplicity distribution;Number of Cascades;Events", 200, 0, 4000); //HERE
225         fListHistCascade->Add(fHistCascadeMultiplicity);
226 }
227
228
229
230 if(! fHistVtxStatus ){
231         fHistVtxStatus   = new TH1F( "fHistVtxStatus" , "Does a Trckg Prim.vtx exist ?; true=1 or false=0; Nb of Events" , 4, -1.0, 3.0 );
232         fListHistCascade->Add(fHistVtxStatus);
233 }
234
235
236
237
238
239
240         // - Vertex Positions
241   
242 if(! fHistPosTrkgPrimaryVtxX ){
243         fHistPosTrkgPrimaryVtxX   = new TH1F( "fHistPosTrkgPrimaryVtxX" , "Trkg Prim. Vertex Position in x; x (cm); Events" , 200, -0.5, 0.5 );
244         fListHistCascade->Add(fHistPosTrkgPrimaryVtxX);
245 }
246
247
248 if(! fHistPosTrkgPrimaryVtxY){
249         fHistPosTrkgPrimaryVtxY   = new TH1F( "fHistPosTrkgPrimaryVtxY" , "Trkg Prim. Vertex Position in y; y (cm); Events" , 200, -0.5, 0.5 );
250         fListHistCascade->Add(fHistPosTrkgPrimaryVtxY);
251 }
252
253 if(! fHistPosTrkgPrimaryVtxZ ){
254         fHistPosTrkgPrimaryVtxZ   = new TH1F( "fHistPosTrkgPrimaryVtxZ" , "Trkg Prim. Vertex Position in z; z (cm); Events" , 100, -15.0, 15.0 );
255         fListHistCascade->Add(fHistPosTrkgPrimaryVtxZ);
256 }
257
258 if(! fHistTrkgPrimaryVtxRadius ){
259         fHistTrkgPrimaryVtxRadius  = new TH1F( "fHistTrkgPrimaryVtxRadius",  "Trkg Prim. Vertex radius; r (cm); Events" , 150, 0., 15.0 );
260         fListHistCascade->Add(fHistTrkgPrimaryVtxRadius);
261 }
262
263
264
265
266 if(! fHistPosBestPrimaryVtxX ){
267         fHistPosBestPrimaryVtxX   = new TH1F( "fHistPosBestPrimaryVtxX" , "Best Prim. Vertex Position in x; x (cm); Events" , 200, -0.5, 0.5 );
268         fListHistCascade->Add(fHistPosBestPrimaryVtxX);
269 }
270
271 if(! fHistPosBestPrimaryVtxY){
272         fHistPosBestPrimaryVtxY   = new TH1F( "fHistPosBestPrimaryVtxY" , "Best Prim. Vertex Position in y; y (cm); Events" , 200, -0.5, 0.5 );
273         fListHistCascade->Add(fHistPosBestPrimaryVtxY);
274 }
275
276 if(! fHistPosBestPrimaryVtxZ ){
277         fHistPosBestPrimaryVtxZ   = new TH1F( "fHistPosBestPrimaryVtxZ" , "Best Prim. Vertex Position in z; z (cm); Events" , 100, -15.0, 15.0 );
278         fListHistCascade->Add(fHistPosBestPrimaryVtxZ);
279 }
280
281 if(! fHistBestPrimaryVtxRadius ){
282         fHistBestPrimaryVtxRadius  = new TH1F( "fHistBestPrimaryVtxRadius",  "Best Prim.  vertex radius; r (cm); Events" , 150, 0., 15.0 );
283         fListHistCascade->Add(fHistBestPrimaryVtxRadius);
284 }
285
286 if(! f2dHistTrkgPrimVtxVsBestPrimVtx) {
287         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.);
288         fListHistCascade->Add(f2dHistTrkgPrimVtxVsBestPrimVtx);
289 }
290
291
292
293
294 // - Typical histos for cascades
295
296
297 if(! fHistEffMassXi) {
298      fHistEffMassXi = new TH1F("fHistEffMassXi", "Cascade candidates ; Invariant Mass (GeV/c^{2}) ; Counts", 200, 1.2, 2.0);
299      fListHistCascade->Add(fHistEffMassXi);
300 }
301    
302 if(! fHistChi2Xi ){
303         fHistChi2Xi = new TH1F("fHistChi2Xi", "Cascade #chi^{2}; #chi^{2}; Number of Cascades", 160, 0, 40);
304         fListHistCascade->Add(fHistChi2Xi);
305 }
306   
307 if(! fHistDcaXiDaughters ){
308         fHistDcaXiDaughters = new TH1F( "fHistDcaXiDaughters",  "DCA between Xi Daughters; DCA (cm) ; Number of Cascades", 100, 0., 0.5);
309         fListHistCascade->Add(fHistDcaXiDaughters);
310 }
311
312 if(! fHistDcaBachToPrimVertex) {
313         fHistDcaBachToPrimVertex = new TH1F("fHistDcaBachToPrimVertex", "DCA of Bach. to Prim. Vertex;DCA (cm);Number of Cascades", 250, 0., 0.25);
314         fListHistCascade->Add(fHistDcaBachToPrimVertex);
315 }
316
317 if(! fHistXiCosineOfPointingAngle) {
318         fHistXiCosineOfPointingAngle = new TH1F("fHistXiCosineOfPointingAngle", "Cosine of Xi Pointing Angle; Cos (Xi Point.Angl);Number of Xis", 200, 0.98, 1.0);
319         fListHistCascade->Add(fHistXiCosineOfPointingAngle);
320 }
321
322 if(! fHistXiRadius ){
323         fHistXiRadius  = new TH1F( "fHistXiRadius",  "Casc. decay transv. radius; r (cm); Counts" , 200, 0., 20.0 );
324         fListHistCascade->Add(fHistXiRadius);
325 }
326
327
328 // - Histos about ~ the "V0 part" of the cascade,  coming by inheritance from AliESDv0
329
330
331
332 if (! fHistMassLambdaAsCascDghter) {
333      fHistMassLambdaAsCascDghter = new TH1F("fHistMassLambdaAsCascDghter","#Lambda associated to Casc. candidates;Eff. Mass (GeV/c^{2});Counts", 160,1.00,1.8);
334     fListHistCascade->Add(fHistMassLambdaAsCascDghter);
335 }
336
337 if (! fHistV0Chi2Xi) {
338         fHistV0Chi2Xi = new TH1F("fHistV0Chi2Xi", "V0 #chi^{2}, in cascade; #chi^{2};Counts", 160, 0, 40);
339         fListHistCascade->Add(fHistV0Chi2Xi);
340 }
341
342 if (! fHistDcaV0DaughtersXi) {
343         fHistDcaV0DaughtersXi = new TH1F("fHistDcaV0DaughtersXi", "DCA between V0 daughters, in cascade;DCA (cm);Number of V0s", 120, 0., 0.6);
344         fListHistCascade->Add(fHistDcaV0DaughtersXi);
345 }
346
347 if (! fHistDcaV0ToPrimVertexXi) {
348         fHistDcaV0ToPrimVertexXi = new TH1F("fHistDcaV0ToPrimVertexXi", "DCA of V0 to Prim. Vertex, in cascade;DCA (cm);Number of Cascades", 200, 0., 1.);
349         fListHistCascade->Add(fHistDcaV0ToPrimVertexXi);
350 }
351
352 if (! fHistV0CosineOfPointingAngleXi) {
353         fHistV0CosineOfPointingAngleXi = new TH1F("fHistV0CosineOfPointingAngleXi", "Cosine of V0 Pointing Angle, in cascade;Cos(V0 Point. Angl); Counts", 200, 0.98, 1.0);
354         fListHistCascade->Add(fHistV0CosineOfPointingAngleXi);
355 }
356
357 if (! fHistV0RadiusXi) {
358         fHistV0RadiusXi = new TH1F("fHistV0RadiusXi", "V0 decay radius, in cascade; radius (cm); Counts", 200, 0, 20);
359         fListHistCascade->Add(fHistV0RadiusXi);
360 }
361
362 if (! fHistDcaPosToPrimVertexXi) {
363         fHistDcaPosToPrimVertexXi = new TH1F("fHistDcaPosToPrimVertexXi", "DCA of V0 pos daughter to Prim. Vertex;DCA (cm);Counts", 300, 0, 3);
364         fListHistCascade->Add(fHistDcaPosToPrimVertexXi);
365 }
366
367 if (! fHistDcaNegToPrimVertexXi) {
368         fHistDcaNegToPrimVertexXi = new TH1F("fHistDcaNegToPrimVertexXi", "DCA of V0 neg daughter to Prim. Vertex;DCA (cm);Counts", 300, 0, 3);
369         fListHistCascade->Add(fHistDcaNegToPrimVertexXi);
370 }
371
372
373
374
375         // - Effective mass histos for cascades.
376   
377 if (! fHistMassXiMinus) {
378     fHistMassXiMinus = new TH1F("fHistMassXiMinus","#Xi^{-} candidates;M( #Lambda , #pi^{-} ) (GeV/c^{2});Counts", 200,1.2,2.0);
379     fListHistCascade->Add(fHistMassXiMinus);
380 }
381   
382 if (! fHistMassXiPlus) {
383     fHistMassXiPlus = new TH1F("fHistMassXiPlus","#Xi^{+} candidates;M( #bar{#Lambda}^{0} , #pi^{+} ) (GeV/c^{2});Counts",200,1.2,2.0);
384     fListHistCascade->Add(fHistMassXiPlus);
385 }
386
387 if (! fHistMassOmegaMinus) {
388     fHistMassOmegaMinus = new TH1F("fHistMassOmegaMinus","#Omega^{-} candidates;M( #Lambda , K^{-} ) (GeV/c^{2});Counts", 250,1.2,2.2);
389     fListHistCascade->Add(fHistMassOmegaMinus);
390 }
391  
392 if (! fHistMassOmegaPlus) {
393     fHistMassOmegaPlus = new TH1F("fHistMassOmegaPlus","#Omega^{+} candidates;M( #bar{#Lambda}^{0} , K^{+} ) (GeV/c^{2});Counts",250,1.2,2.2);
394     fListHistCascade->Add(fHistMassOmegaPlus);
395 }
396
397
398
399         // - Complements for QA
400
401 if(! fHistXiTransvMom ){
402         fHistXiTransvMom  = new TH1F( "fHistXiTransvMom" , "Xi transverse momentum ; p_{t}(#Xi) (GeV/c); Counts", 100, 0.0, 10.0);
403         fListHistCascade->Add(fHistXiTransvMom);
404 }
405
406 if(! fHistXiTotMom ){
407         fHistXiTotMom  = new TH1F( "fHistXiTotMom" , "Xi momentum norm; p_{tot}(#Xi) (GeV/c); Counts", 150, 0.0, 15.0);
408         fListHistCascade->Add(fHistXiTotMom);
409 }
410
411
412 if(! fHistBachTransvMom ){
413         fHistBachTransvMom  = new TH1F( "fHistBachTransvMom" , "Bach. transverse momentum ; p_{t}(Bach.) (GeV/c); Counts", 100, 0.0, 5.0);
414         fListHistCascade->Add(fHistBachTransvMom);
415 }
416
417 if(! fHistBachTotMom ){
418         fHistBachTotMom  = new TH1F( "fHistBachTotMom" , "Bach. momentum norm; p_{tot}(Bach.) (GeV/c); Counts", 100, 0.0, 5.0);
419         fListHistCascade->Add(fHistBachTotMom);
420 }
421
422
423 if(! fHistChargeXi ){
424         fHistChargeXi  = new TH1F( "fHistChargeXi" , "Charge of casc. candidates ; Sign ; Counts", 5, -2.0, 3.0);
425         fListHistCascade->Add(fHistChargeXi);
426 }
427
428
429 if (! fHistV0toXiCosineOfPointingAngle) {
430         fHistV0toXiCosineOfPointingAngle = new TH1F("fHistV0toXiCosineOfPointingAngle", "Cos. of V0 Ptng Angl / Xi vtx ;Cos(V0 Point. Angl / Xi vtx); Counts", 100, 0.99, 1.0);
431         fListHistCascade->Add(fHistV0toXiCosineOfPointingAngle);
432 }
433
434
435 if(! fHistRapXi ){
436         fHistRapXi  = new TH1F( "fHistRapXi" , "Rapidity of Xi candidates ; y ; Counts", 200, -5.0, 5.0);
437         fListHistCascade->Add(fHistRapXi);
438 }
439
440 if(! fHistRapOmega ){
441         fHistRapOmega  = new TH1F( "fHistRapOmega" , "Rapidity of Omega candidates ; y ; Counts", 200, -5.0, 5.0);
442         fListHistCascade->Add(fHistRapOmega);
443 }
444
445 if(! fHistEta ){
446         fHistEta  = new TH1F( "fHistEta" , "Pseudo-rap. of casc. candidates ; #eta ; Counts", 120, -3.0, 3.0);
447         fListHistCascade->Add(fHistEta);
448 }
449
450 if(! fHistTheta ){
451         fHistTheta  = new TH1F( "fHistTheta" , "#theta of casc. candidates ; #theta (deg) ; Counts", 180, 0., 180.0);
452         fListHistCascade->Add(fHistTheta);
453 }
454
455 if(! fHistPhi ){
456         fHistPhi  = new TH1F( "fHistPhi" , "#phi of casc. candidates ; #phi (deg) ; Counts", 360, 0., 360.);
457         fListHistCascade->Add(fHistPhi);
458 }
459
460
461 if(! f2dHistArmenteros) {
462         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);
463         fListHistCascade->Add(f2dHistArmenteros);
464 }
465
466 if(! f2dHistEffMassXiVsEffMassLambda) {
467         f2dHistEffMassXiVsEffMassLambda = new TH2F( "f2dHistEffMassXiVsEffMassLambda", "M_{#Xi^{-} candidates} Vs M_{#Lambda} ; M( #Lambda , #pi^{-} ) (GeV/c^{2}) ; Inv. M_{#Lambda^{0}} (GeV/c^{2})", 200, 1.2, 2.0, 300, 1.1,1.13);
468         fListHistCascade->Add(f2dHistEffMassXiVsEffMassLambda);
469 }
470
471 if(! f2dHistXiRadiusVsEffMassXi) {
472         f2dHistXiRadiusVsEffMassXi = new TH2F( "f2dHistXiRadiusVsEffMassXi", "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(f2dHistXiRadiusVsEffMassXi);
474 }
475
476
477
478
479
480
481
482 }// end CreateOutputObjects
483
484
485
486
487
488
489 //________________________________________________________________________
490 void AliAnalysisTaskESDCheckCascade::Exec(Option_t *) 
491 {
492   // Main loop
493   // Called for each event
494
495   if (!fESD) {
496     Printf("ERROR: fESD not available \n");
497     return;
498   }
499  
500   //  Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
501
502
503
504 //---------------------------------------------------
505 //   fESD->ResetCascades();
506 //      AliCascadeVertexer CascVtxer;
507 //      CascVtxer.V0sTracks2CascadeVertices(fESD);
508 //---------------------------------------------------
509
510
511
512
513
514
515
516   // ---------------------------------------------------------------
517   // I - Initialisation of the part dedicated to cascade vertices
518   
519   Int_t ncascades = -1;
520   ncascades = fESD->GetNumberOfCascades();
521
522         // - General histos (filled for any event)
523
524         fHistTrackMultiplicity->Fill(fESD->GetNumberOfTracks());
525         fHistCascadeMultiplicity->Fill( ncascades );
526
527
528
529
530         Short_t lStatusTrackingPrimVtx = -2;
531   
532         // - 1st part of initialisation : variables needed to store AliESDCascade data members
533         Double_t lEffMassXi  = 0. ;
534         Double_t lChi2Xi = 0. ;
535         Double_t lDcaXiDaughters = 0. ;
536         Double_t lXiCosineOfPointingAngle = 0. ;
537         Double_t lPosXi[3] = { -1000.0, -1000.0, -1000.0 };
538         Double_t lXiRadius= 0. ;
539         
540                 
541         // - 2nd part of initialisation : about V0 part in cascades
542
543         Double_t lInvMassLambdaAsCascDghter = 0.;
544         Double_t lV0Chi2Xi = 0. ;
545         Double_t lDcaV0DaughtersXi = 0.;
546         
547         Double_t lDcaBachToPrimVertexXi = 0. , lDcaV0ToPrimVertexXi = 0.;
548         Double_t lDcaPosToPrimVertexXi = 0.;
549         Double_t lDcaNegToPrimVertexXi = 0.;
550         Double_t lV0CosineOfPointingAngleXi = 0. ;
551         Double_t lPosV0Xi[3] = { -1000. , -1000., -1000. }; // Position of VO coming from cascade
552         Double_t lV0RadiusXi = -1000.0;
553
554         
555         // - 3rd part of initialisation : Effective masses
556         Double_t lV0quality = 0.;
557                 Double_t lInvMassXiMinus = 0.;
558                 Double_t lInvMassXiPlus = 0.;
559                 Double_t lInvMassOmegaMinus = 0.;
560                 Double_t lInvMassOmegaPlus = 0.;
561
562
563         // - 4th part of initialisation : extra info for QA
564         Double_t lXiMomX  = 0.;         Double_t lXiMomY  = 0.;         Double_t lXiMomZ   = 0.;
565         Double_t lXiTransvMom  = 0. ;
566         Double_t lXiTotMom  = 0. ;
567         
568         Double_t lBachMomX  = 0.;       Double_t lBachMomY  = 0.;       Double_t 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
576
577                 
578         
579         
580   
581   
582   // -------------------------------------
583   // II - Calcultaion Part dedicated to Xi vertices
584     
585    
586   for (Int_t iXi = 0; iXi < ncascades; iXi++) 
587     {// This is the begining of the Cascade loop
588       AliESDcascade *xi = fESD->GetCascade(iXi);
589           
590                 if (!xi) continue;
591                 
592       // Just to know which file is currently open : locate the file containing Xi
593       // cout << "Name of the file containing Xi candidate(s) :" <<  fesdH->GetTree()->GetCurrentFile()->GetName() << endl;
594         
595
596                 // - II.Step 1 : Characteristics of the event : prim. Vtx + magnetic field
597                 //-------------
598         const AliESDVertex *lPrimaryTrackingVtx = fESD->GetPrimaryVertexTracks();  
599                 // get the vtx stored in ESD found with tracks
600         Double_t lTrkgPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
601                 lPrimaryTrackingVtx->GetXYZ( lTrkgPrimaryVtxPos );
602                 Double_t lTrkgPrimaryVtxRadius3D = TMath::Sqrt( lTrkgPrimaryVtxPos[0]*lTrkgPrimaryVtxPos[0] +
603                                                 lTrkgPrimaryVtxPos[1] * lTrkgPrimaryVtxPos[1] +
604                                                 lTrkgPrimaryVtxPos[2] * lTrkgPrimaryVtxPos[2] );
605
606                 lStatusTrackingPrimVtx = lPrimaryTrackingVtx->GetStatus();
607
608
609         const AliESDVertex *lPrimaryBestVtx = fESD->GetPrimaryVertex(); 
610                 // get the best primary vertex available for the event
611                 // As done in AliCascadeVertexer, we keep the one which is the best one available.
612                 // between SPD vertex and Tracking vertex.
613                 // This one will be used for next calculations (DCA essentially)
614         Double_t lBestPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
615                 lPrimaryBestVtx->GetXYZ( lBestPrimaryVtxPos );
616                 Double_t lBestPrimaryVtxRadius3D = TMath::Sqrt( lBestPrimaryVtxPos[0]*lBestPrimaryVtxPos[0] +
617                                                 lBestPrimaryVtxPos[1] * lBestPrimaryVtxPos[1] +
618                                                 lBestPrimaryVtxPos[2] * lBestPrimaryVtxPos[2] );
619         
620                         
621         Double_t lMagneticField = fESD->GetMagneticField( );
622         
623         
624         
625                 // - II.Step 2 : Assigning the necessary variables for specific AliESDcascade data members
626                 //-------------
627         lV0quality = 0.;
628         xi->ChangeMassHypothesis(lV0quality , 3312); // default working hypothesis : cascade = Xi- decay
629
630         lEffMassXi                      = xi->GetEffMassXi();
631         lChi2Xi                         = xi->GetChi2Xi();
632         lDcaXiDaughters                 = xi->GetDcaXiDaughters();
633         lXiCosineOfPointingAngle        = xi->GetCascadeCosineOfPointingAngle( lBestPrimaryVtxPos[0], lBestPrimaryVtxPos[1], lBestPrimaryVtxPos[2] );
634                 // Take care : the best available vertex should be used (like in AliCascadeVertexer)
635         
636         xi->GetXYZcascade( lPosXi[0],  lPosXi[1], lPosXi[2] ); 
637                 lXiRadius               = TMath::Sqrt( lPosXi[0]*lPosXi[0]  +  lPosXi[1]*lPosXi[1] );
638                 
639                 
640
641                 // - II.Step 3 : around the tracks : Bach + V0 
642                 // ~ Necessary variables for ESDcascade data members coming from the ESDv0 part (inheritance)
643                 //-------------
644                 
645                 UInt_t lIdxPosXi        = (UInt_t) TMath::Abs( xi->GetPindex() );
646                 UInt_t lIdxNegXi        = (UInt_t) TMath::Abs( xi->GetNindex() );
647                 UInt_t lBachIdx         = (UInt_t) TMath::Abs( xi->GetBindex() );
648                         // Care track label can be negative in MC production (linked with the track quality)
649
650         AliESDtrack *pTrackXi           = fESD->GetTrack( lIdxPosXi );
651         AliESDtrack *nTrackXi           = fESD->GetTrack( lIdxNegXi );
652         AliESDtrack *bachTrackXi        = fESD->GetTrack( lBachIdx );
653         if (!pTrackXi || !nTrackXi || !bachTrackXi ) {
654                 Printf("ERROR: Could not retrieve one of the 3 daughter tracks of the cascade ...");
655                 continue;
656         }
657         
658         lInvMassLambdaAsCascDghter      = xi->GetEffMass();
659                 // This value shouldn't change, whatever the working hyp. is : Xi-, Xi+, Omega-, Omega+
660         lDcaV0DaughtersXi               = xi->GetDcaV0Daughters(); 
661         lV0Chi2Xi                       = xi->GetChi2V0();
662         
663         lV0CosineOfPointingAngleXi      = xi->GetV0CosineOfPointingAngle( lBestPrimaryVtxPos[0], lBestPrimaryVtxPos[1], lBestPrimaryVtxPos[2] );
664
665         lDcaV0ToPrimVertexXi            = xi->GetD( lBestPrimaryVtxPos[0], lBestPrimaryVtxPos[1], lBestPrimaryVtxPos[2] );
666         
667         
668         
669         lDcaBachToPrimVertexXi  = bachTrackXi->GetD( lBestPrimaryVtxPos[0], lBestPrimaryVtxPos[1], lMagneticField  ); // return an algebraic value ...
670         lDcaBachToPrimVertexXi  = TMath::Abs( lDcaBachToPrimVertexXi );
671         
672         
673         
674         
675                 xi->GetXYZ( lPosV0Xi[0],  lPosV0Xi[1], lPosV0Xi[2] ); 
676                 lV0RadiusXi     = TMath::Sqrt( lPosV0Xi[0]*lPosV0Xi[0]  +  lPosV0Xi[1]*lPosV0Xi[1] );
677         
678         Float_t  tDcaPosToPrimVertexXi[2];
679         if(pTrackXi) pTrackXi->GetImpactParameters(tDcaPosToPrimVertexXi[0],tDcaPosToPrimVertexXi[1]);
680                 // void GetImpactParameters(Float_t &xy,Float_t &z) const {xy=fD; z=fZ;}
681         else { 
682                         tDcaPosToPrimVertexXi[0]= -999.;  
683                         tDcaPosToPrimVertexXi[1]= -999.;
684                 }
685         lDcaPosToPrimVertexXi = TMath::Sqrt(tDcaPosToPrimVertexXi[0]*tDcaPosToPrimVertexXi[0]
686                                         +   tDcaPosToPrimVertexXi[1]*tDcaPosToPrimVertexXi[1]);
687
688         Float_t  tDcaNegToPrimVertexXi[2];
689         if(nTrackXi) nTrackXi->GetImpactParameters(tDcaNegToPrimVertexXi[0],tDcaNegToPrimVertexXi[1]);
690         else { 
691                         tDcaNegToPrimVertexXi[0]= -999.;  
692                         tDcaNegToPrimVertexXi[1]= -999.;
693                 }
694         lDcaNegToPrimVertexXi = TMath::Sqrt(tDcaNegToPrimVertexXi[0]*tDcaNegToPrimVertexXi[0]
695                                         +   tDcaNegToPrimVertexXi[1]*tDcaNegToPrimVertexXi[1]);
696         
697                 
698
699
700                 // - II.Step 4 : around effective masses
701                 // ~ change mass hypotheses to cover all the possibilities :  Xi-/+, Omega -/+
702                 //-------------
703
704         lV0quality = 0.;
705         xi->ChangeMassHypothesis(lV0quality , 3312);    // Calculate the effective mass of the Xi- candidate. 
706                                                         // pdg code 3312 = Xi-
707                 if( bachTrackXi->Charge() < 0 )         
708                         lInvMassXiMinus = xi->GetEffMassXi();
709                         
710         lV0quality = 0.;
711         xi->ChangeMassHypothesis(lV0quality , -3312);   // Calculate the effective mass of the Xi+ candidate. 
712                                                         // pdg code -3312 = Xi+
713                 if( bachTrackXi->Charge() >  0 )
714                         lInvMassXiPlus = xi->GetEffMassXi();
715                         
716         lV0quality = 0.;
717         xi->ChangeMassHypothesis(lV0quality , 3334);    // Calculate the effective mass of the Xi- candidate. 
718                                                         // pdg code 3334 = Omega-
719                 if( bachTrackXi->Charge() < 0 )
720                         lInvMassOmegaMinus = xi->GetEffMassXi();
721                 
722         lV0quality = 0.;
723         xi->ChangeMassHypothesis(lV0quality , -3334);   // Calculate the effective mass of the Xi+ candidate. 
724                                                         // pdg code -3334  = Omega+
725                 if( bachTrackXi->Charge() >  0 )
726                         lInvMassOmegaPlus = xi->GetEffMassXi();
727
728         lV0quality = 0.;
729         xi->ChangeMassHypothesis(lV0quality , 3312);    // Back to default hyp.
730                                                         
731
732
733                 // - II.Step 5 : extra info for QA
734                 // miscellaneous pieces onf info that may help regarding data quality assessment.
735                 //-------------
736
737         xi->GetPxPyPz( lXiMomX, lXiMomY, lXiMomZ );
738                 lXiTransvMom    = TMath::Sqrt( lXiMomX*lXiMomX   + lXiMomY*lXiMomY );
739                 lXiTotMom       = TMath::Sqrt( lXiMomX*lXiMomX   + lXiMomY*lXiMomY   + lXiMomZ*lXiMomZ );
740                 
741         xi->GetBPxPyPz(  lBachMomX,  lBachMomY,  lBachMomZ );
742                 lBachTransvMom          = TMath::Sqrt( lBachMomX*lBachMomX   + lBachMomY*lBachMomY );
743                 lBachTotMom             = TMath::Sqrt( lBachMomX*lBachMomX   + lBachMomY*lBachMomY  +  lBachMomZ*lBachMomZ  );
744
745
746         lChargeXi = xi->Charge();
747
748         lV0toXiCosineOfPointingAngle = xi->GetV0CosineOfPointingAngle( lPosXi[0], lPosXi[1], lPosXi[2] );
749
750         // Double_t lRapXi, lRapOmega, lEta, lTheta, lPhi ;
751
752  
753
754
755
756   // -------------------------------------
757   // III - Filling the TH1Fs
758   
759
760                 // - III.Step 1 
761
762                 fHistVtxStatus->Fill( lStatusTrackingPrimVtx );  // 1 if tracking vtx = ok
763
764                 if( lStatusTrackingPrimVtx ){
765                         fHistPosTrkgPrimaryVtxX->Fill( lTrkgPrimaryVtxPos[0]  );
766                         fHistPosTrkgPrimaryVtxY->Fill( lTrkgPrimaryVtxPos[1]  );        
767                         fHistPosTrkgPrimaryVtxZ->Fill( lTrkgPrimaryVtxPos[2]  ); 
768                         fHistTrkgPrimaryVtxRadius->Fill( lTrkgPrimaryVtxRadius3D );
769                 }
770
771                 fHistPosBestPrimaryVtxX->Fill( lBestPrimaryVtxPos[0] );
772                 fHistPosBestPrimaryVtxY->Fill( lBestPrimaryVtxPos[1] ); 
773                 fHistPosBestPrimaryVtxZ->Fill( lBestPrimaryVtxPos[2] ); 
774                 fHistBestPrimaryVtxRadius->Fill( lBestPrimaryVtxRadius3D  );
775                 
776                 f2dHistTrkgPrimVtxVsBestPrimVtx->Fill( lTrkgPrimaryVtxRadius3D, lBestPrimaryVtxRadius3D );
777
778                         
779                 // - III.Step 2
780                 fHistEffMassXi->Fill( lEffMassXi   );
781                 fHistChi2Xi->Fill( lChi2Xi  );                  // Flag CascadeVtxer: Cut Variable a
782                 fHistDcaXiDaughters->Fill( lDcaXiDaughters );   // Flag CascadeVtxer: Cut Variable e 
783                 fHistDcaBachToPrimVertex->Fill( lDcaBachToPrimVertexXi  );      // Flag CascadeVtxer: Cut Variable d
784                 fHistXiCosineOfPointingAngle->Fill( lXiCosineOfPointingAngle ); // Flag CascadeVtxer: Cut Variable f
785                 fHistXiRadius->Fill( lXiRadius );               // Flag CascadeVtxer: Cut Variable g+h
786                 
787                 
788                 // - III.Step 3
789                 fHistMassLambdaAsCascDghter->Fill(  lInvMassLambdaAsCascDghter  ); // Flag CascadeVtxer: Cut Variable c
790                 fHistV0Chi2Xi->Fill( lV0Chi2Xi  );      
791                 fHistDcaV0DaughtersXi->Fill( lDcaV0DaughtersXi );
792                 fHistV0CosineOfPointingAngleXi->Fill( lV0CosineOfPointingAngleXi ); 
793                 fHistV0RadiusXi->Fill( lV0RadiusXi );
794                 
795                 fHistDcaV0ToPrimVertexXi->Fill( lDcaV0ToPrimVertexXi ); // Flag CascadeVtxer: Cut Variable b
796                 fHistDcaPosToPrimVertexXi->Fill( lDcaPosToPrimVertexXi );
797                 fHistDcaNegToPrimVertexXi->Fill( lDcaNegToPrimVertexXi );
798                         
799         
800                 // - III.Step 4
801                 if( lChargeXi < 0 )     fHistMassXiMinus->Fill( lInvMassXiMinus );
802                 if( lChargeXi > 0 )     fHistMassXiPlus->Fill( lInvMassXiPlus );
803                 if( lChargeXi < 0 )     fHistMassOmegaMinus->Fill( lInvMassOmegaMinus );
804                 if( lChargeXi > 0 )     fHistMassOmegaPlus->Fill( lInvMassOmegaPlus );  
805
806
807                 // - III.Step 5
808                 fHistXiTransvMom->Fill( lXiTransvMom );
809                 fHistXiTotMom->Fill( lXiTotMom );
810                 
811                 fHistBachTransvMom->Fill( lBachTransvMom );
812                 fHistBachTotMom->Fill( lBachTotMom );
813
814                 fHistChargeXi->Fill( lChargeXi );
815                 fHistV0toXiCosineOfPointingAngle->Fill( lV0toXiCosineOfPointingAngle );
816         
817                 fHistRapXi->Fill( xi->RapXi() );
818                 fHistRapOmega->Fill( xi->RapOmega()  );
819                 fHistEta->Fill( xi->Eta() );
820                 fHistTheta->Fill( xi->Theta() *180.0/TMath::Pi() );
821                 fHistPhi->Fill( xi->Phi() *180.0/TMath::Pi() );
822
823                 f2dHistArmenteros->Fill( xi->AlphaXi(), xi->PtArmXi()  );
824                 if( lChargeXi < 0 ) f2dHistEffMassXiVsEffMassLambda->Fill( lInvMassXiMinus, lInvMassLambdaAsCascDghter ); 
825                 f2dHistXiRadiusVsEffMassXi->Fill( lXiRadius, lInvMassXiMinus );
826
827       
828
829  
830     }// This is the end of the Cascade loop
831
832  
833   // Post output data.
834  PostData(0, fListHistCascade);
835 }
836
837
838
839
840
841
842
843
844
845
846 //________________________________________________________________________
847 void AliAnalysisTaskESDCheckCascade::Terminate(Option_t *) 
848 {
849   // Draw result to the screen
850   // Called once at the end of the query
851
852   fHistTrackMultiplicity = dynamic_cast<TH1F*> (  ((TList*)GetOutputData(0))->FindObject("fHistTrackMultiplicity")  );
853   if (!fHistTrackMultiplicity) {
854     Printf("ERROR: fHistTrackMultiplicity not available");
855     return;
856   }
857  
858   fHistCascadeMultiplicity = dynamic_cast<TH1F*> (  ((TList*)GetOutputData(0))->FindObject("fHistCascadeMultiplicity")  );
859   if (!fHistCascadeMultiplicity) {
860     Printf("ERROR: fHistCascadeMultiplicity not available");
861     return;
862   }
863    
864   TCanvas *c2 = new TCanvas("AliAnalysisTaskESDCheckCascade","Multiplicity",10,10,510,510);
865   c2->cd(1)->SetLogy();
866
867   fHistTrackMultiplicity->SetMarkerStyle(22);
868   fHistTrackMultiplicity->DrawCopy("E");
869   fHistCascadeMultiplicity->SetMarkerStyle(26);
870   fHistCascadeMultiplicity->DrawCopy("ESAME");
871   
872 }