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