]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMD.h
fix requested by Laurent Aphecetche to clean up the memory properly during simulation...
[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.  These derive from
160       TTask.  
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   /** Initialize this detector */
338   virtual void   Init();
339   /** This member function is called when ever a track deposites
340       energy (or similar) in an FMD tracking medium.  In this base
341       class this member function is pure abstract.   In concrete
342       sub-classes, the member function may make hits or other
343       stuff. */ 
344   virtual void   StepManager() = 0;
345   /** Called at the end of each simulation event.  If the debug level
346       is high enough a list of @e bad hits is printed. */
347   virtual void   FinishEvent();
348   /** @}*/
349   
350   /** @{*/
351   /** @name Graphics and event display */
352   /** Draw a shaded view of the Forward multiplicity detector.  This 
353       isn't really useful anymore. */
354   virtual        void   DrawDetector() {}
355   /** @}*/
356   
357   /** @{ */
358   /** @name Hit and digit management */
359   /* Create Tree branches for the FMD.
360      @param opt  One of 
361      - @c H Make a branch of TClonesArray of AliFMDHit's
362      - @c D Make a branch of TClonesArray of AliFMDDigit's
363      - @c S Make a branch of TClonesArray of AliFMDSDigit's */
364   virtual void          MakeBranch(Option_t *opt=" ");
365   /** Set the TClonesArray to read hits into.
366       @param b The branch to containn the hits */
367   virtual void          SetHitsAddressBranch(TBranch *b);
368   /** Set the TClonesArray to read sdigits into.
369       @param b The branch to containn the sdigits */
370   virtual void          SetSDigitsAddressBranch(TBranch *b);
371   /** Set branch address for the Hits, Digits, and SDigits Tree. */
372   virtual void          SetTreeAddress();
373   /** Get the array of summable digits
374       @return summable digits */
375   virtual TClonesArray* SDigits() { return fSDigits; }        
376   /** Reset the array of summable digits */
377   virtual void          ResetSDigits();
378   /** Add a hit to the hits tree 
379       @param  track  Track #
380       @param  vol Volume parameters, interpreted as 
381       - ivol[0]  [UShort_t ] Detector # 
382       - ivol[1]  [Char_t   ] Ring ID 
383       - ivol[2]  [UShort_t ] Sector #
384       - ivol[3]  [UShort_t ] Strip # 
385       @param hits Hit information 
386       - hits[0]  [Float_t  ] Track's X-coordinate at hit 
387       - hits[1]  [Float_t  ] Track's Y-coordinate at hit
388       - hits[3]  [Float_t  ] Track's Z-coordinate at hit
389       - hits[4]  [Float_t  ] X-component of track's momentum             
390       - hits[5]  [Float_t  ] Y-component of track's momentum             
391       - hits[6]  [Float_t  ] Z-component of track's momentum            
392       - hits[7]  [Float_t  ] Energy deposited by track                  
393       - hits[8]  [Int_t    ] Track's particle Id # 
394       - hits[9]  [Float_t  ] Time when the track hit */
395   virtual void          AddHit(Int_t track, Int_t *vol, Float_t *hits);
396   /** Add a hit to the list
397       @param track     Track #
398       @param detector  Detector # (1, 2, or 3)                      
399       @param ring      Ring ID ('I' or 'O')
400       @param sector    Sector # (For inner/outer rings: 0-19/0-39)
401       @param strip     Strip # (For inner/outer rings: 0-511/0-255)
402       @param x         Track's X-coordinate at hit
403       @param y         Track's Y-coordinate at hit
404       @param z         Track's Z-coordinate at hit
405       @param px        X-component of track's momentum 
406       @param py        Y-component of track's momentum
407       @param pz        Z-component of track's momentum
408       @param edep      Energy deposited by track
409       @param pdg       Track's particle Id #
410       @param t         Time when the track hit 
411       @param len       Track length through the material. 
412       @param stopped   Whether track was stopped or disappeared */
413   virtual AliFMDHit*    AddHitByFields(Int_t    track, 
414                                        UShort_t detector, 
415                                        Char_t   ring, 
416                                        UShort_t sector, 
417                                        UShort_t strip, 
418                                        Float_t  x=0,
419                                        Float_t  y=0, 
420                                        Float_t  z=0,
421                                        Float_t  px=0, 
422                                        Float_t  py=0, 
423                                        Float_t  pz=0,
424                                        Float_t  edep=0,
425                                        Int_t    pdg=0,
426                                        Float_t  t=0, 
427                                        Float_t  len=0, 
428                                        Bool_t   stopped=kFALSE);
429   /** Add a digit to the Digit tree 
430       @param digits
431       - digits[0]  [UShort_t] Detector #
432       - digits[1]  [Char_t]   Ring ID
433       - digits[2]  [UShort_t] Sector #
434       - digits[3]  [UShort_t] Strip #
435       - digits[4]  [UShort_t] ADC Count 
436       - digits[5]  [Short_t]  ADC Count, -1 if not used
437       - digits[6]  [Short_t]  ADC Count, -1 if not used 
438       @param notused Not used */
439   virtual        void   AddDigit(Int_t *digits, Int_t* notused=0);
440   /** add a real digit
441       @param detector  Detector # (1, 2, or 3)                      
442       @param ring         Ring ID ('I' or 'O')
443       @param sector       Sector # (For inner/outer rings: 0-19/0-39)
444       @param strip        Strip # (For inner/outer rings: 0-511/0-255)
445       @param count1    ADC count (a 10-bit word)
446       @param count2    ADC count (a 10-bit word), or -1 if not used
447       @param count3    ADC count (a 10-bit word), or -1 if not used */
448   virtual        void   AddDigitByFields(UShort_t       detector=0, 
449                                          Char_t         ring='\0', 
450                                          UShort_t       sector=0, 
451                                          UShort_t       strip=0, 
452                                          UShort_t       count1=0, 
453                                          Short_t        count2=-1, 
454                                          Short_t        count3=-1, 
455                                          Short_t        count4=-1, 
456                                          UShort_t       nrefs=0,
457                                          Int_t*         refs=0);
458   /** Add a digit to the Digit tree 
459       @param digits
460       - digits[0]  [UShort_t] Detector #
461       - digits[1]  [Char_t]   Ring ID
462       - digits[2]  [UShort_t] Sector #
463       - digits[3]  [UShort_t] Strip #
464       - digits[4]  [Float_t]  Edep
465       - digits[5]  [UShort_t] ADC Count 
466       - digits[6]  [Short_t]  ADC Count, -1 if not used
467       - digits[7]  [Short_t]  ADC Count, -1 if not used  
468       - digits[8]  [Short_t]  ADC Count, -1 if not used
469       - digits[9]  [UShort_t] N total particles
470       - digits[10] [UShort_t] N total primary particles
471   */
472   virtual        void   AddSDigit(Int_t *digits);
473   /** add a summable digit - as coming from data
474       @param detector  Detector # (1, 2, or 3)                      
475       @param ring      Ring ID ('I' or 'O')
476       @param sector    Sector # (For inner/outer rings: 0-19/0-39)
477       @param strip     Strip # (For inner/outer rings: 0-511/0-255)
478       @param edep      Energy deposited   
479       @param count1    ADC count (a 10-bit word)
480       @param count2    ADC count (a 10-bit word), or -1 if not used 
481       @param count3    ADC count (a 10-bit word), or -1 if not used */
482   virtual        void   AddSDigitByFields(UShort_t       detector=0, 
483                                           Char_t         ring='\0', 
484                                           UShort_t       sector=0, 
485                                           UShort_t       strip=0, 
486                                           Float_t        edep=0,
487                                           UShort_t       count1=0, 
488                                           Short_t        count2=-1, 
489                                           Short_t        count3=-1,
490                                           Short_t        count4=-1, 
491                                           UShort_t       ntot=0, 
492                                           UShort_t       nprim=0,
493                                           Int_t*         refs=0);
494   /** @}*/
495
496   /** @{ */
497   /** @name Digitisation */
498   /** Create a digitizer object
499       @param manager Digitization manager
500       @return a newly allocated AliFMDDigitizer */
501   virtual AliDigitizer* CreateDigitizer(AliRunDigitizer* manager) const;
502   /** Create AliFMDDigit's from AliFMDHit's.  This is done by creating
503       an AliFMDDigitizer object, and executing it.  */
504   virtual        void   Hits2Digits();
505   /** Create AliFMDSDigit's from AliFMDHit's.  This is done by creating
506       an AliFMDSDigitizer object, and executing it.  */
507   virtual        void   Hits2SDigits();
508   /** @}*/
509
510   /** @{ */
511   /** @name Raw data */
512   /** 
513    * Turn digits into raw data. This uses the class AliFMDRawWriter to
514    * do the job.  Please refer to that class for more information.
515    */
516   virtual void   Digits2Raw();
517   /** 
518    * Convert raw data to sdigits
519    * 
520    * @param reader Raw reader
521    * 
522    * @return @c true on success
523    */
524   virtual Bool_t Raw2SDigits(AliRawReader* reader);
525   /** @}*/
526
527   /** @{ */
528   /** 
529    * @name Utility 
530    */
531   /** 
532    * Browse this object 
533    *
534    * @param b Browser to show this object in 
535    */
536   void   Browse(TBrowser* b);
537   /** @}*/
538 protected:
539   /** Initialize hit array if not already done, and return pointert. 
540       @return Hit array */
541   TClonesArray*      HitsArray();
542   /** Initialize digit array if not already done, and return pointert. 
543       @return Digit array */
544   TClonesArray*      DigitsArray();
545   /** Initialize summable digit array if not already done, and return
546       pointert.  
547       @return Summable digit array */
548   TClonesArray*      SDigitsArray();
549
550   TClonesArray*      fSDigits;              // Summable digits
551   Int_t              fNsdigits;             // Number of digits  
552   Bool_t             fDetailed;             // Use detailed geometry
553   Bool_t             fUseOld;               // Use old approx geometry
554   Bool_t             fUseAssembly;          // Use divided volumes
555   
556   enum {
557     kSiId,                 // ID index of Si medium
558     kAirId,                // ID index of Air medium
559     kPlasticId,            // ID index of Plastic medium
560     kPcbId,                // ID index of PCB medium
561     kSiChipId,             // ID index of Si Chip medium
562     kAlId,                 // ID index of Al medium
563     kCarbonId,             // ID index of Carbon medium
564     kCopperId,             // ID index of Copper Medium
565     kKaptonId,             // ID index of Kapton Medium
566     kSteelId               // ID index of Steel medium
567   };  
568
569   TObjArray*         fBad;                  //! debugging - bad hits 
570
571 private:  
572   /** Copy constructor 
573       @param other Object to copy from */
574   AliFMD(const AliFMD& other);
575   /** Assignment operator 
576       @param other Object to assign from
577       @return Reference to this object  */
578   AliFMD& operator=(const AliFMD& other);
579
580   ClassDef(AliFMD,11)     // Base class FMD entry point
581 };
582
583 #endif
584 //____________________________________________________________________
585 //
586 // Local Variables:
587 //   mode: C++
588 // End:
589 //
590 // EOF
591 //