]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMD.h
updating ChangeLog
[u/mrichter/AliRoot.git] / FMD / AliFMD.h
1 #ifndef ALIFMD_H
2 #define ALIFMD_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights
4  * reserved. 
5  *
6  * Latest changes by Christian Holm Christensen <cholm@nbi.dk>
7  *
8  * See cxx source for full Copyright notice                               
9  */
10 /** @file    AliFMD.h
11     @author  Christian Holm Christensen <cholm@nbi.dk>
12     @date    Sun Mar 26 17:59:37 2006
13     @brief   Declaration of AliFMD detector driver 
14 */
15 /** @mainpage ALICE FMD Off-line code 
16     
17     @b Contents 
18     - @ref intro 
19     - @ref structure 
20       - @ref base     (see also @ref FMD_base)
21       - @ref sim      (see also @ref FMD_sim)
22       - @ref rec      (see also @ref FMD_rec)
23       - @ref flow     (see also @ref FMD_flow)
24       - @ref ana      (see also @ref FMD_ana)
25       - @ref util     (see also @ref FMD_util)
26     - @ref script     (see also @ref FMD_script)
27     - @ref quick
28     - @ref authors
29     
30     @section intro Introduction:
31
32     This file contains a short overview of the FMD code.   It is by no 
33     means authoritative  - only the code is good for that.  However,
34     I'll try to explain things a bit here. 
35
36     @section structure Structure:
37
38     There are 4 libraries build for the FMD.  These are 
39
40     - libFMDbase: This contains the basic stuff, like data classes and
41       managers.
42     - libFMDsim: This contains code used by the simulation only.  That
43       includes the detector class AliFMD and it's derivatives.  It
44       also contains the code for building the geometry, as well as the
45       digitisers and raw data writer.
46     - libFMDrec: Code needed for the reconstruction.  This include the
47       reconstruction code itself, as well as a raw data reader.
48     - libFMDutil: This is a special library that contains various
49       utility classes for the FMD expert/developer/tester.  It
50       includes code to read all data produced by the FMD, a simple
51       event display, and code to make fake calibration and alignment
52       data.  This library is @e not loaded by aliroot
53       automatically. The user has to load it herself:
54       @code 
55       gSystem->Load("libFMDutil.so");
56       @endcode
57     The content of these libraries is detailed more below. 
58     
59     @subsection base libFMDbase:
60
61     This currently (18th or March 2006) contains the classes 
62
63     - AliFMDBaseDigit, AliFMDDigit, AliFMDSDigit: Base class for
64       digits, real digits, and summable digits.  The are all data
65       classes that holds the ADC value(s) for a single strip.
66     - AliFMDBoolMap: A map of booleans, one for each strip.
67     - AliFMDUShortMap: A map of unsigned short integers, one for each
68       strip.
69     - AliFMDDetector, AliFMD1, AliFMD2, AliFMD3: 3 classes that holds
70       the parameters for each of the 3 FMD sub-detectors.  These
71       derive from AliFMDDetector, and are managed by AliFMDGeometry.
72       Each of these also contains the translation from sensor
73       reference frame to global reference frame.
74     - AliFMDRing: Parameters for the FMD rings (inner and outer type).
75       These are managed by AliFMDGeometry.
76     - AliFMDGeometry: Manager of FMD geometry data. All code queries
77       this manager for geometry parameters, so that the data is always
78       consistent.  
79     - AliFMDParameters: Manager of FMD parameters, like calibration
80       parameters.  This class fetches data from CDB if necessary.
81       All code queries this manager for parameters, so that the data
82       is always consistent. 
83     - AliFMDCalibPedestal, AliFMDCalibGain, AliFMDCalibSampleRate,
84       AliFMDAltroMapping: Containers of calibration parameters.  These
85       correspond to the pedestal and its width, the gain of each
86       strip, the oversampling rate used in read-out, and the map from
87       hardware address to detector.
88     - AliFMDAltroIO, AliFMDAltroReader, AliFMDAltroWriter: Low-level
89       classes to do input/output on ALTRO formated buffers. 
90
91   
92
93     @subsection sim libFMDsim:
94
95     This currently (18th or March 2006) contains the classes 
96
97     - AliFMDEdepMap: Cache of energy deposited and total number of
98       hits for each strip.  The digitiser AliFMDDigitizer uses this to
99       store simulation data before making digits.
100     - AliFMDHit: A hit in the FMD active elements, as described by the
101       simulation back-end during transport.
102     - AliFMD, AliFMDv0, AliFMDv1: Simulation drivers for the FMD.
103       AliFMD is the base class. AliFMDv0 corresponds to a simulation
104       where no hits are created, but the material distribution is
105       right.  AliFMDv1 is like AliFMDv0, except that hits are
106       produced.
107     - AliFMDGeometryBuilder: Build the FMD geometry in terms of TGeo
108       objects.  The information for building the geometry is retrieved
109       from AliFMDGeometry.
110     - AliFMDBaseDigitizer, AliFMDDigitizer, AliFMDSDigitizer: Base
111       class for the digitisers.  AliFMDDigitizer makes `real' digits
112       (AliFMDDigit) from hits, and AliFMDSDigitizer makes summable
113       digits from hits.
114     - AliFMDRawWriter: Writes a pseudo raw data file from the digits
115       created by the digitisers.  It uses the AliFMDAltroMapping from
116       AliFMDParameters to make the mapping from detector coordinates
117       to hardware addresses.
118
119     @subsection rec libFMDrec:
120
121     This currently (18th or March 2006) contains the classes 
122
123     - AliFMDReconstructor: Reconstruct (in a naiive way) the charged
124       particle multiplicity in the FMD strips.  This also writes an
125       AliESDFMD object to the ESD files (that class is in libESD).
126
127     - AliFMDRecPoint: Reconstructed point in the FMD.  These objects
128       are made AliFMDReconstructor. 
129
130     - AliFMDRawReader: Classes to read raw data files. 
131
132     @subsection flow libFMDflow:
133     
134     This library contains flow analysis code that works similar to 
135     histograms. 
136
137     @subsection ana libFMDanalysis:
138     
139     This library contains analysis code.
140
141
142     @subsection util libFMDutil:
143
144     This currently (18th or March 2006) contains the classes 
145
146     - AliFMDInput, AliFMDInputHits, AliFMDInputDigits,
147       AliFMDInputSDigits, AliFMDInputRecPoints: Base class, and
148       concrete classes to read in FMD generated data.  These provide a
149       simple and unified way of getting the data.  Hooks are defined
150       to process hits, tracks, digits, and reconstructed points, as
151       well as geometry and ESD data.   See for example the scripts
152       @c DrawHits.C, @c DrawHitsDigits.C, @c DrawHitsRecs.C, @c
153       DrawDigitsRecs.C in the @c FMD/scripts sub-directory.
154
155     - AliFMDDisplay: Simple event display for FMD data, including
156       hits, digits, reconstructed points and ESD data. 
157
158     - AliFMDCalibFaker, AliFMDAlignFaker: Classes to write fake (or
159       dummy) calibration and alignment  data. 
160
161
162     @section script Scripts 
163     
164     Most scripts live in @c FMD/scripts.  The notiable exceptions are 
165     @ref Simulate.C, @ref Reconstruct.C, and @ref Config.C 
166
167     @section quick Quick start 
168
169     First, install ROOT.  Then Install TGeant3: 
170     @verbatim 
171     > cd ~/
172     > mkdir alice
173     > cd alice
174     > cvs -d :pserver:cvs@root.cern.ch:/user/cvs login 
175     Password: cvs
176     > cvs -d :pserver:cvs@root.cern.ch:/user/cvs co geant3
177     > cd geant3
178     > make 
179     @endverbatim 
180
181     Now you can install AliRoot 
182     @verbatim 
183     > cd ../
184     > cvs -d :pserver:cvs@alisoft.cern.ch:/soft/cvsroot login
185     Password: <empty>
186     > cvs -d :pserver:cvs@alisoft.cern.ch:/soft/cvsroot co AliRoot
187     > cd AliRoot
188     > export ALICE_TARGET=`root-config --arch`
189     > export ALICE=${HOME}/alice
190     > export ALICE_ROOT=${ALICE}/AliRoot
191     > export ALICE_LEVEL=new
192     > export LD_LIBRARY_PATH=${ALICE_ROOT}/lib/tgt_${ALICE_TERGET}:${LD_LIBRARY_PATH}
193     > export PATH=${ALICE_ROOT}/bin/tgt_${ALICE_TERGET}:${PATH}
194     > export G3SYS=${ALICE}/geant3
195     > make 
196     @endverbatim 
197     
198     To simulate one event, do something like 
199
200     @verbatim 
201     > aliroot ${ALICE_ROOT}/FMD/Simulate.C
202     @endverbatim 
203
204     To reconstruct the generated event, do 
205     @verbatim 
206     > aliroot ${ALICE_ROOT}/FMD/Reconstruct.C
207     @endverbatim 
208
209     Now, open the file `AliESDs.root' in AliRoot, and browse through  that. 
210
211     @section authors Authors:
212
213     - Alla Maevskaya            <Alla.Maevskaia@cern.ch>        
214     - Christian Holm Christensen        <cholm@nbi.dk>
215 */
216 /** @defgroup FMD_sim Simulation */
217
218 //____________________________________________________________________
219 //
220 //  Manager class for the FMD - Base class.
221 //  AliFMDv1, AliFMDv0, and AliFMDAlla 
222 //  provides concrete implementations. 
223 //  This class is sooooo crowded
224 //
225 #ifndef ALIDETECTOR_H  
226 # include <AliDetector.h>
227 #endif
228 #ifndef ROOT_TArrayI
229 # include <TArrayI.h>
230 #endif
231 class TBranch;
232 class TClonesArray;
233 class TBrowser;
234 class TMarker3DBox;
235 class TArrayI;
236 class AliDigitizer;
237 class AliFMDHit;
238
239 //____________________________________________________________________
240 /** @class AliFMD AliFMD.h <FMD/AliFMD.h>
241     @brief Forward Multiplicity Detector based on Silicon wafers. 
242     This class is the driver for especially simulation. 
243
244     The Forward Multiplicity Detector consists of 3 sub-detectors FMD1,
245     FMD2, and FMD3, each of which has 1 or 2 rings of silicon sensors. 
246                                                           
247     This is the base class for all FMD manager classes. 
248                        
249     The actual code is done by various separate classes.   Below is
250     diagram showing the relationship between the various FMD classes
251     that handles the simulation
252
253     @verbatim
254           +----------+   +----------+   
255           | AliFMDv1 |   | AliFMDv0 |   
256           +----------+   +----------+   
257                |              |                    +-----------------+
258           +----+--------------+                 +--| AliFMDDigitizer |
259           |                                     |  +-----------------+
260           |           +---------------------+   |
261           |        +--| AliFMDBaseDigitizer |<--+
262           V     1  |  +---------------------+   |
263      +--------+<>--+                            |  +------------------+
264      | AliFMD |                                 +--| AliFMDSDigitizer |    
265      +--------+<>--+                               +------------------+       
266                1  |  +---------------------+
267                   +--| AliFMDReconstructor |
268                      +---------------------+
269     @endverbatim
270                      
271     -  AliFMD 
272        This defines the interface for the various parts of AliROOT that
273        uses the FMD, like AliFMDSimulator, AliFMDDigitizer, 
274        AliFMDReconstructor, and so on. 
275     -  AliFMDv0
276        This is a concrete implementation of the AliFMD interface. 
277        It is the responsibility of this class to create the FMD
278        geometry.
279     -  AliFMDv1 
280        This is a concrete implementation of the AliFMD interface. 
281        It is the responsibility of this class to create the FMD
282        geometry, process hits in the FMD, and serve hits and digits to
283        the various clients. 
284     -  AliFMDSimulator
285        This is the base class for the FMD simulation tasks.   The
286        simulator tasks are responsible to implment the geoemtry, and
287        process hits. 
288     -  AliFMDReconstructor
289        This is a concrete implementation of the AliReconstructor that
290        reconstructs pseudo-inclusive-multiplicities from digits (raw or
291        from simulation)
292
293     Calibration and geometry parameters are managed by separate
294     singleton managers.  These are AliFMDGeometry and
295     AliFMDParameters.  Please refer to these classes for more
296     information on these.    
297  */
298 class AliFMD : public AliDetector 
299 {
300 public:
301   /** Default constructor.  Do not use. */
302   AliFMD();
303   /** Normal constructor 
304       @param name  Name of object.
305       @param title Title of object. */
306   AliFMD(const char *name, const char *title);
307   /** Destructor */
308   virtual ~AliFMD(); 
309   /** Wheter to make a detailed geometry
310       @param use If true, make detailed geometry  */
311   void UseDetailed(Bool_t use=kTRUE) { fDetailed = use; }
312   
313   /** @{*/
314   /** @name GEometry ANd Tracking (GEANT :-) */
315   /** Define the geometry.  This is done by asking the manager
316       AliFMDGeometry to construct the geometry.  This in turn calls
317       AliFMDGeometryBuilder.   */
318   virtual void   CreateGeometry();
319   /** Create entries for alignable volumes associating the symbolic volume
320       name with the corresponding volume path. Needs to be syncronized with
321       eventual changes in the geometry.   */
322   virtual void  AddAlignableVolumes() const;
323   /** Create the tracking mediums used by the FMD.  This associates
324       the tracking mediums defined with the FMD in the
325       TVirtualMCApplication (AliMC). 
326       
327       The defined mediums are 
328       - @c FMD @c Si$   Silicon (active medium in sensors)
329       - @c FMD @c C$    Carbon fibre (support cone for FMD3 and vacuum pipe)
330       - @c FMD @c Al$   Aluminium (honeycomb support plates)
331       - @c FMD @c PCB$  Printed Circuit Board (FEE board with VA1_3)
332       - @c FMD @c Chip$ Electronics chips (currently not used)
333       - @c FMD @c Air$  Air (Air in the FMD)
334       - @c FMD @c Plastic$ Plastic (Support legs for the hybrid cards)
335   */
336   virtual void   CreateMaterials(); 
337 #if 0
338   /** 
339    * Declare tracking parameters for a medium 
340    * 
341    * Cut offs are in GeV. 
342    * @param imed                     Medium identifier
343    * @param gamma                    Cut off for tracking photons
344    * @param electron                 Cut off for tracking electrons
345    * @param neutral_hadron           Cut off for tracking neutral hadrons
346    * @param charged_hadron           Cut off for tracking charged hadrons
347    * @param muon                     Cut off for tracking muons
348    * @param electron_bremstrahlung   Cut off for tracking electron brehmstralung
349    * @param muon__bremstrahlung      Cut off for tracking muon brehmstralung
350    * @param electron_delta           Cut off for tracking delta electrons
351    * @param muon_delta               Cut off for tracking delta muons
352    * @param muon_pair                Cut off for muon->ee pair production
353    * @param annihilation             Enable annihilation
354    * @param bremstrahlung            Enable brehmstralung
355    * @param compton_scattering       Enable Compton scattering
356    * @param decay                    Enable decays
357    * @param delta_ray                Enable delta rays
358    * @param hadronic                 Enable hadronic interactions
359    * @param energy_loss              Enable energy loss
360    * @param multiple_scattering      Enable multiple scattering
361    * @param pair_production          Enable pair production
362    * @param photon_production        Enable cherenkov photon production
363    * @param rayleigh_scattering      Enable rayleigh scattering
364    */ 
365   void SetTrackingParameters(Int_t imed, 
366                              Float_t gamma,                 
367                              Float_t electron, 
368                              Float_t neutral_hadron, 
369                              Float_t charged_hadron, 
370                              Float_t muon,
371                              Float_t electron_bremstrahlung, 
372                              Float_t muon__bremstrahlung, 
373                              Float_t electron_delta,
374                              Float_t muon_delta,
375                              Float_t muon_pair,
376                              Int_t   annihilation, 
377                              Int_t   bremstrahlung, 
378                              Int_t   compton_scattering, 
379                              Int_t   decay,
380                              Int_t   delta_ray, 
381                              Int_t   hadronic, 
382                              Int_t   energy_loss, 
383                              Int_t   multiple_scattering, 
384                              Int_t   pair_production, 
385                              Int_t   photon_production, 
386                              Int_t   rayleigh_scattering);
387 #endif
388   /** Initialize this detector */
389   virtual void   Init();
390   /** This member function is called when ever a track deposites
391       energy (or similar) in an FMD tracking medium.  In this base
392       class this member function is pure abstract.   In concrete
393       sub-classes, the member function may make hits or other
394       stuff. */ 
395   virtual void   StepManager() = 0;
396   /** Called at the end of each simulation event.  If the debug level
397       is high enough a list of @e bad hits is printed. */
398   virtual void   FinishEvent();
399   /** @}*/
400   
401   /** @{ */
402   /** @name Hit and digit management */
403   /* Create Tree branches for the FMD.
404      @param opt  One of 
405      - @c H Make a branch of TClonesArray of AliFMDHit's
406      - @c D Make a branch of TClonesArray of AliFMDDigit's
407      - @c S Make a branch of TClonesArray of AliFMDSDigit's */
408   virtual void          MakeBranch(Option_t *opt=" ");
409   /** Set the TClonesArray to read hits into.
410       @param b The branch to containn the hits */
411   virtual void          SetHitsAddressBranch(TBranch *b);
412   /** Set the TClonesArray to read sdigits into.
413       @param b The branch to containn the sdigits */
414   virtual void          SetSDigitsAddressBranch(TBranch *b);
415   /** Set branch address for the Hits, Digits, and SDigits Tree. */
416   virtual void          SetTreeAddress();
417   /** Get the array of summable digits
418       @return summable digits */
419   virtual TClonesArray* SDigits() { return fSDigits; }        
420   /** Reset the array of summable digits */
421   virtual void          ResetSDigits();
422   /** Add a hit to the hits tree 
423       @param  track  Track #
424       @param  vol Volume parameters, interpreted as 
425       - ivol[0]  [UShort_t ] Detector # 
426       - ivol[1]  [Char_t   ] Ring ID 
427       - ivol[2]  [UShort_t ] Sector #
428       - ivol[3]  [UShort_t ] Strip # 
429       @param hits Hit information 
430       - hits[0]  [Float_t  ] Track's X-coordinate at hit 
431       - hits[1]  [Float_t  ] Track's Y-coordinate at hit
432       - hits[3]  [Float_t  ] Track's Z-coordinate at hit
433       - hits[4]  [Float_t  ] X-component of track's momentum             
434       - hits[5]  [Float_t  ] Y-component of track's momentum             
435       - hits[6]  [Float_t  ] Z-component of track's momentum            
436       - hits[7]  [Float_t  ] Energy deposited by track                  
437       - hits[8]  [Int_t    ] Track's particle Id # 
438       - hits[9]  [Float_t  ] Time when the track hit */
439   virtual void          AddHit(Int_t track, Int_t *vol, Float_t *hits);
440   /** Add a hit to the list
441       @param track     Track #
442       @param detector  Detector # (1, 2, or 3)                      
443       @param ring      Ring ID ('I' or 'O')
444       @param sector    Sector # (For inner/outer rings: 0-19/0-39)
445       @param strip     Strip # (For inner/outer rings: 0-511/0-255)
446       @param x         Track's X-coordinate at hit
447       @param y         Track's Y-coordinate at hit
448       @param z         Track's Z-coordinate at hit
449       @param px        X-component of track's momentum 
450       @param py        Y-component of track's momentum
451       @param pz        Z-component of track's momentum
452       @param edep      Energy deposited by track
453       @param pdg       Track's particle Id #
454       @param t         Time when the track hit 
455       @param len       Track length through the material. 
456       @param stopped   Whether track was stopped or disappeared */
457   virtual AliFMDHit*    AddHitByFields(Int_t    track, 
458                                        UShort_t detector, 
459                                        Char_t   ring, 
460                                        UShort_t sector, 
461                                        UShort_t strip, 
462                                        Float_t  x=0,
463                                        Float_t  y=0, 
464                                        Float_t  z=0,
465                                        Float_t  px=0, 
466                                        Float_t  py=0, 
467                                        Float_t  pz=0,
468                                        Float_t  edep=0,
469                                        Int_t    pdg=0,
470                                        Float_t  t=0, 
471                                        Float_t  len=0, 
472                                        Bool_t   stopped=kFALSE);
473   /** Add a digit to the Digit tree 
474       @param digits
475       - digits[0]  [UShort_t] Detector #
476       - digits[1]  [Char_t]   Ring ID
477       - digits[2]  [UShort_t] Sector #
478       - digits[3]  [UShort_t] Strip #
479       - digits[4]  [UShort_t] ADC Count 
480       - digits[5]  [Short_t]  ADC Count, -1 if not used
481       - digits[6]  [Short_t]  ADC Count, -1 if not used 
482       @param notused Not used */
483   virtual        void   AddDigit(Int_t *digits, Int_t* notused=0);
484   /** add a real digit
485       @param detector  Detector # (1, 2, or 3)                      
486       @param ring         Ring ID ('I' or 'O')
487       @param sector       Sector # (For inner/outer rings: 0-19/0-39)
488       @param strip        Strip # (For inner/outer rings: 0-511/0-255)
489       @param count1    ADC count (a 10-bit word)
490       @param count2    ADC count (a 10-bit word), or -1 if not used
491       @param count3    ADC count (a 10-bit word), or -1 if not used */
492   virtual        void   AddDigitByFields(UShort_t       detector=0, 
493                                          Char_t         ring='\0', 
494                                          UShort_t       sector=0, 
495                                          UShort_t       strip=0, 
496                                          UShort_t       count1=0, 
497                                          Short_t        count2=-1, 
498                                          Short_t        count3=-1, 
499                                          Short_t        count4=-1, 
500                                          UShort_t       nrefs=0,
501                                          Int_t*         refs=0);
502   /** Add a digit to the Digit tree 
503       @param digits
504       - digits[0]  [UShort_t] Detector #
505       - digits[1]  [Char_t]   Ring ID
506       - digits[2]  [UShort_t] Sector #
507       - digits[3]  [UShort_t] Strip #
508       - digits[4]  [Float_t]  Edep
509       - digits[5]  [UShort_t] ADC Count 
510       - digits[6]  [Short_t]  ADC Count, -1 if not used
511       - digits[7]  [Short_t]  ADC Count, -1 if not used  
512       - digits[8]  [Short_t]  ADC Count, -1 if not used
513       - digits[9]  [UShort_t] N total particles
514       - digits[10] [UShort_t] N total primary particles
515   */
516   virtual        void   AddSDigit(Int_t *digits);
517   /** add a summable digit - as coming from data
518       @param detector  Detector # (1, 2, or 3)                      
519       @param ring      Ring ID ('I' or 'O')
520       @param sector    Sector # (For inner/outer rings: 0-19/0-39)
521       @param strip     Strip # (For inner/outer rings: 0-511/0-255)
522       @param edep      Energy deposited   
523       @param count1    ADC count (a 10-bit word)
524       @param count2    ADC count (a 10-bit word), or -1 if not used 
525       @param count3    ADC count (a 10-bit word), or -1 if not used */
526   virtual        void   AddSDigitByFields(UShort_t       detector=0, 
527                                           Char_t         ring='\0', 
528                                           UShort_t       sector=0, 
529                                           UShort_t       strip=0, 
530                                           Float_t        edep=0,
531                                           UShort_t       count1=0, 
532                                           Short_t        count2=-1, 
533                                           Short_t        count3=-1,
534                                           Short_t        count4=-1, 
535                                           UShort_t       ntot=0, 
536                                           UShort_t       nprim=0,
537                                           Int_t*         refs=0);
538   /** @}*/
539
540   /** @{ */
541   /** @name Digitisation */
542   /** Create a digitizer object
543       @param manager Digitization manager
544       @return a newly allocated AliFMDDigitizer */
545   virtual AliDigitizer* CreateDigitizer(AliDigitizationInput* digInput) const;
546   /** Create AliFMDDigit's from AliFMDHit's.  This is done by creating
547       an AliFMDDigitizer object, and executing it.  */
548   virtual        void   Hits2Digits();
549   /** Create AliFMDSDigit's from AliFMDHit's.  This is done by creating
550       an AliFMDSDigitizer object, and executing it.  */
551   virtual        void   Hits2SDigits();
552   /** @}*/
553
554   /** @{ */
555   /** @name Raw data */
556   /** 
557    * Turn digits into raw data. This uses the class AliFMDRawWriter to
558    * do the job.  Please refer to that class for more information.
559    */
560   virtual void   Digits2Raw();
561   /** 
562    * Convert raw data to sdigits
563    * 
564    * @param reader Raw reader
565    * 
566    * @return @c true on success
567    */
568   virtual Bool_t Raw2SDigits(AliRawReader* reader);
569   /** @}*/
570
571   /** @{ */
572   /** 
573    * @name Utility 
574    */
575   /** 
576    * Browse this object 
577    *
578    * @param b Browser to show this object in 
579    */
580   void   Browse(TBrowser* b);
581   /** @}*/
582 protected:
583   /** Initialize hit array if not already done, and return pointert. 
584       @return Hit array */
585   TClonesArray*      HitsArray();
586   /** Initialize digit array if not already done, and return pointert. 
587       @return Digit array */
588   TClonesArray*      DigitsArray();
589   /** Initialize summable digit array if not already done, and return
590       pointert.  
591       @return Summable digit array */
592   TClonesArray*      SDigitsArray();
593
594   TClonesArray*      fSDigits;              // Summable digits
595   Int_t              fNsdigits;             // Number of digits  
596   Bool_t             fDetailed;             // Use detailed geometry
597   Bool_t             fUseOld;               // Use old approx geometry
598   Bool_t             fUseAssembly;          // Use divided volumes
599   
600   enum {
601     kSiId,                 // ID index of Si medium
602     kAirId,                // ID index of Air medium
603     kPlasticId,            // ID index of Plastic medium
604     kPcbId,                // ID index of PCB medium
605     kSiChipId,             // ID index of Si Chip medium
606     kAlId,                 // ID index of Al medium
607     kCarbonId,             // ID index of Carbon medium
608     kCopperId,             // ID index of Copper Medium
609     kKaptonId,             // ID index of Kapton Medium
610     kSteelId               // ID index of Steel medium
611   };  
612
613   TObjArray*         fBad;                  //! debugging - bad hits 
614
615 private:  
616   /** Copy constructor 
617       @param other Object to copy from */
618   AliFMD(const AliFMD& other);
619   /** Assignment operator 
620       @param other Object to assign from
621       @return Reference to this object  */
622   AliFMD& operator=(const AliFMD& other);
623
624   ClassDef(AliFMD,11)     // Base class FMD entry point
625 };
626
627 #endif
628 //____________________________________________________________________
629 //
630 // Local Variables:
631 //   mode: C++
632 // End:
633 //
634 // EOF
635 //