]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/OnlineDisplay/AliHLTTPCDisplay.cxx
get rid of coding violations and warnings (Jochen)
[u/mrichter/AliRoot.git] / HLT / TPCLib / OnlineDisplay / AliHLTTPCDisplay.cxx
1 // @(#) $Id$
2 // Original: AliHLTDisplay.cxx,v 1.26 2005/06/14 10:55:21 cvetan 
3
4 //*************************************************************************
5 // This file is property of and copyright by the ALICE HLT Project        * 
6 // ALICE Experiment at CERN, All rights reserved.                         *
7 //                                                                        *
8 // Primary Authors: Jochen Thaeder <mailto:thaeder@kip.uni-heidelberg.de> *
9 //                  Anders Vestbo                                         *
10 //                  for The ALICE HLT Project.                            *
11 //                                                                        *
12 // Permission to use, copy, modify and distribute this software and its   *
13 // documentation strictly for non-commercial purposes is hereby granted   *
14 // without fee, provided that the above copyright notice appears in all   *
15 // copies and that both the copyright notice and this permission notice   *
16 // appear in the supporting documentation. The authors make no claims     *
17 // about the suitability of this software for any purpose. It is          *
18 // provided "as is" without express or implied warranty.                  *
19 //*************************************************************************/
20
21 /** @file   AliHLTTPCDisplay.cxx
22     @author Jochen Thaeder, Anders Vestbo
23     @date   
24     @brief  Display class for the HLT TPC events.
25 */
26
27 #define TRACKHELIX 0
28 #define TRACKPOLYMARKER 0
29 #define BACKWARD 0
30 #define FIRSTLASTPOINT 0
31
32 #define TRACKCOLOR 
33 #define USEDCLUSTERCOLOR
34 #define UNUSEDCLUSTERCOLOR
35
36 #include "AliHLTStdIncludes.h"
37 #if defined(HAVE_TVIEW3D_H)
38 #include <TView3D.h>
39 #else
40 #include <TView.h>
41 #endif
42 #include <TPolyMarker3D.h>
43 #include <TPolyLine3D.h>
44 #include <TH2.h>
45 #include <TTree.h>
46 #include <TNode.h>
47 #include <TGeometry.h>
48 #include <TShape.h>
49 #include <TParticle.h>
50 #include <TFile.h>
51 #include <THelix.h>
52 #include <TStyle.h>
53 #include <TGraph.h>
54 #include <TMultiGraph.h>
55 #include <TAttText.h>
56 #include <TAxis.h>
57
58 #if TRACKHELIX
59 #include <THelix.h>
60 #endif
61
62 #ifdef use_aliroot
63 #include <TClonesArray.h>
64 #include <AliRun.h>
65 #include <AliSimDigits.h>
66 #include <AliTPCParam.h>
67 #endif
68
69 #include "AliHLTTPCLogging.h"
70 #include "AliHLTTPCDisplay.h"
71 #include "AliHLTTPCTransform.h"
72 #include "AliHLTTPCTrack.h"
73 #include "AliHLTTPCTrackArray.h"
74 #include "AliHLTTPCSpacePointData.h"
75 #include "AliHLTTPCMemHandler.h"
76 #include "AliHLTTPCDigitReaderPacked.h"
77
78 #if __GNUC__ == 3
79 using namespace std;
80 #endif
81
82 ClassImp(AliHLTTPCDisplay);
83
84 AliHLTTPCDisplay::AliHLTTPCDisplay( Char_t * gfile ) :
85   fTrackParam(),
86   fClusters(),
87   fTracks(),
88   fNcl(),
89   fHistrawcl(NULL),
90   fHistraw(NULL),
91   fHistpad1(NULL),
92   fHistpad2(NULL),
93   fHistpad3(NULL),
94   fHistallresidualsY(NULL),
95   fHistallresidualsZ(NULL),
96   fHistcharge(NULL),
97   fGraphresidualsY(NULL),
98   fGraphresidualsZ(NULL),
99   fGraphresidualsYLength(NULL),
100   fGraphresidualsZLength(NULL),
101   fGeom(NULL),
102   fBackColor(1),
103   fLineColor(0),
104   fKeepView(kFALSE),
105   fPad(-1),
106   fPadRow(0),
107   fSlicePadRow(0),
108   fNPads(0),
109   fNTimes(0),
110   fMinHits(0),
111   fPtThreshold(0),
112   fSelectTrackSwitch(kFALSE),
113   fSelectTrack(-1),
114   fSelectTrackSlice(0),
115   fSelectCluster(0),
116   fMinSlice(0),
117   fMaxSlice(25),
118   fSlicePair(kFALSE),
119   fSliceArray(),
120   fDrawGeo(kFALSE),
121   fcolorbin(),
122   fbinct(),
123   fpmarr(),
124   fSwitch3DCluster(kFALSE),
125   fSwitch3DTracks(kFALSE),
126   fSwitch3DPadRow(kFALSE),
127   fSwitch3DGeometry(kFALSE) 
128 {
129   // constructor    
130   InitDisplay(gfile);
131 }
132
133 // #############################################################################
134 void AliHLTTPCDisplay::InitDisplay(Char_t *gfile) {
135     //constructor
136     memset(fClusters,0,36*6*sizeof(AliHLTTPCSpacePointData*));
137     memset(fNcl, 0, 36*6*sizeof(UInt_t)); 
138    
139 // ---------------------------------------------------
140 // In order to be backward compatible
141 // ---------------------------------------------------
142 #if BACKWARD
143     //fc1 = NULL;
144 #endif 
145 // ---------------------------------------------------
146    
147     SetSliceArray();
148    
149     AliHLTTPCTransform::SetBField(0.4);
150     
151     LoadGeometrie(gfile);
152 }
153
154
155 // #############################################################################
156 AliHLTTPCDisplay::~AliHLTTPCDisplay() {
157     //destructor
158     if(fTracks) delete fTracks;
159     fTracks = NULL;
160 }
161
162 // #############################################################################
163 Bool_t AliHLTTPCDisplay::LoadGeometrie(Char_t *gfile) {
164     if (gfile) {
165         TFile *file = TFile::Open(gfile);
166         if(!file) {
167             LOG(AliHLTTPCLog::kError,"AliHLTTPCDisplay::AliHLTTPCDisplay","File Open") <<"Geometry file " << gfile << " does not exist!"<<ENDLOG;
168             return kFALSE;
169         }
170         
171         fGeom = (TGeometry*)file->Get("AliceGeom");
172
173         file->Close();
174         delete file;
175     }
176     return kTRUE;
177 }
178
179 // #############################################################################
180 //                 EXECUTER
181 // #############################################################################
182 void AliHLTTPCDisplay::ExecPadRow(){
183    int event = gPad->GetEvent();
184    if (event != 11) return;
185
186    printf("TEST !!!!!!!!!!!!!!!");
187 /*   int px = gPad->GetEventX();
188    TObject *select = gPad->GetSelected();
189    if (!select) return;
190    if (select->InheritsFrom("TH1")) {
191       TH1 *h = (TH1*)select;
192       Float_t xx = gPad->AbsPixeltoX(px);
193       Float_t x  = gPad->PadtoX(xx);
194       Int_t binx = h->GetXaxis()->FindBin(x);
195       printf("event=%d, hist:%s, bin=%d, content=%f\n",event,h->GetName(),binx,h->GetBinContent(binx));
196    }
197
198 */
199
200 }
201
202 // #############################################################################
203 //                 SETTER
204 // #############################################################################
205 void AliHLTTPCDisplay::SetHistPadRowAxis() {
206     // Set Axis range of Histogramm, due to variable NPads per padrow
207
208     fNPads = AliHLTTPCTransform::GetNPads(fPadRow);
209     fHistrawcl->SetAxisRange(0,fNPads);
210     fHistraw->SetAxisRange(0,fNPads);
211     fHistrawcl->SetAxisRange(0,fNTimes,"Y");
212     fHistraw->SetAxisRange(0,fNTimes,"Y");
213 }
214
215 void AliHLTTPCDisplay::SetSliceArray() {
216     Int_t slice=0;
217     Int_t minSlice = fMinSlice; 
218     Int_t maxSlice = fMaxSlice; 
219     Int_t realslice = 0;
220
221     for (slice=0;slice<=35;slice++){
222         fSliceArray[slice] = kFALSE;
223     }
224
225     // Single Slice, or Range
226     if (minSlice > maxSlice) maxSlice += 17;
227         
228     for (slice=minSlice;slice<=maxSlice;slice++){
229         realslice = slice % 18;
230         fSliceArray[realslice] = kTRUE;
231         fSliceArray[realslice+18] = kTRUE;
232     }
233
234     // Pair of Slices
235     if (fSlicePair) {
236         minSlice = fMinSlice + 9;
237         maxSlice = fMaxSlice + 9;
238         
239         if (minSlice > maxSlice) maxSlice += 17;
240         
241         for (slice=minSlice;slice<=maxSlice;slice++){
242             realslice = slice % 18;
243             fSliceArray[realslice] = kTRUE;     
244             fSliceArray[realslice+18] = kTRUE;
245         }
246     }
247 }
248
249 // #############################################################################
250 //                 SETUP
251 // #############################################################################
252 void AliHLTTPCDisplay::SetupCluster(Int_t slice, Int_t patch, UInt_t nofClusters, AliHLTTPCSpacePointData* data)  {  
253
254     if (data && slice>=0 && slice<36 && patch>=0 && patch<AliHLTTPCTransform::GetNPatches()) {
255         if (fClusters[slice][patch]!=NULL) {
256             delete(fClusters[slice][patch]);
257             fClusters[slice][patch]=NULL;
258         }
259         Int_t arraysize=nofClusters*sizeof(AliHLTTPCSpacePointData);
260         fClusters[slice][patch] = (AliHLTTPCSpacePointData*)new Byte_t[arraysize];
261         if (fClusters[slice][patch]) {
262             memcpy(fClusters[slice][patch], data, arraysize);
263             fNcl[slice][patch]=nofClusters;
264         } else {
265             fNcl[slice][patch]=nofClusters;
266             LOG(AliHLTTPCLog::kError,"AliHLTTPCDisplay::SetupCluster","memory allocation") << "memory allocation failed "<<ENDLOG; 
267         }
268     } else LOG(AliHLTTPCLog::kError,"AliHLTTPCDisplay::SetupCluster","argument check") << "invalid argument "<<ENDLOG; 
269 }
270
271 // #############################################################################
272 void AliHLTTPCDisplay::SetupTracks(AliHLTTPCTrackArray *tracks) {
273     fTracks=tracks;
274
275     // Set USED cluster
276     Int_t ntracks = fTracks->GetNTracks();
277
278     for(Int_t j=0; j<ntracks; j++) {    
279         AliHLTTPCTrack *gtrack = fTracks->GetCheckedTrack(j); 
280         if(!gtrack) continue;
281
282         Int_t nHits = gtrack->GetNHits();
283         UInt_t *hitnum = gtrack->GetHitNumbers();
284
285         for(Int_t h=0; h<nHits; h++){
286           
287             UInt_t id=hitnum[h];
288             Int_t slice = (id>>25) & 0x7f;
289             Int_t patch = (id>>22) & 0x7;
290             UInt_t pos = id&0x3fffff;       
291                 
292             AliHLTTPCSpacePointData *points = fClusters[slice][patch];
293             
294             if(!points) {
295                 LOG(AliHLTTPCLog::kError,"AliHLTTPCDisplay::Draw3D","Clusterarray") <<"No points at slice "<<slice<<" patch "<<patch<<" pos "<<pos<<ENDLOG;
296                 continue;
297             }
298  
299             if(pos>=fNcl[slice][patch]) {
300                 LOG(AliHLTTPCLog::kError,"AliHLTTPCDisplay::Draw3D","Clusterarray") <<"Pos is too large: pos "<<pos <<" ncl "<<fNcl[slice][patch]<<ENDLOG;
301                 continue;
302             }
303             points[pos].fUsed = kTRUE;
304             points[pos].fTrackN = j;
305         }
306     }
307 }
308
309 // #############################################################################
310 void AliHLTTPCDisplay::SetupHist(){
311
312     Int_t maxpads = 150;
313     fNTimes = AliHLTTPCTransform::GetNTimeBins();
314
315     if ( fHistraw ){
316         delete fHistraw;
317         fHistraw = NULL;
318     }
319     if ( fHistrawcl ){
320         delete fHistrawcl;
321         fHistrawcl = NULL;
322     }
323
324     if ( fHistpad1 ){
325         delete fHistpad1;
326         fHistpad1 = NULL;
327     }
328
329     if ( fHistpad2 ){
330         delete fHistpad2;
331         fHistpad2 = NULL;
332     }
333
334     if ( fHistpad3 ){
335         delete fHistpad3;
336         fHistpad3 = NULL;
337     }
338
339     if ( fHistallresidualsY){
340         delete fHistallresidualsY;
341         fHistallresidualsY = NULL;
342     }
343
344     if ( fHistallresidualsZ){
345         delete fHistallresidualsZ;
346         fHistallresidualsZ = NULL;
347     }
348     if ( fHistcharge){
349         delete fHistcharge;
350         fHistcharge = NULL;
351     }
352
353     // Setup the histograms
354     Int_t padbinning = maxpads*10;
355     fHistraw = new TH2F("fHistraw","Selected PadRow with found Clusters;Pad #;Timebin #",maxpads,0,maxpads-1,fNTimes,0,fNTimes-1);
356     fHistrawcl = new TH1F("fHistrawcl","",padbinning,0,maxpads-1);
357     fHistpad1 = new TH1F ("fHistpad1","Selected Pad -1;Timebin #",fNTimes,0,fNTimes-1);
358     fHistpad2 = new TH1F ("fHistpad2","Selected Pad;Timebin #",fNTimes,0,fNTimes-1); 
359     fHistpad3 = new TH1F ("fHistpad3","Selected Pad +1;Timebin #",fNTimes,0,fNTimes-1);
360     fHistallresidualsY = new TH1F ("fHistallresiduals","Y Residuals of all Tracks in selected slices;residuals",5000,-100,100);
361     fHistallresidualsZ = new TH1F ("fHistallresiduals","Z Residuals of all Tracks in selected slices;residuals",5000,-100,100);
362     fHistcharge = new TH1F ("fHistcharge","Cluster distribution per charge;charge;#cluster",5000,0,30000);
363
364     fHistraw->SetOption("COLZ"); 
365
366     fHistallresidualsY->SetTitleSize(0.03);
367     fHistallresidualsY->GetXaxis()->SetLabelSize(0.03);
368     fHistallresidualsY->GetXaxis()->SetTitleSize(0.03);
369     fHistallresidualsY->GetYaxis()->SetLabelSize(0.03);
370     fHistallresidualsY->GetYaxis()->SetTitleSize(0.03);
371
372     fHistallresidualsZ->SetTitleSize(0.03);
373     fHistallresidualsZ->GetXaxis()->SetLabelSize(0.03);
374     fHistallresidualsZ->GetXaxis()->SetTitleSize(0.03);
375     fHistallresidualsZ->GetYaxis()->SetLabelSize(0.03);
376     fHistallresidualsZ->GetYaxis()->SetTitleSize(0.03);
377
378     fHistcharge->SetTitleSize(0.03);
379     fHistcharge->GetXaxis()->SetLabelSize(0.03);
380     fHistcharge->GetXaxis()->SetTitleSize(0.03);
381     fHistcharge->GetYaxis()->SetLabelSize(0.03);
382     fHistcharge->GetYaxis()->SetTitleSize(0.03);
383
384     fHistraw->SetTitleSize(0.03);
385     fHistraw->GetXaxis()->SetLabelSize(0.03);
386     fHistraw->GetXaxis()->SetTitleSize(0.03);
387     fHistraw->GetYaxis()->SetLabelSize(0.03);
388     fHistraw->GetYaxis()->SetTitleSize(0.03);
389
390     fHistpad1->SetTitleSize(0.03);
391     fHistpad1->GetXaxis()->SetLabelSize(0.03);
392     fHistpad1->GetXaxis()->SetTitleSize(0.03);
393     fHistpad1->GetYaxis()->SetLabelSize(0.03);
394     fHistpad1->GetYaxis()->SetTitleSize(0.03);
395
396     fHistpad2->SetTitleSize(0.03);
397     fHistpad2->GetXaxis()->SetLabelSize(0.03);
398     fHistpad2->GetXaxis()->SetTitleSize(0.03);
399     fHistpad2->GetYaxis()->SetLabelSize(0.03);
400     fHistpad2->GetYaxis()->SetTitleSize(0.03);
401
402     fHistpad3->SetTitleSize(0.03);
403     fHistpad3->GetXaxis()->SetLabelSize(0.03);
404     fHistpad3->GetXaxis()->SetTitleSize(0.03);
405     fHistpad3->GetYaxis()->SetLabelSize(0.03);
406     fHistpad3->GetYaxis()->SetTitleSize(0.03);
407
408     gStyle->SetPalette(1);
409     
410     SetHistPadRowAxis();
411 }
412
413 // ####################################################################################################
414 void AliHLTTPCDisplay::FillPadRow(Int_t patch, ULong_t dataBlock, ULong_t dataLen){
415     AliHLTTPCDigitReader* digitReader = new AliHLTTPCDigitReaderPacked();
416     bool readValue = true;
417     Int_t rowOffset = 0;
418
419     // Initialize RAW DATA
420     Int_t firstRow = AliHLTTPCTransform::GetFirstRow(patch);
421     Int_t lastRow = AliHLTTPCTransform::GetLastRow(patch);
422
423     // Outer sector, patches 2, 3, 4, 5 -  start counting in patch 2 with row 0
424     if ( patch >= 2 ) rowOffset = AliHLTTPCTransform::GetFirstRow( 2 );
425
426     // Initialize block for reading packed data
427     void* tmpdataBlock = (void*) dataBlock;
428     digitReader->InitBlock(tmpdataBlock,dataLen,firstRow,lastRow,patch,0);
429
430     readValue = digitReader->Next();
431
432     if (!readValue){    
433         LOG(AliHLTTPCLog::kError,"AliHLTTPCDisplay::FillPadRow","Read first value") << "No value in data block" << ENDLOG;
434         return;
435     }
436     
437     // FILL PADROW 3D --- Initialize the colorbins
438     if (fSwitch3DPadRow){
439         for (UInt_t ii=0;ii < 20;ii++){
440             fbinct[ii] = 0;
441             fcolorbin[ii] = 0;
442         }
443
444         // read number of entries in colorbin
445         while ( readValue ){ 
446
447             Int_t row = digitReader->GetRow() + rowOffset;
448             
449             if (row == fPadRow){    
450                 UInt_t charge = digitReader->GetSignal();
451                 
452                 for (UInt_t ii=0;ii < 19;ii++){
453                     if ( charge > (ii*15) && charge <= ((ii*15) + 15) ) fcolorbin[ii]++;
454                 }
455                 // larger than 19 * 15  
456                 if (charge > 285 ) fcolorbin[19]++;
457             }
458
459             // read next value
460             readValue = digitReader->Next();
461       
462             if(!readValue) break; //No more value
463         } 
464         //Initialize fpmarr[color][3*colorbin[ii]]  
465         fpmarr[0] = new Float_t[fcolorbin[0]*3]; 
466         fpmarr[1] = new Float_t[fcolorbin[1]*3]; 
467         fpmarr[2] = new Float_t[fcolorbin[2]*3]; 
468         fpmarr[3] = new Float_t[fcolorbin[3]*3]; 
469         fpmarr[4] = new Float_t[fcolorbin[4]*3];  
470         fpmarr[5] = new Float_t[fcolorbin[5]*3]; 
471         fpmarr[6] = new Float_t[fcolorbin[6]*3]; 
472         fpmarr[7] = new Float_t[fcolorbin[7]*3]; 
473         fpmarr[8] = new Float_t[fcolorbin[8]*3]; 
474         fpmarr[9] = new Float_t[fcolorbin[9]*3]; 
475         fpmarr[10] = new Float_t[fcolorbin[10]*3]; 
476         fpmarr[11] = new Float_t[fcolorbin[11]*3]; 
477         fpmarr[12] = new Float_t[fcolorbin[12]*3]; 
478         fpmarr[13] = new Float_t[fcolorbin[13]*3]; 
479         fpmarr[14] = new Float_t[fcolorbin[14]*3]; 
480         fpmarr[15] = new Float_t[fcolorbin[15]*3]; 
481         fpmarr[16] = new Float_t[fcolorbin[16]*3]; 
482         fpmarr[17] = new Float_t[fcolorbin[17]*3]; 
483         fpmarr[18] = new Float_t[fcolorbin[18]*3]; 
484         fpmarr[19] = new Float_t[fcolorbin[19]*3]; 
485         
486         // Rewind the raw reader and fill the polymarker3D
487         digitReader->InitBlock(tmpdataBlock,dataLen,firstRow,lastRow,patch,0);
488         
489         readValue = digitReader->Next();
490     } // END if (fSwitch3DPadRow)
491
492     // -- Fill Raw Data
493     while ( readValue ){ 
494
495         Int_t row = digitReader->GetRow() + rowOffset;
496
497         // select padrow to fill in histogramm
498         if (row == fPadRow){    
499             UChar_t pad = digitReader->GetPad();
500             UShort_t time = digitReader->GetTime();
501             UInt_t charge = digitReader->GetSignal();
502             Float_t xyz[3];
503             fHistraw->Fill(pad,time,charge);
504
505             if (pad == (fPad-1) ) fHistpad1->Fill(time,charge);
506             if (pad == fPad) fHistpad2->Fill(time,charge);
507             if (pad == (fPad+1) ) fHistpad3->Fill(time,charge);
508
509             if (fSwitch3DPadRow) {
510                 // Transform raw coordinates to local coordinates
511                 AliHLTTPCTransform::RawHLT2Global(xyz, fSlicePadRow, fPadRow, pad, time);
512
513                 for (UInt_t ii=0;ii < 19;ii++){
514                     if ( charge > (ii*15) && charge <= ((ii*15) + 15) ){
515                         fpmarr[ii][fbinct[ii]] = xyz[0];
516                         fpmarr[ii][fbinct[ii]+1] = xyz[1];
517                         fpmarr[ii][fbinct[ii]+2] = xyz[2];
518                         fbinct[ii] += 3;
519                     }
520                 }
521                 // larger than 19 * 15
522                 if (charge > 285 ) {
523                     fpmarr[19][fbinct[19]] = xyz[0];
524                     fpmarr[19][fbinct[19]+1] = xyz[1];
525                     fpmarr[19][fbinct[19]+2] = xyz[2];
526                     fbinct[19] += 3;
527                 }
528             } // END if (fSwitch3DPadRow)
529         
530         }
531         
532         // read next value
533         readValue = digitReader->Next();
534       
535         //Check where to stop:
536         if(!readValue) break; //No more value
537     } 
538     
539     if ( digitReader )
540         delete digitReader;
541     digitReader = NULL;
542
543     AliHLTTPCSpacePointData *points = fClusters[fSlicePadRow][patch];
544     if(!points) return;
545     Int_t npoints = fNcl[fSlicePadRow][patch];
546     
547     Float_t xyz[3];
548     for(Int_t i=0; i<npoints; i++){
549         xyz[0] = points[i].fX;
550         xyz[1] = points[i].fY;
551         xyz[2] = points[i].fZ;
552         
553         Int_t clrow = AliHLTTPCTransform::GetPadRow(xyz[0]);
554         // select padrow to fill in histogramm
555         if (clrow == fPadRow){
556             AliHLTTPCTransform::LocHLT2Raw(xyz, fSlicePadRow, fPadRow);
557             fHistrawcl->Fill(xyz[1],xyz[2]);
558         }
559     }
560 }
561
562 // #############################################################################
563 void AliHLTTPCDisplay::ResetHistPadRow(){  
564     fHistraw->Reset();   
565     fHistrawcl->Reset(); 
566     fHistpad1->Reset(); 
567     fHistpad2->Reset();  
568     fHistpad3->Reset(); 
569 }
570
571 // #############################################################################
572 void AliHLTTPCDisplay::ResetHistResiduals(){  
573     fHistallresidualsY->Reset();
574     fHistallresidualsZ->Reset();
575 }
576
577 // #############################################################################
578 void AliHLTTPCDisplay::ResetHistCharge(){  
579     fHistcharge->Reset();
580 }
581
582
583 // #############################################################################
584 //                 DRAWER
585 // #############################################################################
586 void AliHLTTPCDisplay::DrawGeomSector(Int_t sector) {  
587   Char_t fname[256];
588   Int_t realsector = sector;// % 18;
589   
590   if (realsector < 10){
591     sprintf(fname,"LS0%d",realsector);
592     fGeom->GetNode(fname)->SetLineColor(fLineColor);
593     fGeom->GetNode(fname)->Draw("same");
594     sprintf(fname,"US0%d",realsector);
595     fGeom->GetNode(fname)->SetLineColor(fLineColor); 
596     fGeom->GetNode(fname)->Draw("same");
597   }
598   else {
599     sprintf(fname,"LS%d",realsector);
600     fGeom->GetNode(fname)->SetLineColor(fLineColor);
601     fGeom->GetNode(fname)->Draw("same");
602     sprintf(fname,"US%d",realsector);
603     fGeom->GetNode(fname)->SetLineColor(fLineColor); 
604     fGeom->GetNode(fname)->Draw("same");
605   }   
606 }
607 // #############################################################################
608 void AliHLTTPCDisplay::DrawHistPadRow(){  
609     Char_t title[256];
610     sprintf(title,"Selected PadRow %d with found Clusters",fPadRow);
611
612     fHistraw->SetTitle(title);
613     fHistraw->SetStats(kFALSE);
614     fHistraw->Draw("COLZ");
615
616     fHistrawcl->SetStats(kFALSE);
617     fHistrawcl->SetMarkerStyle(28);
618     fHistrawcl->SetMarkerSize(2);
619     fHistrawcl->SetMarkerColor(1);
620     fHistrawcl->Draw("psame");
621 }
622
623 // #############################################################################
624 void AliHLTTPCDisplay::DrawHistPad1(){  
625     Char_t title[256];
626     sprintf(title,"Selected Pad %d",fPad -1);
627     fHistpad1->SetStats(kFALSE);
628     fHistpad1->SetTitle(title);
629     fHistpad1->Draw();
630 }
631
632 // #############################################################################
633 void AliHLTTPCDisplay::DrawHistPad2(){  
634     Char_t title[256];
635     sprintf(title,"Selected Pad %d",fPad);
636
637     fHistpad2->SetStats(kFALSE);
638     fHistpad2->SetTitle(title);
639     fHistpad2->Draw();
640 }
641
642 // #############################################################################
643 void AliHLTTPCDisplay::DrawHistPad3(){  
644     Char_t title[256];
645     sprintf(title,"Selected Pad %d",fPad +1);
646
647     fHistpad3->SetStats(kFALSE);
648     fHistpad3->SetTitle(title);
649     fHistpad3->Draw();
650 }
651
652 // #############################################################################
653 void AliHLTTPCDisplay::DrawHistResiduals(Bool_t ySwitch){  
654     if (fSwitch3DTracks){
655         if (ySwitch){
656             // Y Residual histogram for 1 track
657
658             if (fSelectTrackSwitch){
659                 Char_t title[256];
660                 sprintf(title,"Y Residuals of Track %d in Slice %d",fSelectTrack, fSelectTrackSlice );
661
662                 TMultiGraph *mgY = new TMultiGraph();
663
664
665 //              fGraphresidualsY->SetTitle(title);      
666                 fGraphresidualsY->GetXaxis()->SetTitle("padrow");       
667                 fGraphresidualsY->GetYaxis()->SetTitle("residuals");
668 //              fGraphresidualsY->Draw("A*");
669                 fGraphresidualsY->GetXaxis()->SetLabelSize(0.02);
670                 fGraphresidualsY->GetXaxis()->SetTitleSize(0.02);
671                 fGraphresidualsY->GetYaxis()->SetLabelSize(0.02);
672                 fGraphresidualsY->GetYaxis()->SetTitleSize(0.02);
673                 fGraphresidualsYLength->SetMarkerColor(2);
674                 fGraphresidualsYLength->SetMarkerStyle(5);
675                 fGraphresidualsY->SetMarkerColor(1);
676                 fGraphresidualsY->SetMarkerStyle(3);
677
678 //              fGraphresidualsY->Draw("A*");
679 //              fGraphresidualsYLength->Draw("*");
680
681                 mgY->Add(fGraphresidualsY);
682                 mgY->Add(fGraphresidualsYLength);
683                 mgY->SetTitle(title);
684 //              mgY->GetXaxis()->SetTitle("padrow");    
685 //              mgY->GetYaxis()->SetTitle("residuals");
686                 mgY->Draw("AP");
687             }
688             // Global residuals histogram
689             else{
690                 fHistallresidualsY->SetStats(kFALSE);
691                 fHistallresidualsY->Draw();
692             }
693         }
694         else {
695             // Z Residual histogram for 1 track
696             if (fSelectTrackSwitch){
697                 Char_t title[256];
698                 sprintf(title,"Z Residuals of Track %d in Slice %d",fSelectTrack, fSelectTrackSlice );
699
700                 TMultiGraph *mgZ = new TMultiGraph();
701
702                 fGraphresidualsZ->SetTitle(title);      
703                 fGraphresidualsZ->GetXaxis()->SetTitle("padrow");       
704                 fGraphresidualsZ->GetYaxis()->SetTitle("residuals");
705                 fGraphresidualsZ->GetXaxis()->SetLabelSize(0.02);
706                 fGraphresidualsZ->GetXaxis()->SetTitleSize(0.02);
707                 fGraphresidualsZ->GetYaxis()->SetLabelSize(0.02);
708                 fGraphresidualsZ->GetYaxis()->SetTitleSize(0.02);
709 //              fGraphresidualsZLength->Draw("F*");
710 //              fGraphresidualsZ->Draw("A*");
711         
712                 mgZ->Add(fGraphresidualsZ);
713 //              mgZ->Add(fGraphresidualsZLength);
714                 mgZ->SetTitle(title);
715                 mgZ->Draw("A*");
716             }
717             // Global residuals histogram
718             else{
719                 fHistallresidualsZ->SetStats(kFALSE);
720                 fHistallresidualsZ->Draw();
721             }
722         }
723     }
724 }
725
726 // #############################################################################
727 void AliHLTTPCDisplay::DrawHistCharge(){  
728     if (fSwitch3DCluster){
729 //      fHistcharge->SetStats(kFALSE);
730         fHistcharge->Draw();
731     }
732 }
733
734 // #############################################################################
735 void AliHLTTPCDisplay::Draw3D(){        
736     
737 #if defined(HAVE_TVIEW3D_H)
738     TView3D *v = new TView3D();
739     if (v) v->SetSystem(1);
740 #else
741     TView *v = new TView(1);
742 #endif
743     if (v==NULL) {
744       HLTFatal("can not create viewer");
745       return;
746     }
747     v->SetRange(-800,-800,-800,800,800,800);
748
749     Float_t* etaRange = NULL;   // ------  STILL TO FIX
750     
751
752     //--------------------------------------------------------------------------------------------
753     // DRAW 3D CLUSTER
754     //--------------------------------------------------------------------------------------------
755     if (fSwitch3DCluster){
756         Int_t maxCharge = 0;
757
758         for (Int_t slice=0; slice <= 35; slice++){
759
760             Int_t currenttrack = -1;
761
762             if (fSelectCluster == 1 && fSelectTrackSwitch && slice == fSelectTrackSlice ){
763
764                 Int_t trackcounter = 0;
765                 Int_t ntracks = fTracks->GetNTracks();
766         
767                 for(Int_t j=0; j<ntracks; j++) {        
768
769                     AliHLTTPCTrack *gtrack = fTracks->GetCheckedTrack(j); 
770                     if(!gtrack) continue;
771
772                     Int_t nHits = gtrack->GetNHits();  // Number of associated hits to track
773                     Int_t tmpslice = gtrack->GetSector();
774
775                     // --- CHECK if track is should be drawn
776                     // select Single Track
777                     if(tmpslice != fSelectTrackSlice) continue;
778                         
779                     if (trackcounter != fSelectTrack){
780                         trackcounter++;  
781                         continue;
782                     }
783                     trackcounter++;
784
785                     if((fPtThreshold > 0) && (gtrack->GetPt()< fPtThreshold)) continue;
786                     if(nHits < fMinHits) continue;
787
788                     currenttrack = j;
789                     break;
790                 }
791             }
792             
793             if (!fSliceArray[slice]) continue;
794             
795             for(Int_t p=0;p<6;p++){
796
797                 AliHLTTPCSpacePointData *points = fClusters[slice][p];
798                 if(!points) continue;
799                 Int_t npoints = fNcl[slice][p];
800                 TPolyMarker3D *pmUsed = new TPolyMarker3D(1,6);
801                 TPolyMarker3D *pmUnused = new TPolyMarker3D(1,6);
802                 Int_t nUsedCluster = 0;
803                 Int_t nUnusedCluster = 0;
804
805                 Float_t xyz[3];
806                 for(Int_t i=0; i<npoints; i++){
807                     // Used  cluster only
808                     if (fSelectCluster == 1  && points[i].fUsed == kFALSE) continue; 
809                     // Unused cluster only
810                     if (fSelectCluster == 2  && points[i].fUsed == kTRUE) continue; 
811
812                     // if single track is selcted draw only cluster for this track
813                     if (fSelectCluster == 1 && fSelectTrackSwitch && points[i].fTrackN != currenttrack) continue;
814                     
815                     xyz[0] = points[i].fX;
816                     xyz[1] = points[i].fY;
817                     xyz[2] = points[i].fZ;
818                     
819                     if ( etaRange ){              
820                         // Do this before the transform, because the tracker also uses
821                         // local coordinates when using this limit to determine 
822                         // which clusters to use for tracking
823                         Double_t pointEta = AliHLTTPCTransform::GetEta( xyz );
824                         if ( pointEta<etaRange[0] || pointEta>etaRange[1] )
825                             continue;
826                     }
827
828                     AliHLTTPCTransform::Local2Global(xyz,slice);
829                  
830                     if (points[i].fUsed == kTRUE){
831                         pmUsed->SetPoint(nUsedCluster,xyz[0],xyz[1],xyz[2]);
832                         nUsedCluster++;
833                     }
834                     else {
835                         pmUnused->SetPoint(nUnusedCluster,xyz[0],xyz[1],xyz[2]);
836                         nUnusedCluster++;
837                     }
838
839                     // Fill Charge Histogram
840                     fHistcharge->Fill(points[i].fCharge);
841                     if ((Int_t)points[i].fCharge > maxCharge ) maxCharge = (Int_t) points[i].fCharge; 
842                 }
843                 pmUsed->SetMarkerSize(1);
844                 pmUsed->SetMarkerColor(3); 
845                 pmUsed->Draw("");
846
847                 pmUnused->SetMarkerSize(1);
848                 pmUnused->SetMarkerColor(2); 
849                 pmUnused->Draw("");
850             } // END - PATCH LOOP           
851         }  // END - SLICE LOOP
852         fHistcharge->SetAxisRange(0,maxCharge);
853     }   // END - DRAW 3D CLUSTER 
854
855     //--------------------------------------------------------------------------------------------
856     // DRAW 3D TRACKS
857     //--------------------------------------------------------------------------------------------
858     if (fSwitch3DTracks){
859
860         Int_t trackcounter = 0;
861         Int_t ntracks = fTracks->GetNTracks();
862         Double_t drawStep = 0.2;
863
864         Double_t maxResidualY = 0.;
865         Double_t maxResidualZ = 0.;
866
867         TPolyLine3D *line = new TPolyLine3D[ntracks];
868 #if TRACKHELIX
869         THelix *helix = new THelix[ntracks];
870 #endif
871         for(Int_t j=0; j<ntracks; j++) {        
872
873             AliHLTTPCTrack *gtrack = fTracks->GetCheckedTrack(j); 
874             if(!gtrack) continue;
875
876             Int_t nHits = gtrack->GetNHits();  // Number of associated hits to track
877             Int_t slice = gtrack->GetSector();
878
879             // --- CHECK if track is should be drawn
880             // select if slice should be displayed or not
881             if (!fSliceArray[slice]) continue;  
882             
883             // select Single Track
884             if (fSelectTrackSwitch){
885                 if(slice != fSelectTrackSlice) continue;
886         
887                 if (trackcounter != fSelectTrack){
888                     trackcounter++;  
889                     continue;
890                 }
891                 trackcounter++;
892             }
893     
894             if((fPtThreshold > 0) && (gtrack->GetPt()< fPtThreshold)) continue;
895             if(nHits < fMinHits) continue;
896             
897             TPolyMarker3D *pmL = new TPolyMarker3D(1,2);
898             TPolyMarker3D *pmF = new TPolyMarker3D(1,2);
899
900             Double_t radius = gtrack->GetRadius();      // radius
901             Double_t kappa = gtrack->GetKappa();        // curvature = 1/R , signed
902             Double_t lambda = atan( gtrack->GetTgl() ); // dipAngle lambda
903             Double_t phi0 = gtrack->GetPsi() + (gtrack->GetCharge() * AliHLTTPCTransform::PiHalf() ); // azimuthal angle of startingpoint, with respect to helix axis
904
905             Double_t xyzL[3];      // lastpoint of track
906             Double_t xyzF[3];      // firstpoint of track
907
908             xyzF[0] = gtrack->GetFirstPointX();
909             xyzF[1] = gtrack->GetFirstPointY();
910             xyzF[2] = gtrack->GetFirstPointZ();
911             pmF->SetPoint(0,xyzF[0],xyzF[1],xyzF[2]);
912
913             xyzL[0] = gtrack->GetLastPointX();
914             xyzL[1] = gtrack->GetLastPointY();
915             xyzL[2] = gtrack->GetLastPointZ();
916             pmL->SetPoint(0,xyzL[0],xyzL[1],xyzL[2]);
917
918             Double_t s = 0.;       // length of the track
919
920             // Calculate the length of the track. If it is to flat in in s,z plane use sxy, otherwise use sz
921             if (fabs(lambda) > 0.05){
922                 // length of track calculated out of z
923                 s = fabs( (xyzL[2] - xyzF[2]) / sin(lambda) ); // length of track calculated out of z
924             }
925             else {
926                 Double_t d =  (xyzL[0] - xyzF[0])*(xyzL[0] - xyzF[0]) + (xyzL[1] - xyzF[1])*(xyzL[1] - xyzF[1]);
927                 // length of track calculated out of xy
928                 s = fabs ( acos( 0.5 * (2 - (d / (radius*radius)))) / ( kappa * cos(lambda) ) );                
929             }
930             
931             Int_t nTrackPoints = 2 + (Int_t) floor(s / drawStep);
932
933 #if TRACKPOLYMARKER
934             TPolyMarker3D *pmT = new TPolyMarker3D(nTrackPoints,6);
935 #endif
936
937             Double_t *xT = new Double_t[nTrackPoints];
938             Double_t *yT = new Double_t[nTrackPoints];
939             Double_t *zT = new Double_t[nTrackPoints];
940
941             //Write Track Parameters for single track
942             if (fSelectTrackSwitch){
943                 fTrackParam.id = trackcounter - 1;
944                 fTrackParam.nHits = nHits;
945                 fTrackParam.charge = gtrack->GetCharge();
946                 fTrackParam.lambda = lambda;
947                 fTrackParam.kappa = kappa;
948                 fTrackParam.radius = radius;
949                 fTrackParam.slice = slice;
950                 fTrackParam.phi0 = phi0;
951                 fTrackParam.pt = gtrack->GetPt();
952                 fTrackParam.bfield = AliHLTTPCTransform::GetBFieldValue();
953                 fTrackParam.xyzF[0] = gtrack->GetFirstPointX();
954                 fTrackParam.xyzF[1] = gtrack->GetFirstPointY();
955                 fTrackParam.xyzF[2] = gtrack->GetFirstPointZ();
956                 fTrackParam.xyzL[0] = gtrack->GetLastPointX();
957                 fTrackParam.xyzL[1] = gtrack->GetLastPointY();
958                 fTrackParam.xyzL[2] = gtrack->GetLastPointZ();
959                 fTrackParam.psi = gtrack->GetPsi();
960                 fTrackParam.s = s;
961             }
962
963             Int_t trackPointCounter = 0;
964
965             for (Double_t ds = 0.; ds < s; ds = ds + drawStep){
966                 // FILL ARRAYS IN ORDER TO DRAW THE TRACKPOINTS, OUT OF THE PARAMETER
967                 xT[trackPointCounter] = xyzF[0] + radius * ( cos( phi0 + (ds*kappa*cos(lambda)) ) - cos(phi0) );
968                 yT[trackPointCounter] = xyzF[1] + radius * ( sin( phi0 + (ds*kappa*cos(lambda)) ) - sin(phi0) );
969                 zT[trackPointCounter] = xyzF[2] + ds * sin(lambda);
970 #if TRACKPOLYMARKER
971                 pmT->SetPoint(trackPointCounter,xT[trackPointCounter],yT[trackPointCounter],zT[trackPointCounter]);
972 #endif
973                 trackPointCounter++;
974             }
975
976             xT[trackPointCounter] = xyzF[0] + radius * ( cos( phi0 + (s*kappa*cos(lambda)) ) - cos(phi0) );
977             yT[trackPointCounter] = xyzF[1] + radius * ( sin( phi0 + (s*kappa*cos(lambda)) ) - sin(phi0) );
978             zT[trackPointCounter] = xyzF[2] + s * sin(lambda);
979 #if TRACKPOLYMARKER       
980             pmT->SetPoint(trackPointCounter,xT[trackPointCounter],yT[trackPointCounter],zT[trackPointCounter]);
981 #endif
982             // --- RESIDUALS ---
983             gtrack->Rotate(slice,kTRUE);
984             Int_t nRes = 0;  // number of resiudals
985             
986             UInt_t *hitnum = gtrack->GetHitNumbers();
987
988             Double_t *resY= new Double_t[nHits];
989             Double_t *resZ= new Double_t[nHits];
990
991             Double_t *resYLength= new Double_t[2*nHits];
992             Double_t *resZLength= new Double_t[2*nHits];
993
994             Double_t *padrows = new Double_t[nHits];
995             Double_t *padrowsLength = new Double_t[2*nHits];
996
997             for(Int_t h=0; h<nHits; h++){
998                 UInt_t id=hitnum[h];
999                 Int_t patch = (id>>22) & 0x7;
1000                 UInt_t pos = id&0x3fffff; 
1001
1002                 AliHLTTPCSpacePointData *points = fClusters[slice][patch];
1003
1004                 Float_t xyzCtmp[3];    // cluster tmp
1005                 Float_t xyzTtmp[3];    // track tmp
1006
1007                 xyzCtmp[0] = points[pos].fX;
1008                 xyzCtmp[1] = points[pos].fY;
1009                 xyzCtmp[2] = points[pos].fZ;
1010
1011                 Int_t padrow = AliHLTTPCTransform::GetPadRow(points[pos].fX);
1012                 xyzTtmp[0] = gtrack->GetFirstPointX();
1013                 if(gtrack->GetCrossingPoint(padrow,xyzTtmp)) {
1014
1015                     Float_t deltaY = ( xyzCtmp[1] - xyzTtmp[1] );
1016                     Float_t deltaZ = ( xyzCtmp[2] - xyzTtmp[2] );
1017 //                  Float_t residual = sqrt( deltaY*deltaY + deltaZ*deltaZ );
1018                     
1019                     padrows[nRes] = (Double_t) padrow;
1020                     resY[nRes] = (Double_t) deltaY;
1021                     resZ[nRes] = (Double_t) deltaZ;
1022
1023                     resYLength[(2*nRes)] = 0.5 * AliHLTTPCTransform::GetPadLength(padrow);
1024                     resYLength[(2*nRes)+1] = -0.5 * AliHLTTPCTransform::GetPadLength(padrow);
1025                     resZLength[nRes] = AliHLTTPCTransform::GetZLength();
1026                     padrowsLength[(2*nRes)] = (Double_t) padrow;
1027                     padrowsLength[(2*nRes)+1] = (Double_t) padrow;
1028
1029                     // FILL RESIDUALS HISTOGRAM
1030                     fHistallresidualsY->Fill(resY[nRes]);
1031                     fHistallresidualsZ->Fill(resZ[nRes]);
1032                     if (resY[nRes] > maxResidualY ) maxResidualY = resY[nRes];
1033                     if (resZ[nRes] > maxResidualZ ) maxResidualZ = resZ[nRes];
1034                     nRes++;
1035                 }
1036             }
1037
1038             gtrack->Rotate(slice,kFALSE);
1039             // --- RESIDUALS ---
1040
1041             // Draw last point of Track
1042             pmL->SetMarkerSize(3);
1043             pmL->SetMarkerColor(4); 
1044 //          pmL->Draw();
1045
1046             // Draw first point of Track
1047             pmF->SetMarkerSize(3);
1048             pmF->SetMarkerColor(5); 
1049 //          pmF->Draw();
1050
1051 #if TRACKPOLYMARKER
1052             // Draw Track -- as polymarker
1053             pmT->SetMarkerSize(3);
1054             pmT->SetMarkerColor(3); 
1055             pmT->Draw();
1056 #endif
1057             // Draw Track -- as line
1058             TPolyLine3D *currentline = &(line[j]);
1059             currentline = new TPolyLine3D(nTrackPoints,xT,yT,zT,"");
1060             currentline->SetLineColor(4);   
1061             currentline->SetLineWidth(2);
1062             currentline->Draw("same");
1063
1064 #if TRACKHELIX
1065             // Draw Track -- as helix
1066             // works ok, execpt for very small dipangles -> track almost horizontal
1067             Double_t hrange[2];
1068             Double_t v0[3];
1069             Double_t omega;
1070             hrange[0] = xyzF[2];
1071             hrange[1] = xyzL[2];
1072             v0[0] = gtrack->GetPx();
1073             v0[1] = gtrack->GetPy();
1074             v0[2] = gtrack->GetPz();
1075             omega = AliHLTTPCTransform::GetBFieldValue() * gtrack->GetCharge();
1076
1077             THelix *currenthelix = &(helix[j]);
1078             currenthelix = new THelix(xyzF,v0,omega,hrange,kHelixZ,0);
1079             currenthelix->SetLineColor(6);   
1080             currenthelix->SetLineWidth(1);
1081             currenthelix->Draw("same");             
1082 #endif
1083
1084             //Residuals
1085             if ( fGraphresidualsY){
1086                 delete fGraphresidualsY;
1087                 fGraphresidualsY = NULL;
1088             }
1089
1090             if ( fGraphresidualsZ){
1091                 delete fGraphresidualsZ;
1092                 fGraphresidualsZ = NULL;
1093             }
1094             //Residuals
1095             if ( fGraphresidualsYLength){
1096                 delete fGraphresidualsYLength;
1097                 fGraphresidualsYLength = NULL;
1098             }
1099
1100             if ( fGraphresidualsZLength){
1101                 delete fGraphresidualsZLength;
1102                 fGraphresidualsZLength = NULL;
1103             }
1104
1105
1106
1107             // FILL Y RESIDUALS GRAPH
1108             fGraphresidualsY = new TGraph(nRes-1,padrows,resY);
1109             fGraphresidualsYLength = new TGraph((2*nRes)-2,padrowsLength,resYLength);
1110             // FILL Z RESIDUALS GRAPH
1111             fGraphresidualsZ = new TGraph(nRes-1,padrows,resZ);
1112             fGraphresidualsZLength = new TGraph(nRes-1,padrows,resZLength);
1113
1114             if (xT) delete xT;
1115             if (yT) delete yT;
1116             if (zT) delete zT;
1117
1118         } // END for tracks
1119
1120         fHistallresidualsY->SetAxisRange(-maxResidualY,maxResidualY);
1121         fHistallresidualsZ->SetAxisRange(-maxResidualZ,maxResidualZ);
1122
1123     }   // END - DRAW 3D Tracks
1124
1125     //--------------------------------------------------------------------------------------------
1126     // DRAW 3D GEOMETRY
1127     //--------------------------------------------------------------------------------------------
1128     if (fSwitch3DGeometry){
1129
1130         for (Int_t slice=0; slice <= 17; slice++){
1131             if (!fSliceArray[slice]) continue;
1132             DrawGeomSector(slice);
1133         }
1134     }   // END - DRAW 3D GEOMETRY
1135     
1136     //--------------------------------------------------------------------------------------------
1137     // DRAW 3D PadRow
1138     //--------------------------------------------------------------------------------------------
1139     if (fSwitch3DPadRow && fSliceArray[fSlicePadRow]){
1140         Int_t markercolor = 51;
1141
1142         for (UInt_t ii=0;ii < 20;ii++){
1143             if (fcolorbin[ii]> 0){
1144                 
1145                 TPolyMarker3D *pm = new TPolyMarker3D(fcolorbin[ii], fpmarr[ii], 7 );
1146
1147                 pm->SetMarkerColor(markercolor); 
1148                 pm->Draw(""); 
1149             }
1150
1151             // in order to have the SetPalette(1), so called "pretty"
1152             if (ii % 2 == 0 ) markercolor += 2;
1153             else  markercolor += 3;
1154         }
1155     }
1156
1157     //--------------------------------------------------------------------------------------------
1158     // DRAW 3D 
1159     //--------------------------------------------------------------------------------------------
1160     v->ZoomView(0,4);
1161     v->Draw();   
1162 }
1163
1164 // ---------------------------------------------------
1165 // In order to be backward compatible
1166 // ---------------------------------------------------
1167 #if BACKWARD
1168 void AliHLTTPCDisplay::DisplayClusters(Bool_t x3don,Float_t* etaRange) {
1169     if (!fc1){
1170         fc1 = new TCanvas("c1","",900,900);
1171         fc1->cd();
1172     }
1173
1174     fSwitch3DTracks = kFALSE; 
1175     fSwitch3DCluster = kTRUE; 
1176     fSwitch3DPadRow = kFALSE; 
1177     fSwitch3DGeometry = kFALSE;
1178
1179     Draw3D();
1180 }
1181 // ---------------------------------------------------
1182 void AliHLTTPCDisplay::DisplayTracks(Int_t minhits,Bool_t x3don,Float_t thr) {
1183     if (!fc1){
1184         fc1 = new TCanvas("c1","",900,900);
1185         fc1->cd();
1186     }
1187
1188     fMinHits = minhits; 
1189     fPtThreshold = thr;
1190     fSwitch3DTracks = kTRUE; 
1191     fSwitch3DCluster = kFALSE; 
1192     fSwitch3DPadRow = kFALSE; 
1193     fSwitch3DGeometry = kFALSE;
1194
1195     Draw3D();
1196 }
1197 // ---------------------------------------------------
1198 void AliHLTTPCDisplay::DisplayAll(Int_t minhits,Bool_t clusterswitch,Bool_t trackswitch,Bool_t x3don, Float_t thr, Float_t* etaRange){
1199     if (!fc1){
1200         fc1 = new TCanvas("c1","",900,900);
1201         fc1->cd();
1202     }
1203
1204     fMinHits = minhits; 
1205     fPtThreshold = thr;
1206     fSwitch3DTracks = trackswitch; 
1207     fSwitch3DCluster = clusterswitch; 
1208     fSwitch3DPadRow = kFALSE; 
1209     fSwitch3DGeometry = kFALSE;
1210
1211     Draw3D();
1212 }
1213 #endif
1214 // ---------------------------------------------------