]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDDisplay.cxx
Coding rule violation corrected.
[u/mrichter/AliRoot.git] / FMD / AliFMDDisplay.cxx
1 /**************************************************************************
2  * Copyright(c) 2004, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15 /* $Id$ */
16 /** 
17  * @file    AliFMDDisplay.cxx
18  * @author  Christian Holm Christensen <cholm@nbi.dk>
19  * @date    Mon Mar 27 12:39:09 2006
20  * @brief   FMD Event display 
21  */
22 //___________________________________________________________________
23 //
24 // The classes defined here, are utility classes for reading in data
25 // for the FMD.  They are  put in a seperate library to not polute the
26 // normal libraries.  The classes are intended to be used as base
27 // classes for customized class that do some sort of analysis on the
28 // various types of data produced by the FMD. 
29 //
30 // Latest changes by Christian Holm Christensen
31 //
32
33 #include <TSystem.h>
34 #include <TApplication.h>
35 #include <TButton.h>
36 #include <TCanvas.h>
37 #include <TGeoManager.h>
38 #include <TH1D.h>
39 #include <TMarker3DBox.h>
40 #include <TMath.h>
41 #include <TSlider.h>
42 #include <TSliderBox.h>
43 #include <TStyle.h>
44 #include <TView.h>
45 #include <TVirtualX.h>
46 #include <TVirtualViewer3D.h>
47 #include <TList.h>
48 // #include <TArrayF.h>
49 // #include <TParticle.h>
50
51 #include "AliFMDDisplay.h"      // ALIFMDDISPLAY_H
52 #include "AliFMDHit.h"          // ALIFMDHIT_H
53 #include "AliFMDDigit.h"        // ALIFMDDIGIT_H
54 #include "AliFMDSDigit.h"       // ALIFMDSDIGIT_H
55 #include "AliFMDRecPoint.h"     // ALIFMDRECPOINT_H
56 #include "AliFMDGeometry.h"     // ALIFMDGEOMETRY_H
57 #include "AliFMDParameters.h"   // ALIFMDPARAMETERS_H
58 #include "AliFMDRawReader.h"    // ALIFMDRAWREADER_H
59 #include <AliESDFMD.h>          // ALIESDFMD_H
60 // #include <AliLog.h>
61 #include "AliFMDDebug.h" // Better debug macros
62
63 //____________________________________________________________________
64 ClassImp(AliFMDDisplay)
65 #if 0
66   ; // This is here to keep Emacs for indenting the next line
67 #endif
68
69 //____________________________________________________________________
70 AliFMDDisplay* AliFMDDisplay::fgInstance = 0;
71
72 //____________________________________________________________________
73 const AliFMDDisplay::Range_t AliFMDDisplay::fgkEdepRange = {  100, 0.,    2. };
74 const AliFMDDisplay::Range_t AliFMDDisplay::fgkAdcRange  = { 1024, 0., 1023. };
75 const AliFMDDisplay::Range_t AliFMDDisplay::fgkMultRange = {  500, 0.,   20. };
76
77   
78 //____________________________________________________________________
79 AliFMDDisplay* 
80 AliFMDDisplay::Instance()
81 {
82   // Return static instance 
83   // If the instance does not exist
84   // it is not created!
85   return fgInstance;
86 }
87
88 //____________________________________________________________________
89 AliFMDDisplay::~AliFMDDisplay()
90 {
91   // Destructor. 
92   // Cleans 
93   // up 
94   if (fMarkers) {
95     fMarkers->Delete();
96     delete fMarkers;
97   }
98   if (fHits) {
99     fHits->Clear();
100     delete fHits;
101   }
102   if (fPad)     delete fPad;
103   fButtons.Delete();
104   if (fSlider)  delete fSlider;
105   if (fCanvas)  delete fCanvas;
106 }
107   
108 //____________________________________________________________________
109 AliFMDDisplay::AliFMDDisplay(Bool_t onlyFMD, const char* gAliceFile)
110   : AliFMDInput(gAliceFile),
111     fWait(kFALSE),
112     fMarkers(0),
113     fHits(0),
114     fCanvas(0), 
115     fPad(0), 
116     fButtons(0),
117     fSlider(0),
118     fFactor(0),
119     fZoomMode(kFALSE),
120     fX0(0),
121     fY0(0),
122     fX1(0),
123     fY1(0),
124     fXPixel(0),
125     fYPixel(0),
126     fOldXPixel(0),
127     fOldYPixel(0),
128     fLineDrawn(0),
129     fOnlyFMD(onlyFMD),
130     fSpec(0), 
131     fSpecCut(0),
132     fAux(0), 
133     fReturn(kFALSE), 
134     fContinous(kFALSE),
135     fTimeout("gApplication->StopIdleing()", 10), 
136     fInitialMin(0), 
137     fInitialMax(1), 
138     fInitialFactor(3/10.)
139 {
140   // Constructor of an FMD display object. 
141   // Must be called 
142   // before Instance 
143   SetName("AliFMDDisplay");
144   SetTitle("3D Display of various kinds of FMD data");
145   AddLoad(kGeometry);
146   if (fgInstance) delete fgInstance;
147   fgInstance = this;
148 }
149
150 //____________________________________________________________________
151 void           
152 AliFMDDisplay::MakeCanvas(const char** which)
153 {
154   // Make a canvas 
155   // Parameters: 
156   //   which   Which button to put up. 
157   gStyle->SetPalette(1);
158   // gStyle->SetCanvasPreferGL(kTRUE);
159   Double_t y1 = .10;
160   Int_t    w  = 700;
161   fCanvas = new TCanvas(Form("gl%s", GetName()), 
162                         Form("%s - Display", GetTitle()), 
163                         w, Int_t(w / (1-y1)));
164   fCanvas->SetFillColor(1);
165   fCanvas->ToggleEventStatus();
166   fCanvas->cd();
167   fPad = new TPad("glview", "3DView", 0.0, y1, 1.0, 1.0, 1, 0, 0);
168   fPad->Draw();
169   
170   const char** p = which;
171   const char*  m;
172   Int_t        n  = 0;
173   Int_t        j  = 0;
174   while (*(p++)) n++;
175   AliFMDDebug(1, ("Got %d buttons", n));
176   if (n <= 0) return;
177
178   Double_t yb = 0;
179   Double_t xb = 1;
180   fCanvas->cd();
181   if (TESTBIT(fTreeMask, kDigits) || 
182       TESTBIT(fTreeMask, kRawCalib) || 
183       TESTBIT(fTreeMask, kRaw)) { 
184     yb = .05;
185     xb = .66;
186     fFactor = new TSlider("pedFactor", "Pedestal Factor", xb+.01, 0, 1, yb);
187     fFactor->SetMethod("AliFMDDisplay::Instance()->ChangeFactor()");
188     fFactor->SetRange(fInitialFactor, 1);
189     fFactor->Draw();
190     fFactor->SetMinimum(fInitialFactor);
191     TSliderBox *sbox = 
192       static_cast<TSliderBox*>(fFactor->GetListOfPrimitives()->
193                                FindObject("TSliderBox"));
194     if (sbox) { 
195       sbox->SetToolTipText("Adjust the noise suppression factor by moving "
196                            "lower limit");
197     }
198   }
199   if (TESTBIT(fTreeMask, kHits)    || 
200       TESTBIT(fTreeMask, kESD)     || 
201       TESTBIT(fTreeMask, kDigits)  || 
202       TESTBIT(fTreeMask, kSDigits) || 
203       TESTBIT(fTreeMask, kRaw)     ||
204       TESTBIT(fTreeMask, kRawCalib)) {
205     yb = .05;
206     fSlider = new TSlider("genCut", "Multiplicity cut", 0, 0, xb, yb);
207     fSlider->SetMethod("AliFMDDisplay::Instance()->ChangeCut()");
208     fSlider->SetRange(fInitialMin,fInitialMax);
209     fSlider->Draw();
210     TSliderBox *sbox = 
211       static_cast<TSliderBox*>(fSlider->GetListOfPrimitives()->
212                                FindObject("TSliderBox"));
213     if (sbox) { 
214       sbox->SetToolTipText("Adjust lower and upper limit on data signal");
215     }
216   }
217   // fCanvas->Modified();
218   // fCanvas->Update();
219   // fCanvas->cd();
220   Float_t      x0 = 0;
221   Float_t      dx = 1. / n;
222   p               = which;
223   while ((m = *(p++))) {
224     fCanvas->cd();
225     AliFMDDebug(1, ("Adding button %s", m));
226     TButton* b = new TButton(m, Form("AliFMDDisplay::Instance()->%s()", m),
227                              x0, yb, TMath::Min(x0 + dx,.999F), y1);
228     b->Draw();
229     fButtons.Add(b);
230     x0 += dx;
231     j++;
232   }
233 }
234
235 //____________________________________________________________________
236 void           
237 AliFMDDisplay::ShowOnlyFMD()
238 {
239   // Show only the FMD 
240   // Do not show 
241   // other volumes 
242   if (!fGeoManager) return;
243   static bool once = false;
244   if (once) return;
245   once = true;
246   AliInfo("Will only show the FMD");
247   TGeoVolume* top = gGeoManager->GetTopVolume();
248   top->InvisibleAll(kTRUE);
249   TGeoIterator next(top);
250   TGeoNode* node;
251   TGeoVolume* v = 0;
252   // Bool_t hasFMD1 = kFALSE;
253   // Bool_t hasFMD2 = kFALSE;
254   // Bool_t hasFMD3 = kFALSE;
255   AliFMDDebug(1, ("Getting material FMD_Si$"));
256   TGeoMaterial* si   = gGeoManager->GetMaterial("FMD_Si$");      // kRed 
257   AliFMDDebug(1, ("Getting material FMD_Carbon$"));
258   TGeoMaterial* c    = gGeoManager->GetMaterial("FMD_Carbon$");  // kGray
259   AliFMDDebug(1, ("Getting material FMD_Aluminum$"));
260   TGeoMaterial* al   = gGeoManager->GetMaterial("FMD_Aluminum$");// kGray-2
261   AliFMDDebug(1, ("Getting material FMD_Copper$"));
262   TGeoMaterial* cu   = gGeoManager->GetMaterial("FMD_Copper$"); // kGreen-2
263   AliFMDDebug(1, ("Getting material FMD_PCB$"));
264   TGeoMaterial* pcb  = gGeoManager->GetMaterial("FMD_PCB$");    // kGreen+2
265   AliFMDDebug(1, ("Getting material FMD_PCB$"));
266   TGeoMaterial* chip = gGeoManager->GetMaterial("FMD_Si Chip$");// kGreen+2
267   TObjArray     toshow;
268   while ((node = static_cast<TGeoNode*>(next()))) {
269     const char* name = node->GetName();
270     if (!name) continue;
271     if (!(v = node->GetVolume())) continue;
272
273     if (name[0] == 'F') {
274       TGeoMaterial* m   = (v->IsAssembly() ? 0 : v->GetMaterial());
275       Int_t         col = -1;
276       if      (m == si)   col = kRed;
277       else if (m == c)    col = kGray;
278       else if (m == al)   col = kYellow+4;  
279       else if (m == cu)   col = kRed+6;   
280       else if (m == pcb)  col = kGreen+2; 
281       else if (m == chip) col = kGreen+4; 
282       if (col >= 0) { 
283         v->SetLineColor(col);
284         v->SetFillColor(col);
285       }
286       if (name[2] == 'M' && (name[3] == 'T' || name[3] == 'B')) {
287         // Virtual Master half-ring volume - top-level
288         // Int_t det = node->GetNumber();
289         /* switch (det) {
290            case 1: hasFMD1 = true; break;
291            case 2: hasFMD2 = true; break;
292            case 3: hasFMD3 = true; break;
293            default: continue;
294            } */
295         toshow.Add(v);
296       }
297       else if (name[3] == 'V' && (name[2] == 'T' || name[2] == 'B')) 
298         toshow.Add(v); // Virtual Half-ring, bare detectors
299       else if (name[3] == 'H' && (name[2] == 'F' || name[2] == 'B')) 
300         toshow.Add(v); // Virtual Hybrid container 
301       else if (name[2] == 'S' && name[3] == 'U') 
302         toshow.Add(v); // Virtual support structre 
303       // else if (name[3] == 'H' && (name[2] == 'F' || name[2] == 'B')) 
304       //  toshow.Add(v); // Virtual Hybrid container 
305     }
306     v->SetVisibility(kFALSE);
307     v->SetVisDaughters(kFALSE);
308     v->InvisibleAll(kTRUE);
309   }
310   TIter i(&toshow);
311   while ((v = static_cast<TGeoVolume*>(i()))) {
312     if (!v->IsAssembly())
313       v->SetVisibility(kTRUE);
314     v->InvisibleAll(kFALSE);
315     v->SetVisDaughters(kTRUE);
316     
317   }  
318 }
319
320     
321 //____________________________________________________________________
322 void           
323 AliFMDDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py) 
324 {
325   // Execute an event on canvas 
326   // Parameters: 
327   //   event   What happened
328   //   px, py  Pixel coordinates 
329   if (px == 0 && py == 0) return;
330   if (!fZoomMode && fPad->GetView()) {
331     fPad->GetView()->ExecuteRotateView(event, px, py);
332     return;
333   }
334   fPad->SetCursor(kCross);
335   switch (event) {
336   case kButton1Down: 
337     fPad->TAttLine::Modify();
338     fX0        = fPad->AbsPixeltoX(px);
339     fY0        = fPad->AbsPixeltoY(py);
340     fXPixel    = fOldXPixel = px;
341     fYPixel    = fOldYPixel = py;
342     fLineDrawn = kFALSE;
343     return;
344   case kButton1Motion:
345     if (fLineDrawn) 
346       gVirtualX->DrawBox(fXPixel, fYPixel, fOldXPixel, fOldYPixel, 
347                          TVirtualX::kHollow);
348     fOldXPixel = px;
349     fOldYPixel = py;
350     fLineDrawn = kTRUE;
351     gVirtualX->DrawBox(fXPixel, fYPixel, fOldXPixel, fOldYPixel, 
352                        TVirtualX::kHollow);
353     return;
354   case kButton1Up:
355     fPad->GetCanvas()->FeedbackMode(kFALSE);
356     if (px == fXPixel || py == fYPixel) return;
357     fX1 = fPad->AbsPixeltoX(px);
358     fY1 = fPad->AbsPixeltoY(py);
359     if (fX1 < fX0) std::swap(fX0, fX1); 
360     if (fY1 < fY0) std::swap(fY0, fY1); 
361     fPad->Range(fX0, fY0, fX1, fY1);
362     fPad->Modified();
363     return;
364   }
365 }
366
367 //____________________________________________________________________
368 Bool_t 
369 AliFMDDisplay::Init()
370 {
371   // Initialize.  GEt transforms and such,
372   // so that we can draw thins properly 
373   // Returns true on success 
374   if (!AliFMDInput::Init()) return kFALSE;
375   AliFMDGeometry* geom = AliFMDGeometry::Instance();
376   geom->Init();
377   geom->InitTransformations();
378   if (TESTBIT(fTreeMask, kDigits) || TESTBIT(fTreeMask, kRaw)) 
379     AliFMDParameters::Instance()->Init();
380
381   fMarkers = new TObjArray;
382   fHits    = new TObjArray;
383   fMarkers->SetOwner(kTRUE);
384   fHits->SetOwner(kFALSE);
385   return kTRUE;
386 }
387
388 //____________________________________________________________________
389 void
390 AliFMDDisplay::MakeAux()
391 {
392   // MAke the aux canvas 
393   // This is used to display spectra
394   // etc, 
395   const Range_t* range = 0;
396   if      (TESTBIT(fTreeMask, kESD))      range = &fgkMultRange;
397   else if (TESTBIT(fTreeMask, kRawCalib)) range = &fgkMultRange;
398   else if (TESTBIT(fTreeMask, kDigits))   range = &fgkAdcRange;
399   else if (TESTBIT(fTreeMask, kSDigits))  range = &fgkAdcRange;
400   else if (TESTBIT(fTreeMask, kRaw))      range = &fgkAdcRange;
401   else if (TESTBIT(fTreeMask, kHits))     range = &fgkEdepRange;
402   if (!range) return;
403   
404   if (!fAux) {
405     fAux = new TCanvas(Form("aux_%s", GetName()), 
406                        Form("Aux - %s", GetTitle()));
407     fAux->SetLogy();
408     fAux->SetFillColor(kWhite);
409     fAux->SetBorderMode(0);
410     fAux->SetBorderSize(0);
411     Float_t dBin = (range->fHigh - range->fLow) / range->fNbins;
412     fSpec = new TH1D("spec", "Spectra", range->fNbins, 
413                      range->fLow-dBin/2, range->fHigh+dBin/2);
414     fSpecCut = static_cast<TH1*>(fSpec->Clone("specCut"));
415     TString xTitle((TESTBIT(fTreeMask, kRawCalib) || 
416                     TESTBIT(fTreeMask, kESD)) ? "#Delta E/#Delta E_{mip}" :
417                    (TESTBIT(fTreeMask, kDigits) ||
418                     TESTBIT(fTreeMask, kSDigits) || 
419                     TESTBIT(fTreeMask, kRaw)) ? "ADC [counts]" :
420                    TESTBIT(fTreeMask, kHits) ? "Hits" : "signal");
421     fSpec->SetXTitle(xTitle.Data());
422     fSpec->SetYTitle("events");
423     fSpec->SetFillColor(2);
424     fSpec->SetFillStyle(3001);
425     fSpecCut->SetXTitle(xTitle.Data());
426     fSpecCut->SetYTitle("events");
427     fSpecCut->SetFillColor(4);
428     fSpecCut->SetFillStyle(3001);
429   }
430   else {
431     fSpec->Reset();
432     fSpecCut->Reset();
433   }  
434 }
435
436 //____________________________________________________________________
437 void
438 AliFMDDisplay::DrawAux()
439 {
440   // Draw in the Aux the canvas 
441   // For example draw the spectra 
442   // or such stuff 
443   if (!fAux) return;
444   fAux->cd();
445   fAux->Clear();
446   fAux->SetLogy(fSpec->GetMaximum() > 10);
447   fSpec->Draw();
448   fSpecCut->Draw("same");
449   fAux->Modified();
450   fAux->Update();
451   fAux->cd();
452 }
453
454 //____________________________________________________________________
455 Bool_t 
456 AliFMDDisplay::Begin(Int_t event) 
457 {
458   // Begin of event.  Make canvas is not already done 
459   // Parameters: 
460   //   event   The event number 
461   if (!fCanvas) {
462     const char* m[] = { "Continue", 
463                         "Break", 
464                         "Zoom", 
465                         "Pick", 
466                         "Redisplay", 
467                         "Render", 
468                         0 }; 
469     MakeCanvas(m);
470   }
471   MakeAux();
472   fReturn = kFALSE;
473   
474   // AliInfo("Clearing canvas");
475   // fCanvas->Clear();
476   if (!fGeoManager) {
477     Warning("End", "No geometry manager");
478     return kFALSE;
479   }
480   AliFMDDebug(1, ("Drawing geometry"));
481   fPad->cd();
482   fGeoManager->GetTopVolume()->Draw();
483   if (fOnlyFMD) ShowOnlyFMD();
484   AliFMDDebug(1, ("Adjusting view"));
485   Int_t irep;
486   if (fPad->GetView()) {
487     fPad->GetView()->SetView(-200, -40, 80, irep);
488     fPad->GetView()->Zoom();
489     fPad->Modified();
490     fPad->cd();
491   }
492   return AliFMDInput::Begin(event);
493 }
494
495 //____________________________________________________________________
496 void
497 AliFMDDisplay::AtEnd() 
498 {
499   // Called at of the event. 
500   // Draw stuff. 
501   // Draw spectrum. 
502   fPad->cd();
503   fMarkers->Draw();
504   fPad->cd();
505   AppendPad();
506   fPad->cd();
507   DrawAux();
508 }
509
510 //____________________________________________________________________
511 void
512 AliFMDDisplay::Idle() 
513 {
514   // Idle loop. 
515   // Sends the ROOT loop into the idle loop,
516   // so that we can go on. 
517   fWait = kTRUE;
518   if (fContinous) fTimeout.Start(10, kTRUE);
519   while (fWait) {
520     gApplication->StartIdleing();
521     gSystem->InnerLoop();
522     gApplication->StopIdleing();
523     if (fContinous) break;
524   }
525   AliFMDDebug(3, ("After idle loop"));
526   if (fMarkers) fMarkers->Delete();
527   if (fHits)    fHits->Clear();
528   AliFMDDebug(3, ("After clearing caches"));
529 }
530
531 //____________________________________________________________________
532 Bool_t 
533 AliFMDDisplay::End()
534 {
535   // End of event.  Draw everything 
536   AtEnd();
537   Idle();
538   if (fReturn) return kFALSE;
539   return AliFMDInput::End();
540 }
541
542 //____________________________________________________________________
543 Int_t
544 AliFMDDisplay::LookupColor(Float_t x, Float_t min, Float_t max) const
545 {
546   // Look-up color.  
547   // Get a colour from the  current palette depending 
548   // on the ratio x/max 
549   Float_t range  = (max-min);
550   Float_t l      = fSlider->GetMinimum();
551   Float_t h      = fSlider->GetMaximum();
552   if (l == h) { l = 0; h = 1; }
553   Float_t cmin   = range * l;
554   Float_t cmax   = range * h;
555   Float_t crange = (cmax-cmin);
556   Int_t   idx    = Int_t((x-cmin) / crange * gStyle->GetNumberOfColors());
557   return gStyle->GetColorPalette(idx);
558
559
560 //____________________________________________________________________
561 void
562 AliFMDDisplay::SetCut(Float_t l, Float_t h) 
563 {
564   // Change the cut on the slider. 
565   fInitialMin = l;
566   fInitialMax = h;
567   if (!fSlider) return;
568   fSlider->SetMinimum(fInitialMin);
569   fSlider->SetMaximum(fInitialMax);
570   ChangeCut();
571 }
572
573 //____________________________________________________________________
574 void
575 AliFMDDisplay::SetFactor(Float_t f) 
576 {
577   // Change the cut on the slider. 
578   fInitialFactor = f / 10;
579   if (!fFactor) return;
580   fFactor->SetMinimum(fInitialFactor);
581   ChangeFactor();
582 }
583
584 //____________________________________________________________________
585 void
586 AliFMDDisplay::ChangeCut() 
587 {
588   // Change the cut on the slider. 
589   // The factor depends on what is 
590   // drawn in the AUX canvas
591   AliInfo(Form("Range is now %7.5f - %7.5f", fSlider->GetMinimum(), 
592                fSlider->GetMaximum()));
593   if ((TESTBIT(fTreeMask, kESD)     || 
594        TESTBIT(fTreeMask, kDigits)  || 
595        TESTBIT(fTreeMask, kSDigits) || 
596        TESTBIT(fTreeMask, kRaw)     ||
597        TESTBIT(fTreeMask, kRawCalib))) {
598     Float_t l = fSlider->GetMinimum();
599     Float_t h = fSlider->GetMaximum();
600     l         = 1024 * l + 0;
601     h         = 1024 * h + 0;
602     AliInfo(Form("ADC range is now %4d - %4d", int(l), int(h)));
603   }
604   Redisplay();
605 }
606 //____________________________________________________________________
607 void
608 AliFMDDisplay::ChangeFactor() 
609 {
610   // Change the cut on the slider. 
611   // The factor depends on what is 
612   // drawn in the AUX canvas
613   AliInfo(Form("Noise factor is now %4.1f", 10 * fFactor->GetMinimum()));
614   Redisplay();
615 }
616
617 //____________________________________________________________________
618 void
619 AliFMDDisplay::Redisplay()
620 {
621   // Redisplay stuff. 
622   // Redraw markers, hits, 
623   // spectra 
624   if (fMarkers) fMarkers->Delete();
625   if (fHits)    fHits->Clear();
626   if (fSpec)    fSpec->Reset();
627   if (fSpecCut) fSpecCut->Reset();
628   Event();
629   AtEnd();
630 }
631 //____________________________________________________________________
632 void
633 AliFMDDisplay::Break()
634 {
635   // Redisplay stuff. 
636   // Redraw markers, hits, 
637   // spectra 
638   if (fMarkers) fMarkers->Delete();
639   if (fHits)    fHits->Clear();
640   if (fSpec)    fSpec->Reset();
641   if (fSpecCut) fSpecCut->Reset();
642   fReturn = kTRUE;
643   fWait   = kFALSE;
644 }
645 //____________________________________________________________________
646 void
647 AliFMDDisplay::Render()
648 {
649   fPad->cd();
650   TVirtualViewer3D* viewer = fPad->GetViewer3D("ogl");
651   if (!viewer) return;
652 }
653
654 //____________________________________________________________________
655 void
656 AliFMDDisplay::AddMarker(Float_t x, Float_t y, Float_t z,
657                          TObject* o, Float_t s, Float_t min, Float_t max)
658 {
659   // Add a marker to the display
660   //
661   //    det     Detector
662   //    rng     Ring
663   //    sec     Sector 
664   //    str     Strip
665   //    o       Object to refer to
666   //    s       Signal 
667   //    max     Maximum of signal 
668   //
669   Float_t  size  = .1;
670   Float_t  zsize = (s - min) / (max-min) * 10;
671   Float_t  r     = TMath::Sqrt(x * x + y * y);
672   Float_t  theta = TMath::ATan2(r, z);
673   Float_t  phi   = TMath::ATan2(y, x);
674   Float_t  rz    = z + (z < 0 ? 1 : -1) * zsize;
675   TMarker3DBox* marker = new  TMarker3DBox(x,y,rz,size,size,zsize,theta,phi);
676   if (o) marker->SetRefObject(o);
677   marker->SetLineColor(LookupColor(s, min, max));
678   fMarkers->Add(marker);
679 }
680 //____________________________________________________________________
681 void
682 AliFMDDisplay::AddMarker(UShort_t det, Char_t rng, UShort_t sec, UShort_t str,
683                          TObject* o, Float_t s, Float_t min, Float_t max)
684 {
685   // Add a marker to the display
686   //
687   //    det     Detector
688   //    rng     Ring
689   //    sec     Sector 
690   //    str     Strip
691   //    o       Object to refer to
692   //    s       Signal 
693   //    max     Maximum of signal 
694   //
695   AliFMDGeometry*   geom = AliFMDGeometry::Instance();
696   Double_t x, y, z;
697   geom->Detector2XYZ(det, rng, sec, str, x, y, z);
698   AddMarker(x,y,z,o,s,min,max);
699 }
700   
701 //____________________________________________________________________
702 Bool_t 
703 AliFMDDisplay::InsideCut(Float_t val, const Float_t& min, 
704                          const Float_t& max) const
705 {
706   // 
707   // Whether a point is inside 
708   // 
709   // Parameters:
710   //    v   Point
711   //    min Minimum
712   //    max Maximum
713   // 
714   // Return:
715   //    true if @a v is inside cut 
716   //
717   Float_t r = max - min;
718   Float_t l = fSlider->GetMinimum();
719   Float_t h = fSlider->GetMaximum();
720   if (l == h) { l = 0; h = 1; }
721   if (val < r * l + min || val > r * h + min) { 
722     AliFMDDebug(2, ("Value %f is outside cut %f - %f (range %f - %f)", 
723                     val, min+r*l, min+r*h, min, max));
724     return kFALSE;
725   }
726   return kTRUE;
727 }
728
729
730 //____________________________________________________________________
731 Bool_t 
732 AliFMDDisplay::ProcessHit(AliFMDHit* hit, TParticle* /* p */) 
733 {
734   // Process a hit. 
735   // Parameters: 
736   //   hit   Hit data
737   static const Float_t rMin  = fgkEdepRange.fLow;
738   static const Float_t rMax  = fgkEdepRange.fHigh;
739
740   if (!hit) { AliError("No hit");   return kFALSE; }
741   // if (!p)   { AliError("No track"); return kFALSE; }
742   Float_t  edep  = hit->Edep();
743
744   if (fHits)                        fHits->Add(hit);
745   if (fSpec)                        fSpec->Fill(edep);
746   if (!InsideCut(edep, rMin, rMax)) return kTRUE;
747   if (fSpecCut)                     fSpecCut->Fill(edep);
748   
749   AddMarker(hit->X(), hit->Y(), hit->Z(), hit, edep, rMin, rMax);
750   return kTRUE;
751 }
752
753 //____________________________________________________________________
754 Bool_t 
755 AliFMDDisplay::ProcessDigit(AliFMDDigit* digit)
756 {
757   // Process a digit 
758   // Parameters: 
759   //   digit Digit information 
760   static const Float_t rMin  = fgkAdcRange.fLow;
761   static const Float_t rMax  = fgkAdcRange.fHigh;
762
763   if (!digit) { AliError("No digit");   return kFALSE; }
764
765   UShort_t det           =  digit->Detector();
766   Char_t   ring          =  digit->Ring();
767   UShort_t sec           =  digit->Sector();
768   UShort_t str           =  digit->Strip();
769   Double_t threshold     =  GetADCThreshold(det, ring, sec, str);
770   Float_t  counts        =  digit->Counts();
771   if (threshold < 0) { counts += -threshold; threshold = 0; }
772
773   AliFMDDebug(10, ("FMD%d%c[%02d,%03d] counts %4d threshold %4d", 
774                    det, ring, sec, str, Int_t(counts), Int_t(threshold)));
775   if (fHits)                                    fHits->Add(digit);
776   if (fSpec)                                    fSpec->Fill(counts);
777   if (!InsideCut(counts-threshold, rMin, rMax)) return kTRUE;
778   if (fSpecCut)                                 fSpecCut->Fill(counts);
779   
780
781   AddMarker(det, ring, sec, str, digit, counts, rMin, rMax);
782   return kTRUE;
783 }
784
785 //____________________________________________________________________
786 Bool_t 
787 AliFMDDisplay::ProcessSDigit(AliFMDSDigit* sdigit)
788 {
789   // Process a sdigit 
790   // Parameters: 
791   //   sdigit Digit information 
792   static const Float_t rMin  = fgkAdcRange.fLow;
793   static const Float_t rMax  = fgkAdcRange.fHigh;
794
795   if (!sdigit) { AliError("No sdigit");   return kFALSE; }
796   
797   UShort_t det           =  sdigit->Detector();
798   Char_t   ring          =  sdigit->Ring();
799   UShort_t sec           =  sdigit->Sector();
800   UShort_t str           =  sdigit->Strip();
801   Float_t  counts        =  sdigit->Counts();
802
803   if (fHits)                          fHits->Add(sdigit);
804   if (fSpec)                          fSpec->Fill(counts);
805   if (!InsideCut(counts, rMin, rMax)) return kTRUE;
806   if (fSpecCut)                       fSpecCut->Fill(counts);
807   
808   AddMarker(det, ring, sec, str, sdigit, counts, rMin, rMax);
809   return kTRUE;
810 }
811
812 //____________________________________________________________________
813 Bool_t 
814 AliFMDDisplay::ProcessRawDigit(AliFMDDigit* digit)
815 {
816   // PRocess raw data 
817   // Parameters: 
818   //   digit Digit information 
819   AliFMDDebug(50, ("Forwarding call of ProcessRaw to ProcessDigit "
820                   "for FMD%d%c[%02d,%03d] %d", 
821                   digit->Detector(), digit->Ring(), digit->Sector(), 
822                   digit->Strip(), digit->Counts()));
823   return ProcessDigit(digit);
824 }
825
826 //____________________________________________________________________
827 Bool_t 
828 AliFMDDisplay::ProcessRawCalibDigit(AliFMDDigit* digit)
829 {
830   // Process a digit 
831   // Parameters: 
832   //   digit Digit information 
833   static const Float_t rMin  = fgkMultRange.fLow;
834   static const Float_t rMax  = fgkMultRange.fHigh;
835   static const Float_t aMin  = fgkAdcRange.fLow;
836   static const Float_t aMax  = fgkAdcRange.fHigh;
837
838   if (!digit) { AliError("No digit");   return kFALSE; }
839
840   AliFMDParameters* parm = AliFMDParameters::Instance();
841   UShort_t det           =  digit->Detector();
842   Char_t   ring          =  digit->Ring();
843   UShort_t sec           =  digit->Sector();
844   UShort_t str           =  digit->Strip();
845   Double_t gain          =  parm->GetPulseGain(det, ring, sec, str);
846   Double_t ped           =  parm->GetPedestal(det, ring, sec, str);
847   Double_t threshold     =  GetADCThreshold(det, ring, sec, str);
848   Float_t  counts        =  digit->Counts();
849   if (threshold < 0) { counts += -threshold; threshold = 0; ped = 0; }
850
851   //   Double_t edep = ((counts * parm->GetEdepMip() / 
852   //                (gain * parm->GetDACPerMIP()));
853   Double_t mult = (counts-ped) / (gain * parm->GetDACPerMIP());
854   if (gain < 0.1 || gain > 10) mult = 0;
855   
856
857   AliFMDDebug(10, ("FMD%d%c[%02d,%03d] adc %4d "
858                    "(threshold %4d, gain %6.3f) -> mult %7.4f", 
859                    det, ring, sec, str, int(counts), int(threshold), 
860                    gain, mult));
861   if (fHits)                                    fHits->Add(digit);
862   if (fSpec)                                    fSpec->Fill(mult);
863   if (!InsideCut(counts-threshold, aMin, aMax)) return kTRUE;
864   if (fSpecCut)                                 fSpecCut->Fill(mult);
865   
866   if (mult >= 0) 
867     AddMarker(det, ring, sec, str, digit, mult, rMin, rMax);
868   return kTRUE;
869 }
870
871 //____________________________________________________________________
872 Bool_t 
873 AliFMDDisplay::ProcessRecPoint(AliFMDRecPoint* recpoint)
874 {
875   // Process reconstructed point 
876   // Parameters: 
877   //  recpoint  Reconstructed multiplicity/energy
878   static const Float_t rMin  = fgkMultRange.fLow;
879   static const Float_t rMax  = fgkMultRange.fHigh;
880
881
882   if (!recpoint) { AliError("No recpoint");   return kFALSE; }
883
884   if (!InsideCut(recpoint->Particles(), rMin, rMax)) return kTRUE;
885
886   if (fHits) fHits->Add(recpoint);
887   AddMarker(recpoint->Detector(), recpoint->Ring(), recpoint->Sector(),  
888             recpoint->Strip(), recpoint, recpoint->Particles(), rMin, rMax);
889   return kTRUE;
890 }
891
892 //____________________________________________________________________
893 Bool_t 
894 AliFMDDisplay::ProcessESD(UShort_t det, Char_t rng, UShort_t sec, UShort_t str,
895                           Float_t, Float_t mult)
896 {
897   // Process data from ESD 
898   // Parameters 
899   //   det,rng,sec,str   Detector coordinates. 
900   //   mult              Multiplicity. 
901   static const Float_t rMin  = fgkMultRange.fLow;
902   static const Float_t rMax  = fgkMultRange.fHigh;
903   
904   Double_t cmult = mult;
905   if (fSpec) fSpec->Fill(cmult);
906   if (!InsideCut(cmult, rMin, rMax) || cmult == AliESDFMD::kInvalidMult) 
907     return kTRUE;
908
909   AddMarker(det,rng,sec,str, 0, cmult, rMin, rMax);
910
911   if (fSpecCut) fSpecCut->Fill(cmult);
912
913   return kTRUE;
914 }
915
916 //____________________________________________________________________
917 Double_t
918 AliFMDDisplay::GetADCThreshold(UShort_t d, Char_t r, 
919                                UShort_t s, UShort_t t) const
920 {
921   // 
922   // Get the ADC threshold
923   // 
924   // Parameters:
925   //    d Detector
926   //    r Ring 
927   //    s Sector 
928   //    t Strip 
929   // 
930   // Return:
931   //    The threshold 
932   //
933   AliFMDParameters* parm = AliFMDParameters::Instance();
934   Double_t ped           =  parm->GetPedestal(d,r, s, t);
935   Double_t pedW          =  parm->GetPedestalWidth(d,r, s, t);
936   Double_t threshold     = 0;
937   if (fFMDReader && fFMDReader->IsZeroSuppressed(d-1))
938     threshold = - fFMDReader->NoiseFactor(d-1) * pedW;
939   else 
940     threshold = ped + pedW * 10 * fFactor->GetMinimum();
941   AliFMDDebug(10, ("FMD%d%c[%2d,%3d] ped: %f +/- %f [factor: %f-%f]", 
942                    d, r, s, t, ped, pedW, 
943                    fFactor->GetMinimum(), fFactor->GetMaximum()));
944   if (threshold > fgkAdcRange.fHigh) threshold = fgkAdcRange.fHigh;
945   return threshold;
946 }
947
948 //____________________________________________________________________
949 //
950 // EOF
951 //
952