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