]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/EVE/AliHLTEveHLT.cxx
coverity 15108, the code needs further revision since the loop does not make sense...
[u/mrichter/AliRoot.git] / HLT / EVE / AliHLTEveHLT.cxx
1 // $Id$
2 /**************************************************************************
3  * This file is property of and copyright by the ALICE HLT Project        *
4  * ALICE Experiment at CERN, All rights reserved.                         *
5  *                                                                        *
6  * Primary Authors: Svein Lindal <slindal@fys.uio.no   >                  *
7  *                  for The ALICE HLT Project.                            *
8  *                                                                        *
9  * Permission to use, copy, modify and distribute this software and its   *
10  * documentation strictly for non-commercial purposes is hereby granted   *
11  * without fee, provided that the above copyright notice appears in all   *
12  * copies and that both the copyright notice and this permission notice   *
13  * appear in the supporting documentation. The authors make no claims     *
14  * about the suitability of this software for any purpose. It is          *
15  * provided "as is" without express or implied warranty.                  *
16  **************************************************************************/
17
18 /// @file   AliHLTEvePhos.cxx
19 /// @author Svein Lindal <slindal@fys.uio.no>
20 /// @brief  HLT class for the HLT EVE display
21
22 #include "AliHLTEveHLT.h"
23 #include "AliHLTHOMERBlockDesc.h"
24 #include "AliHLTEveBase.h"
25 #include "AliEveHLTEventManager.h"
26 #include "AliHLTGlobalTriggerDecision.h"
27 #include "TEveManager.h"
28 #include "TEvePointSet.h"
29 #include "TEveTrack.h"
30 #include "TEveElement.h"
31 #include "TCanvas.h"
32 #include "AliESDEvent.h"
33 #include "AliESDtrack.h"
34 #include "TEveTrackPropagator.h"
35 #include "AliEveTrack.h"
36 #include "TEveVSDStructs.h"
37 #include "TString.h"
38 #include "TPCLib/tracking-ca/AliHLTTPCCATrackParam.h"
39 #include "TPCLib/tracking-ca/AliHLTTPCCATrackConvertor.h"
40 #include "AliEveMagField.h"
41 #include "TH1.h"
42 #include "TH1F.h"
43 #include "TH2F.h"
44 #include "TThread.h"
45
46 ClassImp(AliHLTEveHLT)
47
48 AliHLTEveHLT::AliHLTEveHLT() : 
49   AliHLTEveBase("TPC tracks"), 
50   fTrueField(kFALSE),
51   fUseIpOnFailedITS(kFALSE),
52   fUseRkStepper(kFALSE),
53   fTrackList(NULL),
54   fTrackLists(NULL),
55   fOldTrackList(NULL),
56   fPointSetVertex(NULL),
57   fTrCanvas(NULL),
58   fVertexCanvas(NULL),
59   fHistEta(NULL),
60   fHistPhi(NULL),
61   fHistnClusters(NULL),
62   fHistMult(NULL),
63   fHistDCAr(NULL),
64   fTrCount(0), 
65   fVCount(0),
66   fNTrackBins(10)
67 {
68   // Constructor.
69   CreateHistograms();
70
71 }
72 ///___________________________________________________________________
73 AliHLTEveHLT::~AliHLTEveHLT()
74 {
75   //Destructor, not implemented
76   if(fTrackList)
77     delete fTrackList;
78   fTrackList = NULL;
79
80   if(fTrackLists)
81     delete fTrackLists;
82   fTrackLists = NULL;
83
84
85 }
86
87
88 ///________________________________________________________________________
89 void AliHLTEveHLT::CreateHistograms(){
90   //See header file for documentation
91
92   //fHistPt        = new TH1F("fHistPt",       "transverse momentum",    100, 0, 10); // KK   
93   //fHistP         = new TH1F("fHistP",        "signed momentum",        100,-7,  7);      
94
95   //fHistPt   ->SetXTitle("p_{t} (GeV/c)");   // KK
96   //fHistP    ->SetXTitle("P*charge (GeV/c)");
97
98   // fHistTheta     = new TH1F("fHistTheta",    "polar angle",            180, 0,180);   
99   // fHistTheta->SetXTitle("#theta (degrees)");
100   
101   fHistEta       = new TH1F("fHistEta",      "pseudorapidity",         100,-2,  2);        
102   fHistEta  ->SetXTitle("#eta");
103
104   fHistPhi       = new TH1F("fHistPhi",      "azimuthal angle",        180, 0,360);   
105   fHistPhi  ->SetXTitle("#phi (degrees)");
106
107   fHistnClusters = new TH1F("fHistnClusters","TPC clusters per track", 160, 0,160);
108
109   fHistMult      = new TH1F("fHistMult",     "event track multiplicity",150, 0, 15000);    
110   
111   fHistDCAr = new TH1F("fHistDCAr", "DCA r", 200, -100, 100);
112
113 }
114
115
116 ///_____________________________________________________________________
117 void AliHLTEveHLT::ProcessBlock(AliHLTHOMERBlockDesc * block) {
118   //See header file for documentation
119   if ( ! block->GetDataType().CompareTo("ALIESDV0") ) {
120     ProcessEsdBlock(block);
121   } 
122   
123   else if ( ! block->GetDataType().CompareTo("ROOTTOBJ") ) {
124     //processROOTTOBJ( block, gHLTText );
125   } 
126
127   else if ( ! block->GetDataType().CompareTo("HLTRDLST") ) {
128     cout << "ignoring hlt rdlst"<<endl;
129     //processHLTRDLST( block );
130   } 
131
132   else if ( ! block->GetDataType().CompareTo("GLOBTRIG") ) {
133     ProcessGlobalTrigger( block );
134   } 
135
136   else if ( !block->GetDataType().CompareTo("ROOTHIST") ) {      
137     if( !fCanvas ) { 
138       fCanvas = CreateCanvas("PVtx/Tr QA", "PrimV");
139       fCanvas->Divide(3, 3);
140     }
141     ProcessHistograms( block , fCanvas);
142   }  
143 }
144
145 ///____________________________________________________________________________
146 void AliHLTEveHLT::UpdateElements() {
147   //See header file for documentation
148
149   DrawHistograms();
150
151   if(fTrackLists) {
152     for(Int_t il = 0; il < fNTrackBins; il++) {
153       TEveTrackList * trackList = dynamic_cast<TEveTrackList*>(fTrackLists->FindChild(Form("Tracks_%d", il)));
154       if(trackList) trackList->ElementChanged();
155     } 
156   }
157
158   if(fPointSetVertex) fPointSetVertex->ResetBBox();
159   if(fTrCanvas) fTrCanvas->Update();
160   if(fVertexCanvas) fVertexCanvas->Update();
161   if(fCanvas) fCanvas->Update();
162   
163 }
164
165 ///_________________________________________________________________________________
166 void AliHLTEveHLT::ResetElements(){
167   //See header file for documentation
168   
169   cout << "destroy"<<endl;
170
171   if(fTrackLists) {
172     for(Int_t il = 0; il < fNTrackBins; il++) {
173       TEveTrackList * trackList = dynamic_cast<TEveTrackList*>(fTrackLists->FindChild(Form("Tracks_%d", il)));
174       if(trackList) trackList->DestroyElements();
175     } 
176   }
177
178   
179   fTrCount = 0;
180   fVCount = 0;
181
182   if(fPointSetVertex) fPointSetVertex->Reset();
183   cout<< "reset done"<<endl;
184   fHistoCount = 0;
185
186 }
187
188 ///_____________________________________________________________________________________
189 void * AliHLTEveHLT::DestroyGarbage(void * arg) {
190   AliHLTEveHLT * hlt = reinterpret_cast<AliHLTEveHLT*>(arg);
191   if(hlt) hlt->DestroyOldTrackList();
192   return (void*)0;
193 }
194 ///_____________________________________________________________________________________
195 void AliHLTEveHLT::DestroyOldTrackList() {
196   cout << "Destroying the old tracklist's elements"<<endl;
197   fOldTrackList->DestroyElements();
198   cout << "Destroying the old tracklist itself"<<endl;
199   fOldTrackList->Destroy();
200 }
201
202 ///_____________________________________________________________________________________
203 void AliHLTEveHLT::ProcessHistograms(AliHLTHOMERBlockDesc * block, TCanvas * canvas) {
204   //See header file for documentation
205   if (!fTrCanvas) {
206     fTrCanvas = CreateCanvas("ESD", "ESD");
207     fTrCanvas->Divide(4, 2);
208   }
209
210   if(!fVertexCanvas) {
211     fVertexCanvas = CreateCanvas("Vtx", "Vtx");
212     fVertexCanvas->Divide(4, 2);
213   }
214
215
216
217   if ( ! block->GetClassName().CompareTo("TH1F")) {
218     TH1F* histo = reinterpret_cast<TH1F*>(block->GetTObject());
219     if( histo ){
220       TString name(histo->GetName());
221       cout << "TH1F " <<name << endl;
222       if( !name.CompareTo("primVertexZ") ){
223         canvas->cd(2);
224         histo->Draw();
225       }else if( !name.CompareTo("primVertexX") ){
226         canvas->cd(3);
227         histo->Draw();
228       }else if( !name.CompareTo("primVertexY") ){
229         canvas->cd(4);
230         histo->Draw();
231       } else if ( name.Contains("Track")) {
232         AddHistogramToCanvas(histo, fTrCanvas, fTrCount);
233       } else {
234         AddHistogramToCanvas(histo, fVertexCanvas, fVCount);
235       }
236     }
237   }  else if ( ! block->GetClassName().CompareTo("TH2F")) {
238     TH2F *hista = reinterpret_cast<TH2F*>(block->GetTObject());
239     if (hista ){
240        TString name(hista->GetName());
241        cout << "TH2F " << name << endl;
242        if( !name.CompareTo("primVertexXY")) {      
243          canvas->cd(1);
244          hista->Draw();
245        } else if ( name.Contains("Track")) {
246         AddHistogramToCanvas(hista, fTrCanvas, fTrCount);
247       } else {
248         AddHistogramToCanvas(hista, fVertexCanvas, fVCount);
249       }
250     }
251   }
252   canvas->cd();
253 }
254
255 ///_____________________________________________________________________________________
256 void AliHLTEveHLT::AddHistogramToCanvas(TH1* histogram, TCanvas * canvas, Int_t &cdCount) {
257   canvas->cd(++cdCount);
258   histogram->Draw();
259 }
260
261
262 ///________________________________________________________________________________________
263 void AliHLTEveHLT::CreateTrackList() {
264   //See header file for documentation
265   fTrackList = new TEveTrackList("ESD Tracks");
266   fTrackList->SetMainColor(kOrange);
267   AddElement(fTrackList);
268 }
269
270 ///________________________________________________________________________________________
271 void AliHLTEveHLT::CreateTrackLists() {
272   //See header file for documentation
273   fTrackLists = new TEveElementList("ESD Track lists");
274   AddElement(fTrackLists);
275   for(Int_t i = 0; i < 10; i++) {
276     TEveTrackList * trackList = new TEveTrackList(Form("Tracks_%d", i));
277     trackList->SetMainColor(kOrange-i);
278     fTrackLists->AddElement(trackList);
279   }
280 }
281
282
283 ///_________________________________________________________________________________________
284 void AliHLTEveHLT::CreateVertexPointSet() {
285   //See header file for documentation
286   fPointSetVertex = new TEvePointSet("Primary Vertex");
287   fPointSetVertex->SetMainColor(6);
288   fPointSetVertex->SetMarkerStyle((Style_t)kFullStar);
289
290   AddElement(fPointSetVertex);
291 }
292
293 ///________________________________________________________________________
294 void AliHLTEveHLT::ProcessGlobalTrigger( AliHLTHOMERBlockDesc * block ) {
295   //See header file for documentation
296   AliHLTGlobalTriggerDecision * decision = dynamic_cast<AliHLTGlobalTriggerDecision*>(block->GetTObject());
297   if(decision) decision->Print();
298
299 }
300
301
302 ///______________________________________________________________________
303 void AliHLTEveHLT::ProcessEsdBlock( AliHLTHOMERBlockDesc * block) {
304   //See header file for documentation
305
306   AliESDEvent* esd = (AliESDEvent *) (block->GetTObject());
307   if (!esd) return;
308   
309   ProcessEsdEvent(esd);
310 }
311
312
313
314 ///_________________________________________________________________________________
315 Color_t AliHLTEveHLT::GetColor(Float_t pt) {
316   //See header file
317   Color_t baseColor = kOrange;
318   
319   Float_t binlimit = 0.1;
320   for(Int_t i = 0; i< 10; i++) {
321    
322     if (pt < binlimit)
323       return baseColor - i;
324     
325     binlimit +=0.1;
326   }
327   
328   return baseColor - 9;
329
330 }
331
332 ///_________________________________________________________________________________
333 Int_t AliHLTEveHLT::GetColorBin(Float_t pt) {
334   //See header file
335   
336   Float_t binlimit = 0.1;
337   for(Int_t i = 0; i< 10; i++) {
338    
339     if (pt < binlimit)
340       return i;
341     
342     binlimit +=0.1;
343   }
344   
345   return 9;
346   
347 }
348
349 ///____________________________________________________________________
350 void AliHLTEveHLT::ProcessEsdEvent( AliESDEvent * esd ) {
351   //See header file for documentation
352   if(!fPointSetVertex) CreateVertexPointSet();
353   if(!fTrackLists) CreateTrackLists();
354
355   cout << "ProcessESDEvent() :"<< esd->GetEventNumberInFile()<< "  " << esd->GetNumberOfCaloClusters() << " tracks : " << esd->GetNumberOfTracks() << endl;
356
357   //fEventManager->SetRunNumber(esd->GetRunNumber());
358
359   Double_t vertex[3];
360   const AliESDVertex * esdVertex = esd->GetPrimaryVertex();
361   
362   if(esdVertex) {
363     esdVertex->GetXYZ(vertex);
364     fPointSetVertex->SetNextPoint(vertex[0], vertex[1], vertex[2]);
365   }
366   
367   TEveTrackList * trackList = dynamic_cast<TEveTrackList*>(fTrackLists->FirstChild());
368   if(trackList) SetUpTrackPropagator(trackList->GetPropagator(),-0.1*esd->GetMagneticField(), 520);
369   for (Int_t iter = 0; iter < esd->GetNumberOfTracks(); ++iter) {
370
371    AliESDtrack * esdTrack = dynamic_cast<AliESDtrack*>(esd->GetTrack(iter));
372    FillTrackList(esdTrack);
373    FillHistograms(esdTrack);
374    
375   }
376   fHistMult->Fill(esd->GetNumberOfTracks()); // KK
377   
378   //BALLE hardcoded size
379   for(Int_t il = 0; il < fNTrackBins; il++) {
380     trackList = dynamic_cast<TEveTrackList*>(fTrackLists->FindChild(Form("Tracks_%d", il)));
381     trackList->MakeTracks();
382   } 
383 }
384 ///__________________________________________________________________________
385 void AliHLTEveHLT::FillTrackList(AliESDtrack * esdTrack) {
386   //See header file for documentation
387   TEveTrackList * trackList = dynamic_cast<TEveTrackList*>(fTrackLists->FirstChild());
388   if (!trackList) return;
389   
390   AliEveTrack* track = dynamic_cast<AliEveTrack*>(MakeEsdTrack(esdTrack, trackList));        
391   Int_t bin = GetColorBin(esdTrack->Pt());
392   trackList = dynamic_cast<TEveTrackList*>(fTrackLists->FindChild(Form("Tracks_%d", bin)));
393   if(trackList) {
394     track->SetAttLineAttMarker(trackList);
395     trackList->AddElement(track);
396   } else cout << "BALLE"<<endl;
397
398
399 }
400
401
402 ///____________________________________________________________________________________
403 void AliHLTEveHLT::FillHistograms(AliESDtrack * esdTrack) {
404
405   if(esdTrack->GetTPCNcls() == 0) return;
406   
407   fHistEta->Fill(esdTrack->Eta());
408   // fHistTheta->Fill(esdTrack->Theta()*TMath::RadToDeg());
409   fHistPhi->Fill(esdTrack->Phi()*TMath::RadToDeg());
410   
411   
412   Float_t DCAr, DCAz = -99;
413   esdTrack->GetImpactParametersTPC(DCAr, DCAz);
414   fHistDCAr->Fill(DCAr);
415   
416   
417   if(esdTrack->GetStatus()&AliESDtrack::kTPCin || (esdTrack->GetStatus()&AliESDtrack::kTPCin && esdTrack->GetStatus()&AliESDtrack::kITSin)){
418     fHistnClusters->Fill(esdTrack->GetTPCNcls());  
419   }
420 }
421
422 void AliHLTEveHLT::DrawHistograms(){
423   //See header file for documentation
424   if(!fCanvas) {
425     fCanvas = CreateCanvas("PVtx/Tr QA", "Primary vertex, Track QA");
426     fCanvas->Divide(3, 3);
427   }
428
429   fCanvas->cd(5);
430   fHistEta->Draw();
431
432   fCanvas->cd(6);
433   fHistPhi->Draw();
434
435   fCanvas->cd(7);
436   fHistnClusters->Draw();
437
438   fCanvas->cd(8);
439   fHistMult->Draw();
440
441   fCanvas->cd(9);
442   fHistDCAr->Draw();
443
444   fCanvas->cd();
445
446 }
447
448 AliEveTrack* AliHLTEveHLT::MakeEsdTrack (AliESDtrack *at, TEveTrackList* cont) {
449   //See header file for documentation
450   
451     
452
453
454
455   const double kCLight = 0.000299792458;
456   double bz = - kCLight*10.*( cont->GetPropagator()->GetMagField(0,0,0).fZ);
457
458   Bool_t innerTaken = kFALSE;
459   if ( ! at->IsOn(AliESDtrack::kITSrefit) && fUseIpOnFailedITS)
460   {
461     //tp = at->GetInnerParam();
462     innerTaken = kTRUE;
463   }
464
465   // Add inner/outer track parameters as path-marks.
466
467   Double_t     pbuf[3], vbuf[3];
468
469   AliExternalTrackParam trackParam = *at;
470
471   // take parameters constrained to vertex (if they are)
472
473   if( at->GetConstrainedParam() ){
474     trackParam = *at->GetConstrainedParam();
475   }
476   else if( at->GetInnerParam() ){
477     float x = at->GetX();
478     float y = at->GetY();
479     if(sqrt(x*x+y*y)>.5 ) trackParam = *(at->GetInnerParam());
480   }
481   if( at->GetStatus()&AliESDtrack::kTRDin ){
482     // transport to TRD in
483     trackParam = *at;
484     trackParam.PropagateTo( 290.45, -10.*( cont->GetPropagator()->GetMagField(0,0,0).fZ) );
485   }
486
487   TEveRecTrack rt;
488   {
489     rt.fLabel  = at->GetLabel();
490     rt.fIndex  = (Int_t) at->GetID();
491     rt.fStatus = (Int_t) at->GetStatus();
492     rt.fSign   = (Int_t) trackParam.GetSign();  
493     trackParam.GetXYZ(vbuf);
494     trackParam.GetPxPyPz(pbuf);    
495     rt.fV.Set(vbuf);
496     rt.fP.Set(pbuf);
497     Double_t ep = at->GetP(), mc = at->GetMass();
498     rt.fBeta = ep/TMath::Sqrt(ep*ep + mc*mc);
499   }
500
501   AliEveTrack* track = new AliEveTrack(&rt, cont->GetPropagator());
502   track->SetName(Form("AliEveTrack %d", at->GetID()));
503   track->SetElementTitle(CreateTrackTitle(at));
504   track->SetSourceObject(at);
505
506
507   // Set reference points along the trajectory
508   // and the last point
509
510   { 
511     TEvePathMark startPoint(TEvePathMark::kReference);
512     trackParam.GetXYZ(vbuf);
513     trackParam.GetPxPyPz(pbuf);    
514     startPoint.fV.Set(vbuf);
515     startPoint.fP.Set(pbuf);
516     rt.fV.Set(vbuf);
517     rt.fP.Set(pbuf);
518     Double_t ep = at->GetP(), mc = at->GetMass();
519     rt.fBeta = ep/TMath::Sqrt(ep*ep + mc*mc);
520
521     track->AddPathMark( startPoint );    
522   }
523
524   bool ok = 1;
525   if( at->GetOuterParam() && at->GetTPCPoints(2)>80 ){
526
527     //
528     // use AliHLTTPCCATrackParam propagator 
529     // since AliExternalTrackParam:PropagateTo()
530     // has an offset at big distances
531     //
532     double rot = at->GetOuterParam()->GetAlpha() - trackParam.GetAlpha();
533     double crot = cos(rot), srot = sin(rot);
534     double xEnd = at->GetTPCPoints(2)*crot -  at->GetTPCPoints(3)*srot;
535   // take parameters constrained to vertex (if they are)
536  
537
538     AliHLTTPCCATrackParam t;
539     AliHLTTPCCATrackConvertor::SetExtParam( t, trackParam );
540     
541     Double_t x0 = trackParam.GetX();
542     Double_t dx = xEnd - x0;
543     
544     if( dx<0 ) ok = 0;
545     //
546     // set a reference at the half of trajectory for better drawing
547     //
548     
549     if( ok ){ 
550       double dxx=dx/2; 
551       if( TMath::Abs(dxx)>=1. ){
552
553         if( !t.TransportToX(x0+dxx, bz, .999 ) ){
554           ok = 0;
555         } else {
556           AliExternalTrackParam tt;
557           AliHLTTPCCATrackConvertor::GetExtParam( t, tt, trackParam.GetAlpha() ); 
558           tt.GetXYZ(vbuf);
559           tt.GetPxPyPz(pbuf);
560           TEvePathMark midPoint(TEvePathMark::kReference);
561           midPoint.fV.Set(vbuf);
562           midPoint.fP.Set(pbuf);    
563           track->AddPathMark( midPoint );     
564         }
565       }
566     }
567    
568     //
569     // Set a reference at the end of the trajectory
570     // and a "decay point", to let the event display know where the track ends
571     //
572     
573     for( ; ok; dx*=.9 ){
574
575       // FIXME this loop does not make sense anymore
576       // Matthias 2010-05-02
577       if( !t.TransportToX(x0+dx, bz, .999 ) ){
578         ok = 0; 
579         break;
580         // if( TMath::Abs(dx)<1. ) break;      
581         // continue;
582       }
583       break;
584     }
585
586     {
587       if( !ok ){ 
588         AliHLTTPCCATrackConvertor::SetExtParam( t, trackParam );
589       }
590       AliExternalTrackParam tt;
591       AliHLTTPCCATrackConvertor::GetExtParam( t, tt, trackParam.GetAlpha() ); 
592       tt.GetXYZ(vbuf);
593       tt.GetPxPyPz(pbuf);
594       TEvePathMark endPoint(TEvePathMark::kReference);
595       TEvePathMark decPoint(TEvePathMark::kDecay);
596       endPoint.fV.Set(vbuf);
597       endPoint.fP.Set(pbuf);
598       decPoint.fV.Set(vbuf);
599       decPoint.fP.Set(pbuf);
600       track->AddPathMark( endPoint );
601       track->AddPathMark( decPoint );
602     }  
603   }
604
605   if ( ok &&  at->IsOn(AliESDtrack::kTPCrefit))
606   {
607     if ( ! innerTaken)
608     {
609       AddTrackParamToTrack(track, at->GetInnerParam());
610     }
611     AddTrackParamToTrack(track, at->GetOuterParam());
612   }
613   
614   return track;
615 }
616
617 void AliHLTEveHLT::SetUpTrackPropagator(TEveTrackPropagator* trkProp, Float_t magF, Float_t maxR) {
618   //See header file for documentation
619
620   if (fTrueField) {
621     trkProp->SetMagFieldObj(new AliEveMagField);
622   
623   } else {
624     trkProp->SetMagField(magF);
625   }
626  
627   if (fUseRkStepper) {
628     trkProp->SetStepper(TEveTrackPropagator::kRungeKutta);
629   }
630
631   trkProp->SetMaxR(maxR);
632 }
633
634
635 void AliHLTEveHLT::AddTrackParamToTrack(AliEveTrack* track, const AliExternalTrackParam* tp) {
636   //See header file for documentation
637
638   if (tp == 0)
639     return;
640
641   Double_t pbuf[3], vbuf[3];
642   tp->GetXYZ(vbuf);
643   tp->GetPxPyPz(pbuf);
644
645   TEvePathMark pm(TEvePathMark::kReference);
646   pm.fV.Set(vbuf);
647   pm.fP.Set(pbuf);
648   track->AddPathMark(pm);
649 }
650
651
652
653 TString AliHLTEveHLT::CreateTrackTitle(AliESDtrack* t) {
654   // Add additional track parameters as a path-mark to track.
655
656   TString s;
657
658   Int_t label = t->GetLabel(), index = t->GetID();
659   TString idx(index == kMinInt ? "<undef>" : Form("%d", index));
660   TString lbl(label == kMinInt ? "<undef>" : Form("%d", label));
661
662   Double_t p[3], v[3];
663   t->GetXYZ(v);
664   t->GetPxPyPz(p);
665   Double_t pt    = t->Pt();
666   Double_t ptsig = TMath::Sqrt(t->GetSigma1Pt2());
667   Double_t ptsq  = pt*pt;
668   Double_t ptm   = pt / (1.0 + pt*ptsig);
669   Double_t ptM   = pt / (1.0 - pt*ptsig);
670
671   s = Form("Index=%s, Label=%s\nChg=%d, Pdg=%d\n"
672            "pT = %.3f + %.3f - %.3f [%.3f]\n"
673            "P  = (%.3f, %.3f, %.3f)\n"
674            "V  = (%.3f, %.3f, %.3f)\n",
675            idx.Data(), lbl.Data(), t->Charge(), 0,
676            pt, ptM - pt, pt - ptm, ptsig*ptsq,
677            p[0], p[1], p[2],
678            v[0], v[1], v[2]);
679
680   Int_t   o;
681   s += "Det (in,out,refit,pid):\n";
682   o  = AliESDtrack::kITSin;
683   s += Form("ITS (%d,%d,%d,%d)  ",  t->IsOn(o), t->IsOn(o<<1), t->IsOn(o<<2), t->IsOn(o<<3));
684   o  = AliESDtrack::kTPCin;
685   s += Form("TPC(%d,%d,%d,%d)\n",   t->IsOn(o), t->IsOn(o<<1), t->IsOn(o<<2), t->IsOn(o<<3));
686   o  = AliESDtrack::kTRDin;
687   s += Form("TRD(%d,%d,%d,%d) ",    t->IsOn(o), t->IsOn(o<<1), t->IsOn(o<<2), t->IsOn(o<<3));
688   o  = AliESDtrack::kTOFin;
689   s += Form("TOF(%d,%d,%d,%d)\n",   t->IsOn(o), t->IsOn(o<<1), t->IsOn(o<<2), t->IsOn(o<<3));
690   o  = AliESDtrack::kHMPIDout;
691   s += Form("HMPID(out=%d,pid=%d)\n", t->IsOn(o), t->IsOn(o<<1));
692   s += Form("ESD pid=%d", t->IsOn(AliESDtrack::kESDpid));
693
694   if (t->IsOn(AliESDtrack::kESDpid))
695   {
696     Double_t pid[5];
697     t->GetESDpid(pid);
698     s += Form("\n[%.2f %.2f %.2f %.2f %.2f]", pid[0], pid[1], pid[2], pid[3], pid[4]);
699   }
700
701   return s;
702 }
703
704