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