]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/OnlineDisplay/AliHLTTPCDisplay3D.cxx
- package adapted to new HOMER version
[u/mrichter/AliRoot.git] / HLT / TPCLib / OnlineDisplay / AliHLTTPCDisplay3D.cxx
1 // $Id$
2
3 /** \class AliHLTTPCDisplayPadRow
4 <pre>
5 //_____________________________________________________________
6 // AliHLTTPCDisplay3D
7 //
8 // Display class for the HLT TPC-3D events.
9 </pre>
10 */
11 // Author: Jochen Thaeder <mailto:thaeder@kip.uni-heidelberg.de>
12 //*-- Copyright &copy ALICE HLT Group 
13
14 #define TRACKHELIX 0       // use THelix for tracks
15 #define TRACKPOLYMARKER 0  // use TPolymarker3D for tracks
16 #define FIRSTLASTPOINT 0   // show first / last point of tracks
17 #define DEBUG 1
18
19 #define DRAWSTEP 0.2
20
21 #define UNUSEDCLUSTERCOLOR 2
22 #define USEDCLUSTERCOLOR 3
23 #define TRACKCOLOR 4
24 #define TRACKPOLYMARKERCOLOR 5
25 #define TRACKHELIXCOLOR 6
26
27 #if defined(HAVE_HOMERREADER) 
28 #include "HOMERReader.h"
29 #endif // defined(HAVE_HOMERREADER) 
30
31 #include "AliHLTTPCDisplay3D.h"
32 #include "AliHLTTPCDisplayPadRow.h"
33
34 #include "AliHLTStdIncludes.h"
35 #include <TView.h>
36 #include <TPolyMarker3D.h>
37 #include <TPolyLine3D.h>
38 #include <TH2.h>
39 #include <TTree.h>
40 #include <TNode.h>
41 #include <TGeometry.h>
42 #include <TShape.h>
43 #include <TParticle.h>
44 #include <TFile.h>
45 #include <THelix.h>
46 #include <TStyle.h>
47 #include <TGraph.h>
48 #include <TMultiGraph.h>
49 #include <TAttText.h>
50 #include <TAxis.h>
51 #include <TCanvas.h>
52
53 #ifdef use_aliroot
54 #include <TClonesArray.h>
55 #include <AliRun.h>
56 #include <AliSimDigits.h>
57 #include <AliTPCParam.h>
58 #endif
59
60 #include "AliHLTTPCDefinitions.h"
61 #include "AliHLTDataTypes.h"
62 #include "AliHLTTPCSpacePointData.h"
63 #include "AliHLTTPCClusterDataFormat.h"
64 #include "AliHLTTPCTrackletDataFormat.h"
65
66
67 #include "AliHLTTPCDigitReader.h"
68 #include "AliHLTTPCDigitReaderRaw.h"
69 #include "AliHLT_C_Component_WrapperInterface.h"
70
71 #include "AliHLTTPCDisplayMain.h"
72
73 #include "AliHLTTPCLogging.h"
74 #include "AliHLTTPCDisplay.h"
75 #include "AliHLTTPCTransform.h"
76 #include "AliHLTTPCTrack.h"
77 #include "AliHLTTPCTrackArray.h"
78 #include "AliHLTTPCMemHandler.h"
79 #include "AliHLTTPCDigitReaderPacked.h"
80
81
82 #if __GNUC__ >= 3
83 using namespace std;
84 #endif
85
86 ClassImp(AliHLTTPCDisplay3D)
87
88 //____________________________________________________________________________________________________
89 AliHLTTPCDisplay3D::AliHLTTPCDisplay3D(AliHLTTPCDisplayMain* display, Char_t* gfile ) {
90     // constructor
91     fDisplay = display;
92     
93     fGeom = NULL;
94     LoadGeometrie(gfile);
95 }
96
97 //____________________________________________________________________________________________________
98 AliHLTTPCDisplay3D::~AliHLTTPCDisplay3D() {
99     // destructor   
100
101 }
102
103 //____________________________________________________________________________________________________
104 void AliHLTTPCDisplay3D::Save(){
105   fDisplay->GetCanvas3D()->SaveAs("HLT-3D-View.eps");
106 }
107
108 //____________________________________________________________________________________________________
109 void AliHLTTPCDisplay3D::Draw(){
110     fDisplay->GetCanvas3D()->cd();
111     fDisplay->GetCanvas3D()->Clear();
112
113     TView *v = new TView(1);
114     //    TView v(1);
115     v->SetRange(-800,-800,-800,800,800,800);
116
117     Float_t* etaRange = NULL;   // ------  STILL TO FIX
118
119
120     //--------------------------------------------------------------------------------------------
121     // DRAW 3D GEOMETRY
122     //--------------------------------------------------------------------------------------------
123     if (fDisplay->Get3DSwitchGeometry()){
124
125         TList* masterNodeList = fGeom->GetListOfNodes();
126         TNode* masterNode=0;
127         TIter next(masterNodeList);
128         
129         while ((masterNode = static_cast<TNode*> (next()))) {  
130         
131             TList* nodeList = masterNode->GetListOfNodes();
132             TNode* node=0;
133             TIter next(nodeList);
134             
135             while ((node = static_cast<TNode*> (next()))) {  
136                 
137                 ULong_t tmpslice = atol(node->GetName() + 2);
138
139                 if (fDisplay->GetDisplaySlice(tmpslice)) {
140                     node->SetVisibility(1);
141                     node->SetFillColor(0);
142                     node->SetLineColor(fDisplay->GetLineColor());
143                 }
144                 else node->SetVisibility(0);
145             } // end while son Nodes
146         } // end while master Node
147
148         fGeom->Draw("");
149     }   // END - DRAW 3D GEOMETRY
150
151     
152     //--------------------------------------------------------------------------------------------
153     // DRAW 3D CLUSTER
154     //--------------------------------------------------------------------------------------------
155     if (fDisplay->Get3DSwitchCluster() && fDisplay->ExistsClusterData()){
156
157         for (Int_t slice=0; slice <= 35; slice++){
158
159             Int_t currenttrack = -1;
160
161             if (fDisplay->ExistsTrackData() && fDisplay->GetSelectTrackSwitch()) currenttrack = fDisplay->GetGlobalTrack(slice);
162
163             if (!fDisplay->GetDisplaySlice(slice)) continue;
164             
165             for(Int_t patch=0;patch<6;patch++){
166
167                 AliHLTTPCSpacePointData *points = fDisplay->GetSpacePointDataPointer(slice,patch);
168                 if(!points) continue;
169
170                 TPolyMarker3D *pmUsed = new TPolyMarker3D(1,6);
171                 TPolyMarker3D *pmUnused = new TPolyMarker3D(1,6);
172                 pmUnused->SetBit(kCanDelete);
173                 pmUsed->SetBit(kCanDelete);
174
175
176                 Int_t nUsedCluster = 0;
177                 Int_t nUnusedCluster = 0;
178
179                 Float_t xyz[3];    
180                 for(Int_t i=0; i< fDisplay->GetNumberSpacePoints(slice,patch); i++){
181                     // Used  cluster only
182                     if (fDisplay->GetSelectCluster() == 1  && points[i].fUsed == kFALSE) continue; 
183                     // Unused cluster only
184                     if (fDisplay->GetSelectCluster() == 2  && points[i].fUsed == kTRUE) continue; 
185                     
186                     // if single track is selcted draw only cluster for this track
187                     if (fDisplay->GetSelectCluster() == 1 && fDisplay->GetSelectTrackSwitch() && points[i].fTrackN != currenttrack) continue;
188                     
189                     xyz[0] = points[i].fX;
190                     xyz[1] = points[i].fY;
191                     xyz[2] = points[i].fZ;
192                     
193                     if ( etaRange ){              
194                         // Do this before the transform, because the tracker also uses
195                         // local coordinates when using this limit to determine 
196                         // which clusters to use for tracking
197                         Double_t pointEta = AliHLTTPCTransform::GetEta( xyz );
198                         if ( pointEta<etaRange[0] || pointEta>etaRange[1] )
199                             continue;
200                     }
201                     
202                     AliHLTTPCTransform::Local2Global(xyz,slice);
203                     
204                     if (points[i].fUsed == kTRUE){
205                         pmUsed->SetPoint(nUsedCluster,xyz[0],xyz[1],xyz[2]);
206                         nUsedCluster++;
207                     }
208                     else {
209                         pmUnused->SetPoint(nUnusedCluster,xyz[0],xyz[1],xyz[2]);
210                         nUnusedCluster++;
211                     }
212                     
213                 }
214                 pmUsed->SetMarkerSize(1);
215                 pmUsed->SetMarkerColor(USEDCLUSTERCOLOR); 
216                 pmUsed->Draw("same");
217
218                 pmUnused->SetMarkerSize(1);
219                 pmUnused->SetMarkerColor(UNUSEDCLUSTERCOLOR); 
220                 pmUnused->Draw("same");
221
222                 fDisplay->GetCanvas3D()->Modified();
223                 fDisplay->GetCanvas3D()->Update();
224
225             } // END - PATCH LOOP           
226         }  // END - SLICE LOOP  
227     }   // END - DRAW 3D CLUSTER 
228
229     //--------------------------------------------------------------------------------------------
230     // DRAW 3D TRACKS
231     //--------------------------------------------------------------------------------------------
232     if (fDisplay->Get3DSwitchTracks() && fDisplay->ExistsTrackData()){
233
234       AliHLTTPCTransform::SetBField( 0.5 );  // ++++++
235
236
237         AliHLTTPCTrackArray* tracks = fDisplay->GetTrackArrayPointer();
238         Int_t ntracks = tracks->GetNTracks();
239
240         //      TPolyLine3D **line = new (TPolyLine3D*)[ntracks];
241         //      for(Int_t j=0; j<ntracks; j++) line[j] = 0;
242 #if TRACKHELIX
243         //      THelix **helix = new (THelix*)[ntracks];
244         //      for(Int_t j=0; j<ntracks; j++) helix[j] = 0;
245
246 #endif
247         for(Int_t j=0; j<ntracks; j++) {        
248
249             AliHLTTPCTrack *gtrack = tracks->GetCheckedTrack(j); 
250             if(!gtrack) continue;
251
252             Int_t nHits = gtrack->GetNHits();  // Number of associated hits to track
253             Int_t slice = gtrack->GetSector();
254
255             // --- CHECK if track is should be drawn
256             // select if slice should be displayed or not
257             if (!fDisplay->GetDisplaySlice(slice)) continue;
258
259             if (fDisplay->GetSelectTrackSwitch() && fDisplay->GetGlobalTrack(slice) != j) continue;
260
261             Double_t radius = gtrack->GetRadius();      // radius
262             Double_t kappa = gtrack->GetKappa();        // curvature = 1/R , signed
263             Double_t lambda = atan( gtrack->GetTgl() ); // dipAngle lambda
264             Double_t phi0 = gtrack->GetPsi() + (gtrack->GetCharge() * AliHLTTPCTransform::PiHalf() ); // azimuthal angle of startingpoint, with respect to helix axis
265
266             if (kappa == 0 && AliHLTTPCTransform::GetBFieldValue() > 0.) {
267                 printf("================================KAPPA == 0");
268                 continue;
269             }
270
271             Double_t xyzL[3];      // lastpoint of track
272             Double_t xyzF[3];      // firstpoint of track
273
274             xyzF[0] = gtrack->GetFirstPointX();
275             xyzF[1] = gtrack->GetFirstPointY();
276             xyzF[2] = gtrack->GetFirstPointZ();
277
278             xyzL[0] = gtrack->GetLastPointX();
279             xyzL[1] = gtrack->GetLastPointY();
280             xyzL[2] = gtrack->GetLastPointZ();
281
282 #if FIRSTLASTPOINT          
283             //      TPolyMarker3D *pmL = new TPolyMarker3D(1,2);
284             //TPolyMarker3D *pmF = new TPolyMarker3D(1,2);
285
286             TPolyMarker3D pmL(1,2);
287             TPolyMarker3D pmF(1,2);
288
289
290             pmF.SetPoint(0,xyzF[0],xyzF[1],xyzF[2]);
291             pmL.SetPoint(0,xyzL[0],xyzL[1],xyzL[2]);
292 #endif
293
294             Double_t s = 0.;       // length of the track
295
296
297             if (  AliHLTTPCTransform::GetBFieldValue() == 0.) 
298                 s = sqrt ( (xyzL[0] - xyzF[0])*(xyzL[0] - xyzF[0]) + (xyzL[1] - xyzF[1])*(xyzL[1] - xyzF[1]) ); 
299             else {
300                 // Calculate the length of the track. If it is to flat in in s,z plane use sxy, otherwise use sz
301                 if (fabs(lambda) > 0.05){
302                     // length of track calculated out of z
303                     s = fabs( (xyzL[2] - xyzF[2]) / sin(lambda) ); // length of track calculated out of z
304                 }
305                 else {
306                     Double_t d = (xyzL[0] - xyzF[0])*(xyzL[0] - xyzF[0]) + (xyzL[1] - xyzF[1])*(xyzL[1] - xyzF[1]); 
307                     // length of track calculated out of xy
308                     s = fabs ( acos( 0.5 * (2 - (d / (radius*radius)))) / ( kappa * cos(lambda) ) );            
309                 }
310             }
311
312
313             // CUTS on tracks
314 /*
315             if (nHits < fDisplay->GetCutHits() ) continue;
316             if (s < fDisplay->GetCutS() ) continue;
317             if (gtrack->GetPsi() < fDisplay->GetCutPsi() ) continue;
318             if (lambda < fDisplay->GetCutLambda() ) continue; 
319             if (gtrack->GetPt() < fDisplay->GetCutPt()  &&  AliHLTTPCTransform::GetBFieldValue() != 0. ) continue;
320             if ( AliHLTTPCTransform::GetPadRow((Float_t)xyzF[0]) >   fDisplay->GetIncidentPadrow() ) continue;
321 */    
322             Int_t nTrackPoints = 2 + (Int_t) floor(s / DRAWSTEP);
323                             
324 #if TRACKPOLYMARKER
325             //    TPolyMarker3D *pmT = new TPolyMarker3D(nTrackPoints,6);
326             TPolyMarker3D pmT(nTrackPoints,6);
327
328 #endif
329             Double_t *xT = new Double_t[nTrackPoints];
330             Double_t *yT = new Double_t[nTrackPoints];
331             Double_t *zT = new Double_t[nTrackPoints];
332
333             Int_t trackPointCounter = 0;
334
335             //Write Track Parameters for single track
336             if (fDisplay->GetSelectTrackSwitch() ){
337                 fDisplay->fTrackParam.id = j;
338                 fDisplay->fTrackParam.nHits = nHits;
339                 fDisplay->fTrackParam.charge = gtrack->GetCharge();
340                 fDisplay->fTrackParam.lambda = lambda;
341                 fDisplay->fTrackParam.kappa = kappa;
342                 fDisplay->fTrackParam.radius = radius;
343                 fDisplay->fTrackParam.slice = slice;
344                 fDisplay->fTrackParam.phi0 = phi0;
345                 fDisplay->fTrackParam.pt = gtrack->GetPt();
346                 fDisplay->fTrackParam.bfield = AliHLTTPCTransform::GetBFieldValue();
347                 fDisplay->fTrackParam.psi = gtrack->GetPsi();
348                 fDisplay->fTrackParam.s = s;
349             }
350
351
352             if (  AliHLTTPCTransform::GetBFieldValue() == 0.) {
353
354                 for (Double_t ds = 0.; ds < s; ds = ds + DRAWSTEP){
355                     // FILL ARRAYS IN ORDER TO DRAW THE TRACKPOINTS, OUT OF THE PARAMETER
356                     xT[trackPointCounter] = xyzF[0] + ds * cos(phi0); 
357                     yT[trackPointCounter] = xyzF[1] + ds * sin(phi0);
358                     zT[trackPointCounter] = xyzF[2] + ds * sin(lambda);
359 #if TRACKPOLYMARKER
360                     pmT.SetPoint(trackPointCounter,xT[trackPointCounter],yT[trackPointCounter],zT[trackPointCounter]);
361 #endif
362                     trackPointCounter++;
363                 }
364                 
365                 if (trackPointCounter > nTrackPoints) printf("N=%d  n=%d", nTrackPoints,trackPointCounter);
366                 else {
367                     xT[trackPointCounter] = xyzF[0] + s * cos(phi0); 
368                     yT[trackPointCounter] = xyzF[1] + s * sin(phi0);
369                     zT[trackPointCounter] = xyzF[2] + s * sin(lambda);
370 #if TRACKPOLYMARKER       
371                     pmT.SetPoint(trackPointCounter,xT[trackPointCounter],yT[trackPointCounter],zT[trackPointCounter]);
372 #endif
373                 }
374
375             }
376             else {
377
378                 for (Double_t ds = 0.; ds < s; ds = ds + DRAWSTEP){
379                     // FILL ARRAYS IN ORDER TO DRAW THE TRACKPOINTS, OUT OF THE PARAMETER
380                     xT[trackPointCounter] = xyzF[0] + radius * ( cos( phi0 + (ds*kappa*cos(lambda)) ) - cos(phi0) );
381                     yT[trackPointCounter] = xyzF[1] + radius * ( sin( phi0 + (ds*kappa*cos(lambda)) ) - sin(phi0) );
382                     zT[trackPointCounter] = xyzF[2] + ds * sin(lambda);
383 #if TRACKPOLYMARKER
384                     pmT.SetPoint(trackPointCounter,xT[trackPointCounter],yT[trackPointCounter],zT[trackPointCounter]);
385 #endif
386                     trackPointCounter++;
387                 }
388                 
389                 if (trackPointCounter > nTrackPoints) printf("N=%d  n=%d", nTrackPoints,trackPointCounter);
390                 else {
391                     xT[trackPointCounter] = xyzF[0] + radius * ( cos( phi0 + (s*kappa*cos(lambda)) ) - cos(phi0) );
392                     yT[trackPointCounter] = xyzF[1] + radius * ( sin( phi0 + (s*kappa*cos(lambda)) ) - sin(phi0) );
393                     zT[trackPointCounter] = xyzF[2] + s * sin(lambda);
394 #if TRACKPOLYMARKER       
395                     pmT.SetPoint(trackPointCounter,xT[trackPointCounter],yT[trackPointCounter],zT[trackPointCounter]);
396 #endif
397                 }
398             }
399
400             // Draw Track -- as line
401             //line[j] = new TPolyLine3D(nTrackPoints,xT,yT,zT,"");
402             //TPolyLine3D *currentline = line[j];
403             //* currentline = new TPolyLine3D(nTrackPoints,xT,yT,zT,"");
404             //      TPolyLine3D currentline(nTrackPoints,xT,yT,zT,"");
405
406             TPolyLine3D *currentline = new TPolyLine3D(nTrackPoints,xT,yT,zT,"");
407             currentline->SetBit(kCanDelete);
408             currentline->SetLineColor(TRACKCOLOR);   
409             currentline->SetLineWidth(2);
410             currentline->Draw("same");
411
412
413
414             // ---  ADDITIONAL DRAW OPTIONS
415 #if FIRSTLASTPOINT
416             // Draw last point of Track
417             pmL.SetMarkerSize(3);
418             pmL.SetMarkerColor(4); 
419             pmL.Draw();
420
421             // Draw first point of Track
422             pmF.SetMarkerSize(3);
423             pmF.SetMarkerColor(5); 
424             pmF.Draw();
425 #endif
426 #if TRACKPOLYMARKER
427             // Draw Track -- as polymarker
428             pmT.SetMarkerSize(3);
429             pmT.SetMarkerColor(TRACKPOLYMARKERCOLOR); 
430             pmT.Draw();
431 #endif
432 #if TRACKHELIX
433             // Draw Track -- as helix
434             // works ok, execpt for very small dipangles -> track almost horizontal
435             Double_t hrange[2];
436             Double_t v0[3];
437             Double_t omega;
438             hrange[0] = xyzF[2];
439             hrange[1] = xyzL[2];
440             v0[0] = gtrack->GetPx();
441             v0[1] = gtrack->GetPy();
442             v0[2] = gtrack->GetPz();
443             omega = AliHLTTPCTransform::GetBFieldValue() * gtrack->GetCharge();
444
445             //      helix[j] = new THelix(xyzF,v0,omega,hrange,kHelixZ,0);
446             //      THelix *currenthelix = helix[j];
447             //      currenthelix = new THelix(xyzF,v0,omega,hrange,kHelixZ,0);
448             THelix currenthelix(xyzF,v0,omega,hrange,kHelixZ,0);
449             currenthelix.SetLineColor(TRACKHELIXCOLOR);   
450             currenthelix.SetLineWidth(1);
451             currenthelix.Draw("same");              
452 #endif
453
454
455             // delete[]
456             //Double_t *xT = new Double_t[nTrackPoints];
457             //      Double_t *yT = new Double_t[nTrackPoints];
458             //      Double_t *zT = new Double_t[nTrackPoints];
459             if (xT){ 
460               delete[] xT; 
461               xT = NULL;
462             }
463             if (yT){ 
464               delete[] yT; 
465               yT = NULL;
466             }
467             if (zT){ 
468               delete[] zT; 
469               zT = NULL;
470             }
471
472         } // END for track loop
473
474
475         // NO !!! DELETE line #ifdef helix delete helix
476
477
478     }   // END - DRAW 3D Tracks
479
480     //--------------------------------------------------------------------------------------------
481     // DRAW 3D PadRow
482     //--------------------------------------------------------------------------------------------
483     if (fDisplay->ExistsRawData() &&  fDisplay->Get3DSwitchPadRow() && fDisplay->GetDisplaySlice(fDisplay->GetSlicePadRow())){
484
485       // -- only one padrow
486       if ( fDisplay->Get3DSwitchRaw() == 0 ) {
487
488         fDisplay->GetPadRowPointer()->Draw3D();
489       }
490       // show all padrows 
491       else {
492
493 #if defined(HAVE_HOMERREADER) 
494         HOMERReader* reader = (HOMERReader*)fDisplay->fReader;
495
496         char* rawID = "KPWR_LDD";
497         ULong_t blk;
498         blk = reader->FindBlockNdx( rawID, " CPT",0xFFFFFFFF );
499
500         Int_t NRawPoints = 0;
501         TPolyMarker3D* pm = new  TPolyMarker3D( );
502         pm->SetBit(kCanDelete);
503         pm->SetMarkerColor(51); 
504     
505         while ( blk != ~(ULong_t)0 ) {
506           
507 #if DEBUG
508           printf( "Found raw data block %lu\n", blk );
509 #endif
510           // Check for corrupt data
511 #if HOMER_VERSION >= 2
512           AliHLTUInt64_t corruptFlag = reader->GetBlockStatusFlags( blk );
513           if (corruptFlag & 0x00000001) {
514             LOG(AliHLTTPCLog::kError,"AliHLTTPCDisplayMain::ReadData","Block status flags") << "Data block is corrupt"<<ENDLOG; 
515             continue;
516           }
517 #endif
518
519           unsigned long rawDataBlock = (unsigned long) reader->GetBlockData( blk );
520           unsigned long rawDataLen = reader->GetBlockDataLength( blk );
521
522           ULong_t spec = reader->GetBlockDataSpec( blk );
523           Int_t patch = AliHLTTPCDefinitions::GetMinPatchNr( spec );
524           Int_t slice = AliHLTTPCDefinitions::GetMinSliceNr( spec );
525 #if DEBUG
526           printf( "Raw data found for slice %u - patch %u\n", slice, patch );
527 #endif  
528
529           // slice should(not) be displayed
530           if (!fDisplay->GetDisplaySlice(slice)) continue;
531           
532
533 #if defined(HAVE_TPC_MAPPING)
534           AliHLTTPCDigitReaderRaw digitReader(0);
535
536           bool readValue = true;
537           Int_t rowOffset = 0;
538     
539           // Initialize RAW DATA
540           Int_t firstRow = AliHLTTPCTransform::GetFirstRow(patch);
541           Int_t lastRow = AliHLTTPCTransform::GetLastRow(patch);
542           
543           // Outer sector, patches 2, 3, 4, 5 -  start counting in patch 2 with row 0
544           if ( patch >= 2 ) rowOffset = AliHLTTPCTransform::GetFirstRow( 2 );
545
546           // Initialize block for reading packed data
547           void* tmpDataBlock = (void*) rawDataBlock;
548           digitReader.InitBlock(tmpDataBlock,rawDataLen,firstRow,lastRow,patch,slice);
549
550           readValue = digitReader.Next();
551
552           if (!readValue){      
553             LOG(AliHLTTPCLog::kError,"AliHLTTPCDisplayPadRow::Fill","Read first value") << "No value in data block" << ENDLOG;
554             continue;
555           }
556
557
558           //      blk = reader->FindBlockNdx( rawID, " CPT", 0xFFFFFFFF, blk+1 );
559           //      continue;
560
561           // Fill 3D Raw Data
562           while ( readValue ){ 
563             
564             Int_t row = digitReader.GetRow() + rowOffset;
565
566             UChar_t pad = digitReader.GetPad();
567             UShort_t time = digitReader.GetTime();
568             UInt_t charge = digitReader.GetSignal();
569             if ( charge < 50 ) {
570               // read next value
571               readValue = digitReader.Next();
572       
573               if(!readValue) break; //No more value
574               continue;
575             }
576             Float_t xyz[3];
577
578             // Transform raw coordinates to local coordinates
579             AliHLTTPCTransform::RawHLT2Global(xyz, slice, row, pad, time);
580
581             NRawPoints++;
582             pm->SetNextPoint((Double_t)xyz[0],(Double_t)xyz[1],(Double_t)xyz[2]);
583
584             //  printf("%u points\n",NRawPoints);
585
586             // read next value
587             readValue = digitReader.Next();
588       
589             if(!readValue) break; //No more value
590           } // end  while ( readValue ){ 
591
592
593 #else //! defined(HAVE_TPC_MAPPING)
594           HLTFatal("DigitReaderRaw not available - check your build");
595 #endif //defined(HAVE_TPC_MAPPING)
596           
597           blk = reader->FindBlockNdx( rawID, " CPT", 0xFFFFFFFF, blk+1 );
598
599         } // end while ( blk != ~(ULong_t)0 ) {
600         pm->Draw(); 
601 #else
602     HLTFatal("HOMER reader not available");
603 #endif // defined(HAVE_HOMERREADER) 
604
605
606         //////////////////////////////
607
608
609
610
611       } // end show all padrows
612
613     }
614
615     //--------------------------------------------------------------------------------------------
616     // DRAW 3D 
617     //--------------------------------------------------------------------------------------------
618     // v->ZoomView(0,3);
619     // v->Draw();   
620
621     fDisplay->GetCanvas3D()->SetFillColor(fDisplay->GetBackColor());
622
623     if ( !fDisplay->GetKeepView() ){ 
624         fDisplay->GetCanvas3D()->SetTheta(fDisplay->GetTheta());
625         fDisplay->GetCanvas3D()->SetPhi(fDisplay->GetPhi());
626     }
627
628     fDisplay->GetCanvas3D()->Modified();
629     fDisplay->GetCanvas3D()->Update();
630 }
631
632 //____________________________________________________________________________________________________
633 void AliHLTTPCDisplay3D::DrawGeomSector(Int_t sector) {  
634   /*
635     Char_t fname[256];
636     Int_t realsector = sector;// % 18;
637   
638     if (realsector < 10){
639         sprintf(fname,"LS0%d",realsector);
640         fGeom->GetNode(fname)->SetLineColor(fDisplay->GetLineColor());
641         fGeom->GetNode(fname)->Draw("same");
642         sprintf(fname,"US0%d",realsector);
643         fGeom->GetNode(fname)->SetLineColor(fDisplay->GetLineColor()); 
644         fGeom->GetNode(fname)->Draw("same");
645     }
646     else {
647         sprintf(fname,"LS%d",realsector);
648         fGeom->GetNode(fname)->SetLineColor(fDisplay->GetLineColor());
649         fGeom->GetNode(fname)->Draw("same");
650         sprintf(fname,"US%d",realsector);
651         fGeom->GetNode(fname)->SetLineColor(fDisplay->GetLineColor()); 
652         fGeom->GetNode(fname)->Draw("same");
653     }
654
655   */   
656 }
657
658 //____________________________________________________________________________________________________
659 void AliHLTTPCDisplay3D::LoadGeometrie(Char_t *gfile) {
660     if (gfile) {
661         TFile *file = TFile::Open(gfile);
662         if(!file) {
663             LOG(AliHLTTPCLog::kError,"AliHLTTPCDisplay3D::AliHLTDisplay","File Open") 
664               <<"Geometry file " << gfile << " does not exist!"<<ENDLOG;
665             exit(-1);
666         }
667         
668         fGeom = (TGeometry*)file->Get("AliceGeom");
669
670         file->Close();
671         //delete file;  ####
672     }
673 }