]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveDet/AliEveTRDData.cxx
update for PWG1-EVE matching
[u/mrichter/AliRoot.git] / EVE / EveDet / AliEveTRDData.cxx
1 // $Id$
2 // Main authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007
3
4 /**************************************************************************
5  * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
6  * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for          *
7  * full copyright notice.                                                 *
8  **************************************************************************/
9
10 #include "TROOT.h"
11 #include "TStyle.h"
12 #include "TVector.h"
13 #include "TLinearFitter.h"
14 #include "TCanvas.h"
15 #include "TGeoMatrix.h"
16 #include "TGeoManager.h"
17
18 #include "TEveTrans.h"
19 #include "TEveManager.h"
20
21 #include "EveBase/AliEveEventManager.h"
22
23 #include "AliEveTRDData.h"
24 #include "AliEveTRDModuleImp.h"
25 #include "AliEveTRDLoader.h"
26 #include "AliEveTRDLoaderImp.h"
27
28 #include "AliGeomManager.h"
29 #include "AliESDtrack.h"
30 #include "AliLog.h"
31 #include "AliPID.h"
32 #include "AliTrackPointArray.h"
33 #include "AliRieman.h"
34
35 #include "AliTRDhit.h"
36 #include "AliTRDcluster.h"
37 #include "AliTRDseedV1.h"
38 #include "AliTRDtrackletMCM.h"
39 #include "AliTRDtrackletWord.h"
40 #include "AliTRDmcmSim.h"
41 #include "AliTRDtrackV1.h"
42 #include "AliTRDtrackerV1.h"
43 #include "AliTRDpadPlane.h"
44 #include "AliTRDdigitsManager.h"
45 #include "AliTRDmcmSim.h"
46 #include "AliTRDarrayADC.h"
47 #include "AliTRDSignalIndex.h"
48 #include "AliTRDgeometry.h"
49 #include "AliTRDtransform.h"
50 #include "AliTRDReconstructor.h"
51 #include "AliTRDrecoParam.h"
52
53 ClassImp(AliEveTRDHits)
54 ClassImp(AliEveTRDDigits)
55 ClassImp(AliEveTRDClusters)
56 ClassImp(AliEveTRDTracklet)
57 ClassImp(AliEveTRDTrack)
58 ClassImp(AliEveTRDTrackletOnline)
59 ClassImp(AliEveTRDmcm)
60
61 ///////////////////////////////////////////////////////////
62 /////////////   AliEveTRDDigits       /////////////////////
63 ///////////////////////////////////////////////////////////
64
65 //______________________________________________________________________________
66 AliEveTRDDigits::AliEveTRDDigits(AliEveTRDChamber *p) 
67   :TEveQuadSet("digits", "")
68   ,fParent(p)
69 {
70   // Constructor.
71   SetOwnIds(kTRUE);
72   gStyle->SetPalette(1, 0);
73   SetPalette(new TEveRGBAPalette(0, 512));
74   Reset(TEveQuadSet::kQT_RectangleYZ, kFALSE, 32);
75 }
76
77 //______________________________________________________________________________
78 AliEveTRDDigits::~AliEveTRDDigits()
79 {
80 //  AliInfo(GetTitle());
81 }
82
83 //______________________________________________________________________________
84 void AliEveTRDDigits::SetData(AliTRDdigitsManager *digits)
85 {
86   // Set data source.
87
88   Int_t det(fParent->GetID());
89   AliTRDarrayADC *data = digits->GetDigits(det);
90   if(!data->GetDim()) return;
91   data->Expand();
92
93   AliTRDSignalIndex *indexes = digits->GetIndexes(det);
94   if(!indexes->IsAllocated()) digits->BuildIndexes(det);
95
96   Double_t scale, dy, dz;
97   Int_t ly    = AliTRDgeometry::GetLayer(det),
98         stk   = AliTRDgeometry::GetStack(det),
99         sec   = AliTRDgeometry::GetSector(det),
100         vid   = AliGeomManager::LayerToVolUID(AliGeomManager::kTRD1 + ly, stk + AliTRDgeometry::Nstack() * sec);
101   SetNameTitle(Form("digits%03d", det), Form("D-%03d [%02d_%d_%d]", det, sec, stk, ly));
102   Short_t sig[7]={0,0,0,10,0,0,0};
103
104   AliTRDtransform transform(det);
105   AliTRDpadPlane *pp(fParent->fGeo->GetPadPlane(ly, stk));
106
107   Int_t row, col;
108   AliTRDcluster c;
109   indexes->ResetCounters();
110   while (indexes->NextRCIndex(row, col)){
111     dz = pp->GetRowSize(row);
112     dy = pp->GetColSize(col);
113     Short_t *const adc = data->GetDataAddress(row,col);
114     for (Int_t time(0); time<30; time++){     
115       if(data->IsPadCorrupted(row, col, time)) break;
116       if(adc[time] <= 1) continue;
117       new (&c) AliTRDcluster(det, col, row, time, sig, vid);
118       transform.Transform(&c);
119
120       scale = adc[time] < 512 ? adc[time]/512. : 1.;
121       AddQuad(c.GetY()-0.5*dy, c.GetZ()-0.5*dz*scale, c.GetX(), dy*0.95, dz*scale);
122       QuadValue(Int_t(adc[time]));
123       QuadId(new TNamed(Form("ADC %d", adc[time]), Form("det[%3d(%02d_%d_%d)] col[%3d] row[%2d] tb[%2d]", det, sec, stk, ly, col, row, time)));
124     } 
125   }
126
127   // rotate to global coordinates
128   RefitPlex();
129   TEveTrans& t = RefMainTrans();
130   t.SetRotByAngles((sec+.5)*AliTRDgeometry::GetAlpha(), 0.,0.);
131 }
132
133
134 // //______________________________________________________________________________
135 // void AliEveTRDDigits::Paint(Option_t *option)
136 // {
137 //   // Paint the object.
138 // 
139 //   if(fParent->GetDigitsBox()) fBoxes.Paint(option);
140 //   else TEveQuadSet::Paint(option);
141 // }
142
143 // //______________________________________________________________________________
144 // void AliEveTRDDigits::Reset()
145 // {
146 //   // Reset raw and visual data.
147 // 
148 //   TEveQuadSet::Reset(TEveQuadSet::kQT_RectangleYZ, kTRUE, 64);
149 //   // MT fBoxes.fBoxes.clear();
150 //   fData.Reset();
151 // }
152
153 ///////////////////////////////////////////////////////////
154 /////////////   AliEveTRDHits         /////////////////////
155 ///////////////////////////////////////////////////////////
156
157 //______________________________________________________________________________
158 AliEveTRDHits::AliEveTRDHits() : TEvePointSet("hits", 20)
159 {
160   // Constructor.
161   SetMarkerSize(.1);
162   SetMarkerColor(kGreen);
163   SetOwnIds(kTRUE);
164 }
165
166 //______________________________________________________________________________
167 AliEveTRDHits::~AliEveTRDHits()
168 {
169   //AliInfo(GetTitle());
170 }
171
172 //______________________________________________________________________________
173 void AliEveTRDHits::PointSelected(Int_t n)
174 {
175   // Handle an individual point selection from GL.
176
177   AliTRDhit *h = NULL;
178   if(!(h = dynamic_cast<AliTRDhit*>(GetPointId(n)))) return;
179   printf("Id[%3d] Det[%3d] Reg[%c] TR[%c] Q[%3d] MC[%d] t[%f]\n", 
180     n, h->GetDetector(), 
181     h->FromAmplification() ? 'A' : 'D', 
182     h->FromTRphoton() ? 'y' : 'n', 
183     h->GetCharge(), h->GetTrack(), h->GetTime());
184 }
185
186
187 ///////////////////////////////////////////////////////////
188 /////////////   AliEveTRDClusters         /////////////////////
189 ///////////////////////////////////////////////////////////
190
191 //______________________________________________________________________________
192 AliEveTRDClusters::AliEveTRDClusters():AliEveTRDHits()
193 {
194   // Constructor.
195   SetName("clusters");
196
197   SetMarkerSize(.4);
198   SetMarkerStyle(24);
199   SetMarkerColor(kGray);
200   SetOwnIds(kTRUE);
201 }
202
203 //______________________________________________________________________________
204 void AliEveTRDClusters::PointSelected(Int_t n)
205 {
206   // Handle an individual point selection from GL.
207
208   AliTRDcluster *c = dynamic_cast<AliTRDcluster*>(GetPointId(n));
209   if(!c) return;
210   c->Print();
211   Emit("PointSelected(Int_t)", n);
212   // Bool_t     AliCluster::GetGlobalCov(Float_t* cov) const
213   // Bool_t     AliCluster::GetGlobalXYZ(Float_t* xyz) const
214   // Float_t    AliCluster::GetSigmaY2() const
215   // Float_t    AliCluster::GetSigmaYZ() const
216   // Float_t    AliCluster::GetSigmaZ2() const
217 }
218
219 //______________________________________________________________________________
220 void AliEveTRDClusters::Print(Option_t *o) const
221 {
222   AliTRDcluster *c = NULL;
223
224   for(Int_t n = GetN(); n--;){
225     if(!(c = dynamic_cast<AliTRDcluster*>(GetPointId(n)))) continue;
226     c->Print(o);
227   }
228 }
229
230 //______________________________________________________________________________
231 void AliEveTRDClusters::Load(const Char_t *w) const
232 {
233   Int_t typ = -1;
234   if(strcmp(w, "hit")==0) typ = 0;
235   else if(strcmp(w, "dig")==0) typ = 1;
236   else if(strcmp(w, "cls")==0) typ = 2;
237   else if(strcmp(w, "tlt")==0) typ = 3;
238   else if(strcmp(w, "all")==0) typ = 4;
239   else{
240     AliInfo("The following arguments are accepted:");
241     AliInfo("   \"hit\" : loading of MC hits");
242     AliInfo("   \"dig\" : loading of digits");
243     AliInfo("   \"cls\" : loading of reconstructed clusters");
244     AliInfo("   \"tlt\" : loading of online tracklets");
245     AliInfo("   \"all\" : loading everything available");
246     return;
247   }
248
249   AliTRDcluster *c = NULL;
250   Int_t n = 0;
251   while((n = GetN() && !(c = dynamic_cast<AliTRDcluster*>(GetPointId(n))))) n++;
252   if(!c) return;
253
254   Int_t det = c->GetDetector();
255   AliEveTRDLoader *loader(NULL);
256   switch(typ){
257   case 4:
258     loader = new AliEveTRDLoaderSim("MC");
259     if(!loader->Open("galice.root")){delete loader; loader=NULL;}
260     else{
261       loader->SetDataType(AliEveTRDLoader::kTRDHits | AliEveTRDLoader::kTRDDigits | AliEveTRDLoader::kTRDClusters);
262       break;
263     }
264   case 0:  
265     loader = new AliEveTRDLoader("Hits");
266     if(!loader->Open("TRD.Hits.root")){delete loader; loader=NULL;}
267     else{
268       if(typ!=4) break;
269     }
270   case 1:
271     if(!loader) loader = new AliEveTRDLoader("Digits");
272     if(!loader->Open("TRD.Digits.root")){
273       if(typ==1){delete loader; loader=NULL;}
274     } else {
275       if(typ!=4) break;
276     }
277   case 2:
278     if(!loader) loader = new AliEveTRDLoader("Clusters");
279     if(!loader->Open("TRD.RecPoints.root")){
280       if(typ ==2){delete loader; loader=NULL;}
281     } else {
282       if(typ!=4) break;
283     }
284   case 3:
285     if(!loader) loader = new AliEveTRDLoader("Tracklets");
286     if(!loader->Open("TRD.Tracklets.root")){
287       if(typ ==3) {delete loader; loader=NULL;}
288     } else {
289       break;
290     }
291   default: return;
292   }
293   if(!loader) return;
294
295   loader->AddChambers(AliTRDgeometry::GetSector(det),AliTRDgeometry::GetStack(det), AliTRDgeometry::GetLayer(det));
296   // load first event
297   if(!loader->GoToEvent(AliEveEventManager::GetCurrent()->GetEventId())) AliWarning(Form("No data loaded for event %d", AliEveEventManager::GetCurrent()->GetEventId()));
298   // register loader with alieve
299   gEve->AddElement(loader->GetChamber(det), *(BeginParents()));
300   //loader->SpawnEditor();
301   gEve->Redraw3D();
302 }
303
304 ///////////////////////////////////////////////////////////
305 /////////////   AliEveTRDTracklet         /////////////////////
306 ///////////////////////////////////////////////////////////
307
308 //______________________________________________________________________________
309 AliEveTRDTracklet::AliEveTRDTracklet(AliTRDseedV1 *trklt):TEveLine()
310   ,fClusters(NULL)
311 {
312   // Constructor.
313   SetName("tracklet");
314   
315   if(!gGeoManager) AliEveEventManager::AssertGeometry();
316   SetUserData(trklt);
317   // init tracklet line
318   Int_t sec = AliTRDgeometry::GetSector(trklt->GetDetector());
319   Float_t alpha((0.5+sec)*AliTRDgeometry::GetAlpha()),
320           cphi(TMath::Cos(alpha)),
321           sphi(TMath::Sin(alpha));
322   Float_t dx;
323   Float_t x0   = trklt->GetX0();
324   Float_t y0   = trklt->GetYref(0);
325   Float_t z0   = trklt->GetZref(0);
326   Float_t dydx = trklt->GetYref(1);
327   Float_t dzdx = trklt->GetZref(1);
328   Float_t  tilt = trklt->GetTilt();
329   Float_t g[3];
330   AliTRDcluster *c = NULL;
331   for(Int_t ic=0; ic<AliTRDseedV1::kNclusters; ic++){
332     if(!(c = trklt->GetClusters(ic))) continue;
333     if(!fClusters) AddElement(fClusters = new AliEveTRDClusters());
334     dx = x0 - c->GetX();
335     //Float_t yt = y0 - dx*dydx;
336     Float_t zt = z0 - dx*dzdx;
337     // backup yc - for testing purposes
338     Float_t yc = c->GetY(); 
339     c->SetY(yc-tilt*(c->GetZ()-zt));
340     
341     // Ben: Temporary(?) fix -> Issue with wrongly set flag "fIsMisaligned" of AliCluster
342     // Changes can be undone ofter the issue with the flag has been solved
343     //c->GetGlobalXYZ(g);
344     
345     // The following code for the trasformation local-to-global is just adapted from
346     // AliTRDgeometry::RotateBack(...)
347   
348     // tracking to global coordinates transformation
349     Float_t loc[3] = { c->GetX(), c->GetY(), c->GetZ() }; 
350     g[0] = loc[0] * cphi - loc[1] * sphi;
351     g[1] = loc[0] * sphi + loc[1] * cphi;
352     g[2] = loc[2];
353     // Ben: End of fix
354     Int_t id = fClusters->SetNextPoint(g[0], g[1], g[2]);
355     c->SetY(yc);
356     fClusters->SetPointId(id, new AliTRDcluster(*c));
357   } 
358   if(fClusters){
359     fClusters->SetName("TC clusters");
360     fClusters->SetTitle(Form("N[%d]", trklt->GetN()));
361     fClusters->SetMarkerColor(kMagenta);
362   }
363
364   SetTitle(Form("Det[%d] RC[%c] Layer[%d] P[%7.3f]", trklt->GetDetector(), trklt->IsRowCross()?'y':'n', trklt->GetPlane(), trklt->GetMomentum()));
365   SetLineColor(kRed);
366   //SetOwnIds(kTRUE);
367   
368   //trklt->Fit(kTRUE);
369   y0   = trklt->GetYfit(0);
370   dydx = trklt->GetYfit(1);
371   Double_t xg(x0 * cphi - y0 * sphi),
372            yg(x0 * sphi + y0 * cphi);
373   SetPoint(0, xg, yg, z0); 
374   //SetPoint(0, x0, y0, z0);
375
376
377   dx = .5*AliTRDgeometry::CamHght()+AliTRDgeometry::CdrHght();
378   x0 -= dx; 
379   y0 -= dydx*dx,
380   z0 -= dzdx*dx; 
381   xg = x0 * cphi - y0 * sphi;
382   yg = x0 * sphi + y0 * cphi;
383   SetPoint(1, xg, yg, z0);
384   //SetPoint(1, x0, y0, z0);
385 }
386
387
388 AliEveTRDTracklet::~AliEveTRDTracklet() 
389 {
390
391 }
392
393 //______________________________________________________________________________
394 void AliEveTRDTracklet::Print(Option_t *o) const
395 {
396   AliTRDseedV1 *tracklet = (AliTRDseedV1*)GetUserData();
397   if(!tracklet) return;
398   tracklet->Print(o);
399 }
400
401 ///////////////////////////////////////////////////////////
402 /////////////   AliEveTRDTrack         /////////////////////
403 ///////////////////////////////////////////////////////////
404
405 //______________________________________________________________________________
406 AliEveTRDTrack::AliEveTRDTrack(AliTRDtrackV1 *trk) 
407   :TEveLine()
408   ,fTrackState(0)
409   ,fESDStatus(0)
410   ,fAlpha(0.)
411   ,fPoints(NULL)
412   ,fRim(NULL)
413 {
414   // Constructor.
415   SetUserData(trk);
416   SetName("");
417
418 //  AliTRDtrackerV1::SetNTimeBins(24);
419
420   fRim = new AliRieman(trk->GetNumberOfClusters());
421   AliTRDseedV1 *tracklet = NULL;
422   for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++){
423     if(!(tracklet = trk->GetTracklet(il))) continue;
424     if(!tracklet->IsOK()) continue;
425     AddElement(new AliEveTRDTracklet(tracklet));
426
427 //     tracklet->ResetClusterIter(kFALSE);
428 //     while((c = tracklet->PrevCluster())){
429 //    AliTRDcluster *c(NULL);
430 /*    for(Int_t ic=AliTRDseedV1::kNtb; ic--;){
431       if(!(c=tracklet->GetClusters(ic))) continue;
432       Float_t xc = c->GetX();
433       Float_t yc = c->GetY();
434       Float_t zc = c->GetZ();
435       Float_t zt = tracklet->GetZref(0) - (tracklet->GetX0()-xc)*tracklet->GetZref(1); 
436       yc -= tracklet->GetTilt()*(zc-zt);
437       fRim->AddPoint(xc, yc, zc, .05, 2.3);
438     }*/
439   }
440   if(trk->GetNumberOfTracklets()>1) fRim->Update();
441   SetStatus(fTrackState);
442 }
443
444 //______________________________________________________________________________
445 AliEveTRDTrack::~AliEveTRDTrack()
446 {
447   if(fPoints) delete [] fPoints; fPoints = NULL;
448   //delete dynamic_cast<AliTRDtrackV1*>(GetUserData());
449 }
450
451 //______________________________________________________________________________
452 void AliEveTRDTrack::Print(Option_t *o) const
453 {
454   AliTRDtrackV1 *track = (AliTRDtrackV1*)GetUserData();
455   if(!track) return;
456   track->Print(o);
457 }
458
459 //______________________________________________________________________________
460 void AliEveTRDTrack::SetStatus(UChar_t s)
461 {
462   // nothing to be done
463   if(fPoints && fTrackState == s) return;
464   //return;
465   const Int_t nc = AliTRDtrackV1::kMAXCLUSTERSPERTRACK;
466   AliTRDtrackV1 *trk(NULL);
467   if(!(trk=static_cast<AliTRDtrackV1*>(GetUserData()))) {
468     AliError("Failed casting data to TRD track.");
469     return;
470   }
471
472   Bool_t BUILD = kFALSE;
473   if(!fPoints){ 
474     fPoints = new AliTrackPoint[nc];
475
476     // define the radial span of the track in the TRD
477     Double_t xmin = -1., xmax = -1.;
478     Int_t det = 0;
479     AliTRDseedV1 *trklt = NULL;
480     for(Int_t ily=0; ily<AliTRDgeometry::kNlayer; ily++){
481       if(!(trklt = trk->GetTracklet(ily))) continue;
482       if(xmin<0.) xmin = trklt->GetX0() - AliTRDgeometry::CamHght() - AliTRDgeometry::CdrHght();
483       if(trklt->GetX0()>xmax) xmax = trklt->GetX0();
484       det = trklt->GetDetector();
485     }
486     Int_t sec = det/AliTRDgeometry::kNdets;
487     fAlpha = AliTRDgeometry::GetAlpha() * (sec<9 ? sec + .5 : sec - 17.5); //trk->GetAlpha()
488
489     Double_t dx =(xmax - xmin)/nc;
490     for(Int_t ip=0; ip<nc; ip++){
491       fPoints[ip].SetXYZ(xmin, 0., 0.);
492       xmin+=dx;
493     }
494     BUILD = kTRUE;
495   }
496
497   // select track model
498   if(BUILD || ((s&12) != (fTrackState&12))){
499     if(TESTBIT(s, kTrackCosmics)){
500       //printf("Straight track\n");
501       AliTRDtrackerV1::FitLine(trk, NULL, kFALSE, nc, fPoints);
502     } else {
503       if(TESTBIT(s, kTrackModel)){
504         //printf("Kalman track\n");
505         if(trk->GetNumberOfTracklets() >=4) AliTRDtrackerV1::FitKalman(trk, NULL, kFALSE, nc, fPoints);
506       } else { 
507         //printf("Rieman track rim[%p] nc[%d]\n", (void*)fRim, nc);
508         //if(trk->GetNumberOfTracklets() >=4) AliTRDtrackerV1::FitRiemanTilt(trk, NULL, kTRUE, nc, fPoints);
509         Float_t x = 0.;
510         for(Int_t ip = nc; ip--;){
511           x = fPoints[ip].GetX();
512           //printf("%2d x[%f] y[%f] z[%f]\n", ip, x, fRim->GetYat(x), fRim->GetZat(x));
513           fPoints[ip].SetXYZ(x, fRim->GetYat(x), fRim->GetZat(x));
514         }
515       }
516     }
517   
518     Float_t global[3];
519     for(Int_t ip=0; ip<nc; ip++){
520       fPoints[ip].Rotate(-fAlpha).GetXYZ(global);
521       SetPoint(ip, global[0], global[1], global[2]);
522       //printf("*** %2d x[%f] y[%f] z[%f]\n", ip, global[0], global[1], global[2]);
523     
524     }
525     SetSmooth(kTRUE);
526   }
527
528   // set color
529   if(BUILD || ((s&3) != (fTrackState&3))){
530     if(TESTBIT(s, kSource)){
531       //printf("Source color\n");
532       if(fESDStatus&AliESDtrack::kTRDin){
533         SetMarkerColor(kGreen);
534         SetLineColor(kGreen);
535       } else {
536         SetMarkerColor(kMagenta);
537         SetLineColor(kMagenta);
538       }
539     } else {
540       if(TESTBIT(s, kPID) == AliTRDpidUtil::kLQ){
541         //printf("PID color kLQPID\n");
542         //trk->GetReconstructor()->SetOption("!nn");
543       } else {
544         //printf("PID color kNNPID\n");
545         //trk->GetReconstructor()->SetOption("nn");
546       }
547       //trk->CookPID();
548   
549       Int_t species = 0; Float_t pid = 0.;
550       for(Int_t is=0; is<AliPID::kSPECIES; ++is) 
551         if(trk->GetPID(is) > pid){
552           pid = trk->GetPID(is);
553                   species = (AliPID::EParticleType) is;
554         }
555       switch(species){
556       case AliPID::kElectron:
557         SetMarkerColor(kRed);
558         SetLineColor(kRed);
559         break;
560       default:
561         SetMarkerColor(kBlue);
562         SetLineColor(kBlue);
563         break;
564       }
565     }
566     SetLineWidth(2);
567   }
568   
569   const Char_t *model = "line";
570   if(!TESTBIT(s, kTrackCosmics)){
571     if(TESTBIT(s, kTrackModel)) model = "kalman";
572     else model = "rieman";
573   }
574   Int_t species = 0; Float_t pid = 0.;
575   for(Int_t is=0; is<AliPID::kSPECIES; is++) 
576     if(trk->GetPID(is) > pid){
577       pid = trk->GetPID(is);
578       species = is;
579     }
580
581   SetTitle(Form(
582     "Tracklets[%d] Clusters[%d]\n"
583     "Reconstruction Source[%s]\n"
584     "PID[%4.1f %4.1f %4.1f %4.1f %4.1f]\n"
585     "MC[%d]", trk->GetNumberOfTracklets(), trk->GetNumberOfClusters(), fESDStatus&AliESDtrack::kTRDin ? "barrel" : "sa",
586     1.E2*trk->GetPID(0), 1.E2*trk->GetPID(1),
587     1.E2*trk->GetPID(2), 1.E2*trk->GetPID(3), 1.E2*trk->GetPID(4), trk->GetLabel()));
588
589   if(GetName()){
590     char id[6]; snprintf(id, 6, "%s", GetName());
591     SetName(Form("%s %s", id, AliPID::ParticleName(species)));
592   }
593
594   // save track status
595   fTrackState = s;
596 }
597
598 //______________________________________________________________________________
599 void AliEveTRDTrack::Load(const Char_t *what) const
600 {
601 // Spread downwards to tracklets the command "what"
602
603   const AliEveTRDTracklet* trklt(NULL);
604   TEveElement::List_ci itrklt=BeginChildren();
605   while(itrklt!=EndChildren()){
606     if((trklt = dynamic_cast<const AliEveTRDTracklet*>(*itrklt))) trklt->Load(what);
607     itrklt++;
608   }
609 }
610
611 //______________________________________________________________________________
612 AliEveTRDTrackletOnline::AliEveTRDTrackletOnline(AliTRDtrackletMCM *tracklet) :
613   TEveLine(),
614   fDetector(-1),
615   fROB(-1),
616   fMCM(-1)
617 {
618   AliTRDtrackletMCM *trkl = new AliTRDtrackletMCM(*tracklet);
619   SetUserData(trkl);
620
621   fDetector = trkl->GetDetector();
622   fROB = trkl->GetROB();
623   fMCM = trkl->GetMCM();
624
625   SetName("sim. tracklet");
626   SetTitle(Form("Det: %i, ROB: %i, MCM: %i, Label: %i\n0x%08x", 
627                 trkl->GetDetector(), trkl->GetROB(), trkl->GetMCM(), trkl->GetLabel(),
628                 trkl->GetTrackletWord()));
629   SetLineColor(kGreen);
630   SetLineWidth(3);
631
632   AliTRDgeometry geo;
633   TGeoHMatrix *matrix = geo.GetClusterMatrix(trkl->GetDetector());
634
635   fDetector = trkl->GetDetector();
636   fROB = trkl->GetROB();
637   fMCM = trkl->GetMCM();
638   
639   Float_t length = 3.;
640   Double_t x[3];
641   Double_t p[3];
642   Double_t p2[3];
643   x[0] = AliTRDgeometry::AnodePos(); 
644   x[1] = trkl->GetY();
645   x[2] = trkl->GetLocalZ();
646
647   matrix->LocalToMaster(x, p);
648   geo.RotateBack(trkl->GetDetector(), p, p2);
649   SetPoint(0, p2[0], p2[1], p2[2]);
650
651   x[0] -= length;
652   x[1] -= length * trkl->GetdYdX();
653   matrix->LocalToMaster(x, p);
654   p[2] *= p[0] / (p[0] + length);
655   geo.RotateBack(trkl->GetDetector(), p, p2);
656   SetPoint(1, p2[0], p2[1], p2[2]);
657 }
658
659 AliEveTRDTrackletOnline::AliEveTRDTrackletOnline(AliTRDtrackletWord *tracklet) :
660   TEveLine(),
661   fDetector(-1),
662   fROB(-1),
663   fMCM(-1)
664 {
665   AliTRDtrackletWord *trkl = new AliTRDtrackletWord(*tracklet);
666   SetUserData(trkl);
667
668   fDetector = trkl->GetDetector();
669   fROB = trkl->GetROB(); 
670   fMCM = trkl->GetMCM(); 
671
672   SetName("raw tracklet");
673   SetTitle(Form("Det: %i, ROB: %i, MCM: %i, Label: %i\n0x%08x", 
674                 trkl->GetDetector(), fROB, fMCM, -1,
675                 trkl->GetTrackletWord()));
676   SetLineColor(kRed);
677   SetLineWidth(3);
678
679   AliTRDgeometry geo;
680   TGeoHMatrix *matrix = geo.GetClusterMatrix(trkl->GetDetector());
681
682   Float_t length = 3.;
683   Double_t x[3];
684   Double_t p[3];
685   Double_t p2[3];
686   x[0] = AliTRDgeometry::AnodePos();
687   x[1] = trkl->GetY();
688   x[2] = trkl->GetLocalZ();
689   
690   matrix->LocalToMaster(x, p);
691   geo.RotateBack(trkl->GetDetector(), p, p2);
692   SetPoint(0, p2[0], p2[1], p2[2]);
693
694   x[0] -= length;
695   x[1] -= length * trkl->GetdYdX();
696   matrix->LocalToMaster(x, p);
697   p[2] *= p[0] / (p[0] + length);
698   geo.RotateBack(trkl->GetDetector(), p, p2);
699   SetPoint(1, p2[0], p2[1], p2[2]);
700 }
701
702 AliEveTRDTrackletOnline::~AliEveTRDTrackletOnline() 
703 {
704   delete ((AliTRDtrackletBase*) GetUserData());
705   SetUserData(0x0);
706 }
707
708 void AliEveTRDTrackletOnline::ShowMCM(Option_t *opt) const
709 {
710   if (fDetector < 0 || fROB < 0 || fMCM < 0)
711     return;
712
713   AliEveTRDmcm *evemcm = new AliEveTRDmcm();
714   evemcm->Init(fDetector, fROB, fMCM);
715   evemcm->LoadDigits();
716   evemcm->Draw(opt);
717
718   TEveElementList *mcmlist = NULL;
719   if (gEve->GetCurrentEvent()) 
720     mcmlist = (TEveElementList*) gEve->GetCurrentEvent()->FindChild("TRD MCMs");
721   if (!mcmlist) {
722     mcmlist = new TEveElementList("TRD MCMs");
723     gEve->AddElement(mcmlist);
724   }
725   gEve->AddElement(evemcm, mcmlist);
726 }
727
728
729 AliEveTRDmcm::AliEveTRDmcm() :
730   TEveElement(),
731   TNamed(),
732   fMCM(new AliTRDmcmSim())
733 {
734   SetName("MCM");
735   SetTitle("Unknown MCM");
736 }
737
738 AliEveTRDmcm::~AliEveTRDmcm()
739 {
740   delete fMCM;
741 }
742
743 Bool_t AliEveTRDmcm::Init(Int_t det, Int_t rob, Int_t mcm)
744 {
745   SetName(Form("MCM: Det. %i, ROB %i, MCM %i", det, rob, mcm));
746   SetTitle(Form("MCM: Det. %i, ROB %i, MCM %i", det, rob, mcm));
747   fMCM->Init(det, rob, mcm);
748   fMCM->Reset();
749   return kTRUE;
750 }
751
752 Bool_t AliEveTRDmcm::LoadDigits()
753 {
754   AliRunLoader *rl = AliEveEventManager::AssertRunLoader();
755   return fMCM->LoadMCM(rl, fMCM->GetDetector(), 
756                        fMCM->GetRobPos(), fMCM->GetMcmPos());
757 }
758
759 Bool_t AliEveTRDmcm::Filter()
760 {
761   fMCM->Filter();
762   return kTRUE;
763 }
764
765 Bool_t AliEveTRDmcm::Tracklet()
766 {
767   fMCM->Tracklet();
768   return kTRUE;
769 }
770
771 void AliEveTRDmcm::Draw(Option_t* option)
772 {
773   const char *mcmname = Form("mcm_%i_%i_%i", fMCM->GetDetector(),
774                              fMCM->GetRobPos(), fMCM->GetMcmPos());
775
776   TCanvas *c = dynamic_cast<TCanvas*> (gROOT->FindObject(mcmname));
777   if (!c)
778     c = gEve->AddCanvasTab("TRD MCM");
779   c->SetTitle(Form("MCM %i on ROB %i of det. %i", 
780                    fMCM->GetMcmPos(), fMCM->GetRobPos(), fMCM->GetDetector()));
781   c->SetName(mcmname);
782   c->cd();
783   fMCM->Draw(option);
784 }
785
786 Bool_t AliEveTRDmcm::AssignPointer(const char* ptrname)
787 {
788   gROOT->ProcessLine(Form("AliTRDmcmSim* %s = (AliTRDmcmSim *)%p", ptrname, (void*)fMCM));
789   return kTRUE;
790 }