2 /**************************************************************************
3 * Authors : Antonin Maire, Boris Hippolyte *
4 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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 //-----------------------------------------------------------------
29 class AliMCEventHandler;
47 #include "AliAnalysisManager.h"
48 #include "AliESDInputHandler.h"
50 #include "AliESDEvent.h"
51 //#include "AliCascadeVertexer.h"
53 #include "AliESDcascade.h"
55 #include "AliAnalysisTaskESDCheckCascade.h"
57 ClassImp(AliAnalysisTaskESDCheckCascade)
61 //________________________________________________________________________
62 AliAnalysisTaskESDCheckCascade::AliAnalysisTaskESDCheckCascade()
65 fesdH(0), fMCevent(0), fESD(0),
67 // - Cascade part initialisation
69 fHistTrackMultiplicity(0), fHistCascadeMultiplicity(0),
72 fHistPosTrkgPrimaryVtxX(0), fHistPosTrkgPrimaryVtxY(0), fHistPosTrkgPrimaryVtxZ(0), fHistTrkgPrimaryVtxRadius(0),
73 fHistPosBestPrimaryVtxX(0), fHistPosBestPrimaryVtxY(0), fHistPosBestPrimaryVtxZ(0), fHistBestPrimaryVtxRadius(0),
74 f2dHistTrkgPrimVtxVsBestPrimVtx(0),
76 fHistEffMassXi(0), fHistChi2Xi(0),
77 fHistDcaXiDaughters(0), fHistDcaBachToPrimVertex(0), fHistXiCosineOfPointingAngle(0), fHistXiRadius(0),
79 fHistMassLambdaAsCascDghter(0),
81 fHistDcaV0DaughtersXi(0),
82 fHistDcaV0ToPrimVertexXi(0),
83 fHistV0CosineOfPointingAngleXi(0),
85 fHistDcaPosToPrimVertexXi(0), fHistDcaNegToPrimVertexXi(0),
87 fHistMassXiMinus(0), fHistMassXiPlus(0),
88 fHistMassOmegaMinus(0), fHistMassOmegaPlus(0),
90 fHistXiTransvMom(0), fHistXiTotMom(0),
91 fHistBachTransvMom(0), fHistBachTotMom(0),
94 fHistV0toXiCosineOfPointingAngle(0),
96 fHistRapXi(0), fHistRapOmega(0), fHistEta(0),
97 fHistTheta(0), fHistPhi(0),
100 f2dHistEffMassXiVsEffMassLambda(0),
101 f2dHistXiRadiusVsEffMassXi(0)
114 //________________________________________________________________________
115 AliAnalysisTaskESDCheckCascade::AliAnalysisTaskESDCheckCascade(const char *name)
116 : AliAnalysisTask(name, ""),
118 fesdH(0), fMCevent(0), fESD(0),
120 // - Cascade part initialisation
122 fHistTrackMultiplicity(0), fHistCascadeMultiplicity(0),
125 fHistPosTrkgPrimaryVtxX(0), fHistPosTrkgPrimaryVtxY(0), fHistPosTrkgPrimaryVtxZ(0), fHistTrkgPrimaryVtxRadius(0),
126 fHistPosBestPrimaryVtxX(0), fHistPosBestPrimaryVtxY(0), fHistPosBestPrimaryVtxZ(0), fHistBestPrimaryVtxRadius(0),
127 f2dHistTrkgPrimVtxVsBestPrimVtx(0),
129 fHistEffMassXi(0), fHistChi2Xi(0),
130 fHistDcaXiDaughters(0), fHistDcaBachToPrimVertex(0), fHistXiCosineOfPointingAngle(0), fHistXiRadius(0),
132 fHistMassLambdaAsCascDghter(0),
134 fHistDcaV0DaughtersXi(0),
135 fHistDcaV0ToPrimVertexXi(0),
136 fHistV0CosineOfPointingAngleXi(0),
138 fHistDcaPosToPrimVertexXi(0), fHistDcaNegToPrimVertexXi(0),
140 fHistMassXiMinus(0), fHistMassXiPlus(0),
141 fHistMassOmegaMinus(0), fHistMassOmegaPlus(0),
143 fHistXiTransvMom(0), fHistXiTotMom(0),
144 fHistBachTransvMom(0), fHistBachTotMom(0),
147 fHistV0toXiCosineOfPointingAngle(0),
149 fHistRapXi(0), fHistRapOmega(0), fHistEta(0),
150 fHistTheta(0), fHistPhi(0),
152 f2dHistArmenteros(0),
153 f2dHistEffMassXiVsEffMassLambda(0),
154 f2dHistXiRadiusVsEffMassXi(0)
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());
174 //________________________________________________________________________
175 void AliAnalysisTaskESDCheckCascade::ConnectInputData(Option_t *)
177 // Connect ESD or AOD here
180 AliLog::SetGlobalLogLevel(AliLog::kError); // to suppress the extensive info prompted by a run with MC
182 TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
184 Printf("ERROR: Could not read chain from input slot 0");
187 fesdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
188 //AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
191 Printf("ERROR: Could not get ESDInputHandler");
193 fESD = fesdH->GetEvent();
203 //________________________________________________________________________
204 void AliAnalysisTaskESDCheckCascade::CreateOutputObjects()
211 fListHistCascade = new TList();
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);
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);
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);
240 // - Vertex Positions
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);
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);
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);
258 if(! fHistTrkgPrimaryVtxRadius ){
259 fHistTrkgPrimaryVtxRadius = new TH1F( "fHistTrkgPrimaryVtxRadius", "Trkg Prim. Vertex radius; r (cm); Events" , 150, 0., 15.0 );
260 fListHistCascade->Add(fHistTrkgPrimaryVtxRadius);
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);
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);
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);
281 if(! fHistBestPrimaryVtxRadius ){
282 fHistBestPrimaryVtxRadius = new TH1F( "fHistBestPrimaryVtxRadius", "Best Prim. vertex radius; r (cm); Events" , 150, 0., 15.0 );
283 fListHistCascade->Add(fHistBestPrimaryVtxRadius);
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);
294 // - Typical histos for cascades
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);
303 fHistChi2Xi = new TH1F("fHistChi2Xi", "Cascade #chi^{2}; #chi^{2}; Number of Cascades", 160, 0, 40);
304 fListHistCascade->Add(fHistChi2Xi);
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);
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);
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);
322 if(! fHistXiRadius ){
323 fHistXiRadius = new TH1F( "fHistXiRadius", "Casc. decay transv. radius; r (cm); Counts" , 200, 0., 20.0 );
324 fListHistCascade->Add(fHistXiRadius);
328 // - Histos about ~ the "V0 part" of the cascade, coming by inheritance from AliESDv0
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);
337 if (! fHistV0Chi2Xi) {
338 fHistV0Chi2Xi = new TH1F("fHistV0Chi2Xi", "V0 #chi^{2}, in cascade; #chi^{2};Counts", 160, 0, 40);
339 fListHistCascade->Add(fHistV0Chi2Xi);
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);
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);
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);
357 if (! fHistV0RadiusXi) {
358 fHistV0RadiusXi = new TH1F("fHistV0RadiusXi", "V0 decay radius, in cascade; radius (cm); Counts", 200, 0, 20);
359 fListHistCascade->Add(fHistV0RadiusXi);
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);
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);
375 // - Effective mass histos for cascades.
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);
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);
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);
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);
399 // - Complements for QA
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);
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);
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);
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);
423 if(! fHistChargeXi ){
424 fHistChargeXi = new TH1F( "fHistChargeXi" , "Charge of casc. candidates ; Sign ; Counts", 5, -2.0, 3.0);
425 fListHistCascade->Add(fHistChargeXi);
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);
436 fHistRapXi = new TH1F( "fHistRapXi" , "Rapidity of Xi candidates ; y ; Counts", 200, -5.0, 5.0);
437 fListHistCascade->Add(fHistRapXi);
440 if(! fHistRapOmega ){
441 fHistRapOmega = new TH1F( "fHistRapOmega" , "Rapidity of Omega candidates ; y ; Counts", 200, -5.0, 5.0);
442 fListHistCascade->Add(fHistRapOmega);
446 fHistEta = new TH1F( "fHistEta" , "Pseudo-rap. of casc. candidates ; #eta ; Counts", 120, -3.0, 3.0);
447 fListHistCascade->Add(fHistEta);
451 fHistTheta = new TH1F( "fHistTheta" , "#theta of casc. candidates ; #theta (deg) ; Counts", 180, 0., 180.0);
452 fListHistCascade->Add(fHistTheta);
456 fHistPhi = new TH1F( "fHistPhi" , "#phi of casc. candidates ; #phi (deg) ; Counts", 360, 0., 360.);
457 fListHistCascade->Add(fHistPhi);
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);
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);
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);
482 }// end CreateOutputObjects
489 //________________________________________________________________________
490 void AliAnalysisTaskESDCheckCascade::Exec(Option_t *)
493 // Called for each event
496 Printf("ERROR: fESD not available \n");
500 // Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
504 //---------------------------------------------------
505 // fESD->ResetCascades();
506 // AliCascadeVertexer CascVtxer;
507 // CascVtxer.V0sTracks2CascadeVertices(fESD);
508 //---------------------------------------------------
516 // ---------------------------------------------------------------
517 // I - Initialisation of the part dedicated to cascade vertices
519 Int_t ncascades = -1;
520 ncascades = fESD->GetNumberOfCascades();
522 // - General histos (filled for any event)
524 fHistTrackMultiplicity->Fill(fESD->GetNumberOfTracks());
525 fHistCascadeMultiplicity->Fill( ncascades );
530 Short_t lStatusTrackingPrimVtx = -2;
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. ;
541 // - 2nd part of initialisation : about V0 part in cascades
543 Double_t lInvMassLambdaAsCascDghter = 0.;
544 Double_t lV0Chi2Xi = 0. ;
545 Double_t lDcaV0DaughtersXi = 0.;
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;
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.;
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. ;
568 Double_t lBachMomX = 0.; Double_t lBachMomY = 0.; Double_t lBachMomZ = 0.;
569 Double_t lBachTransvMom = 0.;
570 Double_t lBachTotMom = 0.;
572 Short_t lChargeXi = -2;
573 Double_t lV0toXiCosineOfPointingAngle = 0. ;
582 // -------------------------------------
583 // II - Calcultaion Part dedicated to Xi vertices
586 for (Int_t iXi = 0; iXi < ncascades; iXi++)
587 {// This is the begining of the Cascade loop
588 AliESDcascade *xi = fESD->GetCascade(iXi);
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;
596 // - II.Step 1 : Characteristics of the event : prim. Vtx + magnetic field
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] );
606 lStatusTrackingPrimVtx = lPrimaryTrackingVtx->GetStatus();
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] );
621 Double_t lMagneticField = fESD->GetMagneticField( );
625 // - II.Step 2 : Assigning the necessary variables for specific AliESDcascade data members
628 xi->ChangeMassHypothesis(lV0quality , 3312); // default working hypothesis : cascade = Xi- decay
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)
636 xi->GetXYZcascade( lPosXi[0], lPosXi[1], lPosXi[2] );
637 lXiRadius = TMath::Sqrt( lPosXi[0]*lPosXi[0] + lPosXi[1]*lPosXi[1] );
641 // - II.Step 3 : around the tracks : Bach + V0
642 // ~ Necessary variables for ESDcascade data members coming from the ESDv0 part (inheritance)
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)
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 ...");
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();
663 lV0CosineOfPointingAngleXi = xi->GetV0CosineOfPointingAngle( lBestPrimaryVtxPos[0], lBestPrimaryVtxPos[1], lBestPrimaryVtxPos[2] );
665 lDcaV0ToPrimVertexXi = xi->GetD( lBestPrimaryVtxPos[0], lBestPrimaryVtxPos[1], lBestPrimaryVtxPos[2] );
669 lDcaBachToPrimVertexXi = bachTrackXi->GetD( lBestPrimaryVtxPos[0], lBestPrimaryVtxPos[1], lMagneticField ); // return an algebraic value ...
670 lDcaBachToPrimVertexXi = TMath::Abs( lDcaBachToPrimVertexXi );
675 xi->GetXYZ( lPosV0Xi[0], lPosV0Xi[1], lPosV0Xi[2] );
676 lV0RadiusXi = TMath::Sqrt( lPosV0Xi[0]*lPosV0Xi[0] + lPosV0Xi[1]*lPosV0Xi[1] );
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;}
682 tDcaPosToPrimVertexXi[0]= -999.;
683 tDcaPosToPrimVertexXi[1]= -999.;
685 lDcaPosToPrimVertexXi = TMath::Sqrt(tDcaPosToPrimVertexXi[0]*tDcaPosToPrimVertexXi[0]
686 + tDcaPosToPrimVertexXi[1]*tDcaPosToPrimVertexXi[1]);
688 Float_t tDcaNegToPrimVertexXi[2];
689 if(nTrackXi) nTrackXi->GetImpactParameters(tDcaNegToPrimVertexXi[0],tDcaNegToPrimVertexXi[1]);
691 tDcaNegToPrimVertexXi[0]= -999.;
692 tDcaNegToPrimVertexXi[1]= -999.;
694 lDcaNegToPrimVertexXi = TMath::Sqrt(tDcaNegToPrimVertexXi[0]*tDcaNegToPrimVertexXi[0]
695 + tDcaNegToPrimVertexXi[1]*tDcaNegToPrimVertexXi[1]);
700 // - II.Step 4 : around effective masses
701 // ~ change mass hypotheses to cover all the possibilities : Xi-/+, Omega -/+
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();
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();
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();
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();
729 xi->ChangeMassHypothesis(lV0quality , 3312); // Back to default hyp.
733 // - II.Step 5 : extra info for QA
734 // miscellaneous pieces onf info that may help regarding data quality assessment.
737 xi->GetPxPyPz( lXiMomX, lXiMomY, lXiMomZ );
738 lXiTransvMom = TMath::Sqrt( lXiMomX*lXiMomX + lXiMomY*lXiMomY );
739 lXiTotMom = TMath::Sqrt( lXiMomX*lXiMomX + lXiMomY*lXiMomY + lXiMomZ*lXiMomZ );
741 xi->GetBPxPyPz( lBachMomX, lBachMomY, lBachMomZ );
742 lBachTransvMom = TMath::Sqrt( lBachMomX*lBachMomX + lBachMomY*lBachMomY );
743 lBachTotMom = TMath::Sqrt( lBachMomX*lBachMomX + lBachMomY*lBachMomY + lBachMomZ*lBachMomZ );
746 lChargeXi = xi->Charge();
748 lV0toXiCosineOfPointingAngle = xi->GetV0CosineOfPointingAngle( lPosXi[0], lPosXi[1], lPosXi[2] );
750 // Double_t lRapXi, lRapOmega, lEta, lTheta, lPhi ;
756 // -------------------------------------
757 // III - Filling the TH1Fs
762 fHistVtxStatus->Fill( lStatusTrackingPrimVtx ); // 1 if tracking vtx = ok
764 if( lStatusTrackingPrimVtx ){
765 fHistPosTrkgPrimaryVtxX->Fill( lTrkgPrimaryVtxPos[0] );
766 fHistPosTrkgPrimaryVtxY->Fill( lTrkgPrimaryVtxPos[1] );
767 fHistPosTrkgPrimaryVtxZ->Fill( lTrkgPrimaryVtxPos[2] );
768 fHistTrkgPrimaryVtxRadius->Fill( lTrkgPrimaryVtxRadius3D );
771 fHistPosBestPrimaryVtxX->Fill( lBestPrimaryVtxPos[0] );
772 fHistPosBestPrimaryVtxY->Fill( lBestPrimaryVtxPos[1] );
773 fHistPosBestPrimaryVtxZ->Fill( lBestPrimaryVtxPos[2] );
774 fHistBestPrimaryVtxRadius->Fill( lBestPrimaryVtxRadius3D );
776 f2dHistTrkgPrimVtxVsBestPrimVtx->Fill( lTrkgPrimaryVtxRadius3D, lBestPrimaryVtxRadius3D );
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
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 );
795 fHistDcaV0ToPrimVertexXi->Fill( lDcaV0ToPrimVertexXi ); // Flag CascadeVtxer: Cut Variable b
796 fHistDcaPosToPrimVertexXi->Fill( lDcaPosToPrimVertexXi );
797 fHistDcaNegToPrimVertexXi->Fill( lDcaNegToPrimVertexXi );
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 );
808 fHistXiTransvMom->Fill( lXiTransvMom );
809 fHistXiTotMom->Fill( lXiTotMom );
811 fHistBachTransvMom->Fill( lBachTransvMom );
812 fHistBachTotMom->Fill( lBachTotMom );
814 fHistChargeXi->Fill( lChargeXi );
815 fHistV0toXiCosineOfPointingAngle->Fill( lV0toXiCosineOfPointingAngle );
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() );
823 f2dHistArmenteros->Fill( xi->AlphaXi(), xi->PtArmXi() );
824 if( lChargeXi < 0 ) f2dHistEffMassXiVsEffMassLambda->Fill( lInvMassXiMinus, lInvMassLambdaAsCascDghter );
825 f2dHistXiRadiusVsEffMassXi->Fill( lXiRadius, lInvMassXiMinus );
830 }// This is the end of the Cascade loop
834 PostData(0, fListHistCascade);
846 //________________________________________________________________________
847 void AliAnalysisTaskESDCheckCascade::Terminate(Option_t *)
849 // Draw result to the screen
850 // Called once at the end of the query
852 fHistTrackMultiplicity = dynamic_cast<TH1F*> ( ((TList*)GetOutputData(0))->FindObject("fHistTrackMultiplicity") );
853 if (!fHistTrackMultiplicity) {
854 Printf("ERROR: fHistTrackMultiplicity not available");
858 fHistCascadeMultiplicity = dynamic_cast<TH1F*> ( ((TList*)GetOutputData(0))->FindObject("fHistCascadeMultiplicity") );
859 if (!fHistCascadeMultiplicity) {
860 Printf("ERROR: fHistCascadeMultiplicity not available");
864 TCanvas *c2 = new TCanvas("AliAnalysisTaskESDCheckCascade","Multiplicity",10,10,510,510);
865 c2->cd(1)->SetLogy();
867 fHistTrackMultiplicity->SetMarkerStyle(22);
868 fHistTrackMultiplicity->DrawCopy("E");
869 fHistCascadeMultiplicity->SetMarkerStyle(26);
870 fHistCascadeMultiplicity->DrawCopy("ESAME");