]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDDisplay.h
Fixes to AliFMDDisplay and AliFMDPattern for use in DQM environment.
[u/mrichter/AliRoot.git] / FMD / AliFMDDisplay.h
1 #ifndef AliFMDDISPLAY_H
2 #define AliFMDDISPLAY_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights
4  * reserved. 
5  *
6  * See cxx source for full Copyright notice                               
7  */
8 /** @file    AliFMDDisplay.h
9     @author  Christian Holm Christensen <cholm@nbi.dk>
10     @date    Mon Mar 27 12:39:09 2006
11     @brief   FMD Event display 
12 */
13 //___________________________________________________________________
14 //
15 // The classes defined here, are utility classes for reading in data
16 // for the FMD.  They are  put in a seperate library to not polute the
17 // normal libraries.  The classes are intended to be used as base
18 // classes for customized class that do some sort of analysis on the
19 // various types of data produced by the FMD. 
20 //
21 #include "AliFMDInput.h"
22 #include <TObjArray.h>
23 #include <TTimer.h>
24 class TCanvas;
25 class TPad;
26 class TButton;
27 class TSlider;
28 class TH1;
29
30 //___________________________________________________________________
31 /** @class AliFMDDisplay 
32     @brief Utility class to visualize FMD data in geometry. 
33     @ingroup FMD_util
34  */
35 class AliFMDDisplay : public AliFMDInput
36 {
37 public:
38   /** Constructor
39       @param onlyFMD Only show the FMD
40       @param gAliceFile galice file*/
41   AliFMDDisplay(Bool_t onlyFMD=kTRUE, const char* gAliceFile="galice.root");
42   /** DTOR */
43   virtual ~AliFMDDisplay();
44   /** Singleton access function
45       @return Singleton object. */
46   static AliFMDDisplay* Instance();
47
48   /** Continue to next event */
49   void  Continue() { fWait = kFALSE; }
50   /** Run throug events as fast as possible */ 
51   void Start() { fContinous = kTRUE; fWait = kFALSE; }
52   /** Pause the processing */ 
53   void Pause() { fContinous = kFALSE; fWait = kTRUE; }
54   /** Zoom mode */
55   void  Zoom() { fZoomMode = kTRUE; }
56   /** Pick mode */
57   void  Pick() { fZoomMode = kFALSE; }
58   /** Redisplay the event */ 
59   virtual void Redisplay(); // *MENU*
60   /** Break */
61   virtual void Break();
62   /** Render in 3D */
63   virtual void Render();
64   
65   /** Change cut */
66   virtual void ChangeCut();
67   /** Change cut */
68   virtual void ChangeFactor();
69   /** Called when a mouse or similar event happens in the display. 
70       @param event Event type
71       @param px    where the event happened in pixels along X
72       @param py    where the event happened in pixels along Y */
73   virtual void  ExecuteEvent(Int_t event, Int_t px, Int_t py);
74   /** Paint into canvas 
75       @param option Not used */
76   virtual void  Paint(Option_t* option="") { (void)option; }
77
78   /** Initialize
79       @return  @c false on error */
80   virtual Bool_t Init();
81   /** Called at beginning of an event 
82       @param event Event number
83       @return @c false on error  */
84   virtual Bool_t Begin(Int_t event);
85   /** Called at end of an event 
86       @return @c false on error  */
87   virtual Bool_t End();
88   /** Visualize a hit
89       @param hit Hit
90       @param p   Track
91       @return @c false on error  */
92   virtual Bool_t ProcessHit(AliFMDHit* hit, TParticle* p);
93   /** Visualize a digit
94       @param digit Digit to draw
95       @return @c false on error  */
96   virtual Bool_t ProcessDigit(AliFMDDigit* digit);
97   /** Visualize a summable digit
98       @param sdigit Summable digit to draw
99       @return @c false on error  */
100   virtual Bool_t ProcessSDigit(AliFMDSDigit* sdigit);
101   /** Visualize a raw digit
102       @param digit Raw digit.
103       @return @c false on error  */
104   virtual Bool_t ProcessRawDigit(AliFMDDigit* digit);
105   /** Visualize a raw digit
106       @param digit Raw digit.
107       @return @c false on error  */
108   virtual Bool_t ProcessRawCalibDigit(AliFMDDigit* digit);
109   /** Visualize a reconstructed point.
110       @param recpoint Reconstructed point
111       @return @c false on error  */
112   virtual Bool_t ProcessRecPoint(AliFMDRecPoint* recpoint);
113   /** Process ESD data for the FMD.  Users should overload this to
114       deal with ESD data. 
115       @param d    Detector number (1-3)
116       @param r    Ring identifier ('I' or 'O')
117       @param s    Sector number (0-19, or 0-39)
118       @param t    Strip number (0-511, or 0-255)
119       @param eta  Psuedo-rapidity 
120       @param mult Psuedo-multiplicity 
121       @return  @c false on error  */
122   virtual Bool_t ProcessESD(UShort_t d, Char_t r, UShort_t s, UShort_t t, 
123                             Float_t eta, Float_t mult);
124   /** Look up a color index, based on the value @a x and the maximum
125       value of @a x
126       @param x   Value 
127       @param max Maximum (for example 1023 for digits)
128       @return @c false on error  */
129   virtual Int_t  LookupColor(Float_t x, Float_t min, Float_t max)  const;
130
131   /** Set range of displayed values */
132   virtual void SetCut(Float_t l=0., Float_t h=1.); //*MENU*
133 protected:
134   /** Copy constructor 
135       @param o Object to copy from  */
136   AliFMDDisplay(const AliFMDDisplay& o) 
137     : AliFMDInput(o),
138       fWait(kFALSE),
139       fMarkers(0),
140       fHits(0),
141       fCanvas(0),
142       fPad(0),
143       fButtons(0),
144       fSlider(0),
145       fFactor(0),
146       fZoomMode(0),
147       fX0(0),
148       fY0(0),
149       fX1(0),
150       fY1(0),
151       fXPixel(0),
152       fYPixel(0),
153       fOldXPixel(0),
154       fOldYPixel(0),
155       fLineDrawn(0), 
156       fOnlyFMD(kTRUE),
157       fSpec(0), 
158       fSpecCut(0),
159       fAux(0),
160       fReturn(kFALSE),
161       fContinous(kFALSE), 
162       fTimeout()
163   { } 
164   /** Assignment operator 
165       @return Reference to this object */
166   AliFMDDisplay& operator=(const AliFMDDisplay&) { return *this; } 
167   /** Add a marker to the display
168       @param x   X position
169       @param y   Y position
170       @param z   Z position
171       @param o   Object to refer to
172       @param s   Signal 
173       @param max Maximum of signal */
174   virtual void AddMarker(Float_t x, Float_t y, Float_t z, 
175                          TObject* o, Float_t s, Float_t min, Float_t max);
176   /** Add a marker to the display
177       @param det Detector
178       @param rng Ring
179       @param sec Sector 
180       @param str Strip
181       @param o   Object to refer to
182       @param s   Signal 
183       @param max Maximum of signal */
184   virtual void AddMarker(UShort_t det, Char_t rng, UShort_t sec, UShort_t str, 
185                          TObject* o, Float_t s, Float_t min, Float_t max);
186   
187   /** Show only the FMD detectors. */
188   void ShowOnlyFMD();
189   /** Make base canvas */ 
190   virtual void MakeCanvas(const char** which);
191   virtual void MakeAux();
192   virtual void DrawAux();
193   virtual void Idle();
194   virtual void AtEnd();
195   virtual Bool_t InsideCut(Float_t v, const Float_t& min, 
196                          const Float_t& max) const;
197   virtual Double_t GetADCThreshold(UShort_t d, Char_t r, 
198                                    UShort_t s, UShort_t t) const;
199   
200   static AliFMDDisplay* fgInstance; // Static instance 
201   Bool_t                fWait;      // Wait until user presses `Continue'
202   TObjArray*            fMarkers;   // Cache of markers
203   TObjArray*            fHits;      // Cache of `hits'
204   TCanvas*              fCanvas;    // Canvas to draw in 
205   TPad*                 fPad;       // View pad. 
206   TObjArray             fButtons;   // Continue button
207   TSlider*              fSlider;    // Cut slider
208   TSlider*              fFactor;    // Factor slider
209   Bool_t                fZoomMode;  // Whether we're in Zoom mode
210   Float_t               fX0;        // X at lower left corner or range 
211   Float_t               fY0;        // Y at lower left corner or range 
212   Float_t               fX1;        // X at upper right corner or range 
213   Float_t               fY1;        // Y at upper right corner or range 
214   Int_t                 fXPixel;    // X pixel of mark
215   Int_t                 fYPixel;    // Y pixel of mark
216   Int_t                 fOldXPixel; // Old x pixel of mark
217   Int_t                 fOldYPixel; // Old y pixel of mark
218   Bool_t                fLineDrawn; // Whether we're drawing a box
219   Bool_t                fOnlyFMD;   // Whether to only do FMD
220   TH1*                  fSpec;      // Spectra
221   TH1*                  fSpecCut;   // Cut spectra
222   TCanvas*              fAux;       // Aux canvas.
223   Bool_t                fReturn;    // Stop 
224   Bool_t                fContinous;
225   TTimer                fTimeout;
226   
227
228   struct Range_t { 
229     UInt_t  fNbins;
230     Float_t fLow;
231     Float_t fHigh;
232   };
233   static const Range_t fgkEdepRange;
234   static const Range_t fgkAdcRange;
235   static const Range_t fgkMultRange;
236
237   ClassDef(AliFMDDisplay,0)  // FMD specialised event display
238 };
239
240 #endif
241 //____________________________________________________________________
242 //
243 // Local Variables:
244 //   mode: C++
245 // End:
246 //
247 // EOF
248 //