ROOT logo
#ifndef ALIFMD_H
#define ALIFMD_H
/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights
 * reserved. 
 *
 * Latest changes by Christian Holm Christensen <cholm@nbi.dk>
 *
 * See cxx source for full Copyright notice                               
 */
/** @file    AliFMD.h
    @author  Christian Holm Christensen <cholm@nbi.dk>
    @date    Sun Mar 26 17:59:37 2006
    @brief   Declaration of AliFMD detector driver 
*/
/** @mainpage ALICE FMD Off-line code 
    
    @b Contents 
    - @ref intro 
    - @ref structure 
      - @ref base     (see also @ref FMD_base)
      - @ref sim      (see also @ref FMD_sim)
      - @ref rec      (see also @ref FMD_rec)
      - @ref flow     (see also @ref FMD_flow)
      - @ref ana      (see also @ref FMD_ana)
      - @ref util     (see also @ref FMD_util)
    - @ref script     (see also @ref FMD_script)
    - @ref quick
    - @ref authors
    
    @section intro Introduction:

    This file contains a short overview of the FMD code.   It is by no 
    means authoritative  - only the code is good for that.  However,
    I'll try to explain things a bit here. 

    @section structure Structure:

    There are 4 libraries build for the FMD.  These are 

    - libFMDbase: This contains the basic stuff, like data classes and
      managers.
    - libFMDsim: This contains code used by the simulation only.  That
      includes the detector class AliFMD and it's derivatives.  It
      also contains the code for building the geometry, as well as the
      digitisers and raw data writer.
    - libFMDrec: Code needed for the reconstruction.  This include the
      reconstruction code itself, as well as a raw data reader.
    - libFMDutil: This is a special library that contains various
      utility classes for the FMD expert/developer/tester.  It
      includes code to read all data produced by the FMD, a simple
      event display, and code to make fake calibration and alignment
      data.  This library is @e not loaded by aliroot
      automatically. The user has to load it herself:
      @code 
      gSystem->Load("libFMDutil.so");
      @endcode
    The content of these libraries is detailed more below. 
    
    @subsection base libFMDbase:

    This currently (18th or March 2006) contains the classes 

    - AliFMDBaseDigit, AliFMDDigit, AliFMDSDigit: Base class for
      digits, real digits, and summable digits.  The are all data
      classes that holds the ADC value(s) for a single strip.
    - AliFMDBoolMap: A map of booleans, one for each strip.
    - AliFMDUShortMap: A map of unsigned short integers, one for each
      strip.
    - AliFMDDetector, AliFMD1, AliFMD2, AliFMD3: 3 classes that holds
      the parameters for each of the 3 FMD sub-detectors.  These
      derive from AliFMDDetector, and are managed by AliFMDGeometry.
      Each of these also contains the translation from sensor
      reference frame to global reference frame.
    - AliFMDRing: Parameters for the FMD rings (inner and outer type).
      These are managed by AliFMDGeometry.
    - AliFMDGeometry: Manager of FMD geometry data. All code queries
      this manager for geometry parameters, so that the data is always
      consistent.  
    - AliFMDParameters: Manager of FMD parameters, like calibration
      parameters.  This class fetches data from CDB if necessary.
      All code queries this manager for parameters, so that the data
      is always consistent. 
    - AliFMDCalibPedestal, AliFMDCalibGain, AliFMDCalibSampleRate,
      AliFMDAltroMapping: Containers of calibration parameters.  These
      correspond to the pedestal and its width, the gain of each
      strip, the oversampling rate used in read-out, and the map from
      hardware address to detector.
    - AliFMDAltroIO, AliFMDAltroReader, AliFMDAltroWriter: Low-level
      classes to do input/output on ALTRO formated buffers. 

  

    @subsection sim libFMDsim:

    This currently (18th or March 2006) contains the classes 

    - AliFMDEdepMap: Cache of energy deposited and total number of
      hits for each strip.  The digitiser AliFMDDigitizer uses this to
      store simulation data before making digits.
    - AliFMDHit: A hit in the FMD active elements, as described by the
      simulation back-end during transport.
    - AliFMD, AliFMDv0, AliFMDv1: Simulation drivers for the FMD.
      AliFMD is the base class. AliFMDv0 corresponds to a simulation
      where no hits are created, but the material distribution is
      right.  AliFMDv1 is like AliFMDv0, except that hits are
      produced.
    - AliFMDGeometryBuilder: Build the FMD geometry in terms of TGeo
      objects.  The information for building the geometry is retrieved
      from AliFMDGeometry.
    - AliFMDBaseDigitizer, AliFMDDigitizer, AliFMDSDigitizer: Base
      class for the digitisers.  AliFMDDigitizer makes `real' digits
      (AliFMDDigit) from hits, and AliFMDSDigitizer makes summable
      digits from hits.
    - AliFMDRawWriter: Writes a pseudo raw data file from the digits
      created by the digitisers.  It uses the AliFMDAltroMapping from
      AliFMDParameters to make the mapping from detector coordinates
      to hardware addresses.

    @subsection rec libFMDrec:

    This currently (18th or March 2006) contains the classes 

    - AliFMDReconstructor: Reconstruct (in a naiive way) the charged
      particle multiplicity in the FMD strips.  This also writes an
      AliESDFMD object to the ESD files (that class is in libESD).

    - AliFMDRecPoint: Reconstructed point in the FMD.  These objects
      are made AliFMDReconstructor. 

    - AliFMDRawReader: Classes to read raw data files. 

    @subsection flow libFMDflow:
    
    This library contains flow analysis code that works similar to 
    histograms. 

    @subsection ana libFMDanalysis:
    
    This library contains analysis code.


    @subsection util libFMDutil:

    This currently (18th or March 2006) contains the classes 

    - AliFMDInput, AliFMDInputHits, AliFMDInputDigits,
      AliFMDInputSDigits, AliFMDInputRecPoints: Base class, and
      concrete classes to read in FMD generated data.  These provide a
      simple and unified way of getting the data.  Hooks are defined
      to process hits, tracks, digits, and reconstructed points, as
      well as geometry and ESD data.   See for example the scripts
      @c DrawHits.C, @c DrawHitsDigits.C, @c DrawHitsRecs.C, @c
      DrawDigitsRecs.C in the @c FMD/scripts sub-directory.

    - AliFMDDisplay: Simple event display for FMD data, including
      hits, digits, reconstructed points and ESD data. 

    - AliFMDCalibFaker, AliFMDAlignFaker: Classes to write fake (or
      dummy) calibration and alignment 	data. 


    @section script Scripts 
    
    Most scripts live in @c FMD/scripts.  The notiable exceptions are 
    @ref Simulate.C, @ref Reconstruct.C, and @ref Config.C 

    @section quick Quick start 

    First, install ROOT.  Then Install TGeant3: 
    @verbatim 
    > cd ~/
    > mkdir alice
    > cd alice
    > cvs -d :pserver:cvs@root.cern.ch:/user/cvs login 
    Password: cvs
    > cvs -d :pserver:cvs@root.cern.ch:/user/cvs co geant3
    > cd geant3
    > make 
    @endverbatim 

    Now you can install AliRoot 
    @verbatim 
    > cd ../
    > cvs -d :pserver:cvs@alisoft.cern.ch:/soft/cvsroot login
    Password: <empty>
    > cvs -d :pserver:cvs@alisoft.cern.ch:/soft/cvsroot co AliRoot
    > cd AliRoot
    > export ALICE_TARGET=`root-config --arch`
    > export ALICE=${HOME}/alice
    > export ALICE_ROOT=${ALICE}/AliRoot
    > export ALICE_LEVEL=new
    > export LD_LIBRARY_PATH=${ALICE_ROOT}/lib/tgt_${ALICE_TERGET}:${LD_LIBRARY_PATH}
    > export PATH=${ALICE_ROOT}/bin/tgt_${ALICE_TERGET}:${PATH}
    > export G3SYS=${ALICE}/geant3
    > make 
    @endverbatim 
    
    To simulate one event, do something like 

    @verbatim 
    > aliroot ${ALICE_ROOT}/FMD/Simulate.C
    @endverbatim 

    To reconstruct the generated event, do 
    @verbatim 
    > aliroot ${ALICE_ROOT}/FMD/Reconstruct.C
    @endverbatim 

    Now, open the file `AliESDs.root' in AliRoot, and browse through  that. 

    @section authors Authors:

    - Alla Maevskaya		<Alla.Maevskaia@cern.ch>	
    - Christian Holm Christensen 	<cholm@nbi.dk>
*/
/** @defgroup FMD_sim Simulation */

//____________________________________________________________________
//
//  Manager class for the FMD - Base class.
//  AliFMDv1, AliFMDv0, and AliFMDAlla 
//  provides concrete implementations. 
//  This class is sooooo crowded
//
#ifndef ALIDETECTOR_H  
# include <AliDetector.h>
#endif
#ifndef ROOT_TArrayI
# include <TArrayI.h>
#endif
class TBranch;
class TClonesArray;
class TBrowser;
class TMarker3DBox;
class TArrayI;
class AliDigitizer;
class AliFMDHit;

//____________________________________________________________________
/** @class AliFMD AliFMD.h <FMD/AliFMD.h>
    @brief Forward Multiplicity Detector based on Silicon wafers. 
    This class is the driver for especially simulation. 

    The Forward Multiplicity Detector consists of 3 sub-detectors FMD1,
    FMD2, and FMD3, each of which has 1 or 2 rings of silicon sensors. 
                                                          
    This is the base class for all FMD manager classes. 
                       
    The actual code is done by various separate classes.   Below is
    diagram showing the relationship between the various FMD classes
    that handles the simulation

    @verbatim
          +----------+   +----------+   
          | AliFMDv1 |   | AliFMDv0 |   
          +----------+   +----------+   
               |              |                    +-----------------+
          +----+--------------+                 +--| AliFMDDigitizer |
          |                                     |  +-----------------+
          |           +---------------------+   |
          |        +--| AliFMDBaseDigitizer |<--+
          V     1  |  +---------------------+   |
     +--------+<>--+                            |  +------------------+
     | AliFMD |                                 +--| AliFMDSDigitizer |    
     +--------+<>--+                               +------------------+       
	       1  |  +---------------------+
	          +--| AliFMDReconstructor |
	  	     +---------------------+
    @endverbatim
		     
    -  AliFMD 
       This defines the interface for the various parts of AliROOT that
       uses the FMD, like AliFMDSimulator, AliFMDDigitizer, 
       AliFMDReconstructor, and so on. 
    -  AliFMDv0
       This is a concrete implementation of the AliFMD interface. 
       It is the responsibility of this class to create the FMD
       geometry.
    -  AliFMDv1 
       This is a concrete implementation of the AliFMD interface. 
       It is the responsibility of this class to create the FMD
       geometry, process hits in the FMD, and serve hits and digits to
       the various clients. 
    -  AliFMDSimulator
       This is the base class for the FMD simulation tasks.   The
       simulator tasks are responsible to implment the geoemtry, and
       process hits. 
    -  AliFMDReconstructor
       This is a concrete implementation of the AliReconstructor that
       reconstructs pseudo-inclusive-multiplicities from digits (raw or
       from simulation)

    Calibration and geometry parameters are managed by separate
    singleton managers.  These are AliFMDGeometry and
    AliFMDParameters.  Please refer to these classes for more
    information on these.    
 */
class AliFMD : public AliDetector 
{
public:
  /** Default constructor.  Do not use. */
  AliFMD();
  /** Normal constructor 
      @param name  Name of object.
      @param title Title of object. */
  AliFMD(const char *name, const char *title);
  /** Destructor */
  virtual ~AliFMD(); 
  /** Wheter to make a detailed geometry
      @param use If true, make detailed geometry  */
  void UseDetailed(Bool_t use=kTRUE) { fDetailed = use; }
  
  /** @{*/
  /** @name GEometry ANd Tracking (GEANT :-) */
  /** Define the geometry.  This is done by asking the manager
      AliFMDGeometry to construct the geometry.  This in turn calls
      AliFMDGeometryBuilder.   */
  virtual void   CreateGeometry();
  /** Create entries for alignable volumes associating the symbolic volume
      name with the corresponding volume path. Needs to be syncronized with
      eventual changes in the geometry.   */
  virtual void  AddAlignableVolumes() const;
  /** Create the tracking mediums used by the FMD.  This associates
      the tracking mediums defined with the FMD in the
      TVirtualMCApplication (AliMC). 
      
      The defined mediums are 
      -	@c FMD @c Si$	Silicon (active medium in sensors)
      -	@c FMD @c C$	Carbon fibre (support cone for FMD3 and vacuum pipe)
      -	@c FMD @c Al$	Aluminium (honeycomb support plates)
      -	@c FMD @c PCB$	Printed Circuit Board (FEE board with VA1_3)
      -	@c FMD @c Chip$	Electronics chips (currently not used)
      -	@c FMD @c Air$	Air (Air in the FMD)
      -	@c FMD @c Plastic$ Plastic (Support legs for the hybrid cards)
  */
  virtual void   CreateMaterials(); 
#if 0
  /** 
   * Declare tracking parameters for a medium 
   * 
   * Cut offs are in GeV. 
   * @param imed                     Medium identifier
   * @param gamma                    Cut off for tracking photons
   * @param electron                 Cut off for tracking electrons
   * @param neutral_hadron           Cut off for tracking neutral hadrons
   * @param charged_hadron           Cut off for tracking charged hadrons
   * @param muon                     Cut off for tracking muons
   * @param electron_bremstrahlung   Cut off for tracking electron brehmstralung
   * @param muon__bremstrahlung      Cut off for tracking muon brehmstralung
   * @param electron_delta           Cut off for tracking delta electrons
   * @param muon_delta               Cut off for tracking delta muons
   * @param muon_pair                Cut off for muon->ee pair production
   * @param annihilation             Enable annihilation
   * @param bremstrahlung            Enable brehmstralung
   * @param compton_scattering       Enable Compton scattering
   * @param decay                    Enable decays
   * @param delta_ray                Enable delta rays
   * @param hadronic                 Enable hadronic interactions
   * @param energy_loss              Enable energy loss
   * @param multiple_scattering      Enable multiple scattering
   * @param pair_production          Enable pair production
   * @param photon_production        Enable cherenkov photon production
   * @param rayleigh_scattering      Enable rayleigh scattering
   */ 
  void SetTrackingParameters(Int_t imed, 
			     Float_t gamma,                 
			     Float_t electron, 
			     Float_t neutral_hadron, 
			     Float_t charged_hadron, 
			     Float_t muon,
			     Float_t electron_bremstrahlung, 
			     Float_t muon__bremstrahlung, 
			     Float_t electron_delta,
			     Float_t muon_delta,
			     Float_t muon_pair,
			     Int_t   annihilation, 
			     Int_t   bremstrahlung, 
			     Int_t   compton_scattering, 
			     Int_t   decay,
			     Int_t   delta_ray, 
			     Int_t   hadronic, 
			     Int_t   energy_loss, 
			     Int_t   multiple_scattering, 
			     Int_t   pair_production, 
			     Int_t   photon_production, 
			     Int_t   rayleigh_scattering);
#endif
  /** Initialize this detector */
  virtual void   Init();
  /** This member function is called when ever a track deposites
      energy (or similar) in an FMD tracking medium.  In this base
      class this member function is pure abstract.   In concrete
      sub-classes, the member function may make hits or other
      stuff. */ 
  virtual void   StepManager() = 0;
  /** Called at the end of each simulation event.  If the debug level
      is high enough a list of @e bad hits is printed. */
  virtual void   FinishEvent();
  /** @}*/
  
  /** @{ */
  /** @name Hit and digit management */
  /* Create Tree branches for the FMD.
     @param opt  One of 
     - @c H Make a branch of TClonesArray of AliFMDHit's
     - @c D Make a branch of TClonesArray of AliFMDDigit's
     - @c S Make a branch of TClonesArray of AliFMDSDigit's */
  virtual void          MakeBranch(Option_t *opt=" ");
  /** Set the TClonesArray to read hits into.
      @param b The branch to containn the hits */
  virtual void          SetHitsAddressBranch(TBranch *b);
  /** Set the TClonesArray to read sdigits into.
      @param b The branch to containn the sdigits */
  virtual void          SetSDigitsAddressBranch(TBranch *b);
  /** Set branch address for the Hits, Digits, and SDigits Tree. */
  virtual void          SetTreeAddress();
  /** Get the array of summable digits
      @return summable digits */
  virtual TClonesArray* SDigits() { return fSDigits; }        
  /** Reset the array of summable digits */
  virtual void          ResetSDigits();
  /** Add a hit to the hits tree 
      @param  track  Track #
      @param  vol Volume parameters, interpreted as 
      - ivol[0]  [UShort_t ] Detector # 
      - ivol[1]	 [Char_t   ] Ring ID 
      - ivol[2]	 [UShort_t ] Sector #
      - ivol[3]	 [UShort_t ] Strip # 
      @param hits Hit information 
      - hits[0]	 [Float_t  ] Track's X-coordinate at hit 
      - hits[1]	 [Float_t  ] Track's Y-coordinate at hit
      - hits[3]  [Float_t  ] Track's Z-coordinate at hit
      - hits[4]  [Float_t  ] X-component of track's momentum 	       	 
      - hits[5]	 [Float_t  ] Y-component of track's momentum	       	 
      - hits[6]	 [Float_t  ] Z-component of track's momentum	       	
      - hits[7]	 [Float_t  ] Energy deposited by track		       	
      - hits[8]	 [Int_t    ] Track's particle Id # 
      - hits[9]	 [Float_t  ] Time when the track hit */
  virtual void          AddHit(Int_t track, Int_t *vol, Float_t *hits);
  /** Add a hit to the list
      @param track     Track #
      @param detector  Detector # (1, 2, or 3)                      
      @param ring      Ring ID ('I' or 'O')
      @param sector    Sector # (For inner/outer rings: 0-19/0-39)
      @param strip     Strip # (For inner/outer rings: 0-511/0-255)
      @param x	       Track's X-coordinate at hit
      @param y	       Track's Y-coordinate at hit
      @param z	       Track's Z-coordinate at hit
      @param px	       X-component of track's momentum 
      @param py	       Y-component of track's momentum
      @param pz	       Z-component of track's momentum
      @param edep      Energy deposited by track
      @param pdg       Track's particle Id #
      @param t	       Time when the track hit 
      @param len       Track length through the material. 
      @param stopped   Whether track was stopped or disappeared */
  virtual AliFMDHit*    AddHitByFields(Int_t    track, 
				       UShort_t detector, 
				       Char_t   ring, 
				       UShort_t sector, 
				       UShort_t strip, 
				       Float_t  x=0,
				       Float_t  y=0, 
				       Float_t  z=0,
				       Float_t  px=0, 
				       Float_t  py=0, 
				       Float_t  pz=0,
				       Float_t  edep=0,
				       Int_t    pdg=0,
				       Float_t  t=0, 
				       Float_t  len=0, 
				       Bool_t   stopped=kFALSE);
  /** Add a digit to the Digit tree 
      @param digits
      - digits[0]  [UShort_t] Detector #
      - digits[1]  [Char_t]   Ring ID
      - digits[2]  [UShort_t] Sector #
      - digits[3]  [UShort_t] Strip #
      - digits[4]  [UShort_t] ADC Count 
      - digits[5]  [Short_t]  ADC Count, -1 if not used
      - digits[6]  [Short_t]  ADC Count, -1 if not used 
      @param notused Not used */
  virtual        void   AddDigit(Int_t *digits, Int_t* notused=0);
  /** add a real digit
      @param detector  Detector # (1, 2, or 3)                      
      @param ring	  Ring ID ('I' or 'O')
      @param sector	  Sector # (For inner/outer rings: 0-19/0-39)
      @param strip	  Strip # (For inner/outer rings: 0-511/0-255)
      @param count1    ADC count (a 10-bit word)
      @param count2    ADC count (a 10-bit word), or -1 if not used
      @param count3    ADC count (a 10-bit word), or -1 if not used */
  virtual        void   AddDigitByFields(UShort_t       detector=0, 
					 Char_t         ring='\0', 
					 UShort_t       sector=0, 
					 UShort_t       strip=0, 
					 UShort_t       count1=0, 
					 Short_t        count2=-1, 
					 Short_t        count3=-1, 
					 Short_t        count4=-1, 
					 UShort_t       nrefs=0,
					 Int_t*		refs=0);
  /** Add a digit to the Digit tree 
      @param digits
      - digits[0]  [UShort_t] Detector #
      - digits[1]  [Char_t]   Ring ID
      - digits[2]  [UShort_t] Sector #
      - digits[3]  [UShort_t] Strip #
      - digits[4]  [Float_t]  Edep
      - digits[5]  [UShort_t] ADC Count 
      - digits[6]  [Short_t]  ADC Count, -1 if not used
      - digits[7]  [Short_t]  ADC Count, -1 if not used  
      - digits[8]  [Short_t]  ADC Count, -1 if not used
      - digits[9]  [UShort_t] N total particles
      - digits[10] [UShort_t] N total primary particles
  */
  virtual        void   AddSDigit(Int_t *digits);
  /** add a summable digit - as coming from data
      @param detector  Detector # (1, 2, or 3)                      
      @param ring      Ring ID ('I' or 'O')
      @param sector    Sector # (For inner/outer rings: 0-19/0-39)
      @param strip     Strip # (For inner/outer rings: 0-511/0-255)
      @param edep      Energy deposited   
      @param count1    ADC count (a 10-bit word)
      @param count2    ADC count (a 10-bit word), or -1 if not used 
      @param count3    ADC count (a 10-bit word), or -1 if not used */
  virtual        void   AddSDigitByFields(UShort_t       detector=0, 
					  Char_t         ring='\0', 
					  UShort_t       sector=0, 
					  UShort_t       strip=0, 
					  Float_t        edep=0,
					  UShort_t       count1=0, 
					  Short_t        count2=-1, 
					  Short_t        count3=-1,
					  Short_t        count4=-1, 
					  UShort_t       ntot=0, 
					  UShort_t       nprim=0,
					  Int_t*         refs=0);
  /** @}*/

  /** @{ */
  /** @name Digitisation */
  /** Create a digitizer object
      @param manager Digitization manager
      @return a newly allocated AliFMDDigitizer */
  virtual AliDigitizer* CreateDigitizer(AliDigitizationInput* digInput) const;
  /** Create AliFMDDigit's from AliFMDHit's.  This is done by creating
      an AliFMDDigitizer object, and executing it.  */
  virtual        void   Hits2Digits();
  /** Create AliFMDSDigit's from AliFMDHit's.  This is done by creating
      an AliFMDSDigitizer object, and executing it.  */
  virtual        void   Hits2SDigits();
  /** @}*/

  /** @{ */
  /** @name Raw data */
  /** 
   * Turn digits into raw data. This uses the class AliFMDRawWriter to
   * do the job.  Please refer to that class for more information.
   */
  virtual void   Digits2Raw();
  /** 
   * Convert raw data to sdigits
   * 
   * @param reader Raw reader
   * 
   * @return @c true on success
   */
  virtual Bool_t Raw2SDigits(AliRawReader* reader);
  /** @}*/

  /** @{ */
  /** 
   * @name Utility 
   */
  /** 
   * Browse this object 
   *
   * @param b Browser to show this object in 
   */
  void   Browse(TBrowser* b);
  /** @}*/
protected:
  /** Initialize hit array if not already done, and return pointert. 
      @return Hit array */
  TClonesArray*      HitsArray();
  /** Initialize digit array if not already done, and return pointert. 
      @return Digit array */
  TClonesArray*      DigitsArray();
  /** Initialize summable digit array if not already done, and return
      pointert.  
      @return Summable digit array */
  TClonesArray*      SDigitsArray();

  TClonesArray*      fSDigits;              // Summable digits
  Int_t              fNsdigits;             // Number of digits  
  Bool_t             fDetailed;             // Use detailed geometry
  Bool_t             fUseOld;               // Use old approx geometry
  Bool_t             fUseAssembly;          // Use divided volumes
  
  enum {
    kSiId,                 // ID index of Si medium
    kAirId,                // ID index of Air medium
    kPlasticId,            // ID index of Plastic medium
    kPcbId,                // ID index of PCB medium
    kSiChipId,             // ID index of Si Chip medium
    kAlId,                 // ID index of Al medium
    kCarbonId,             // ID index of Carbon medium
    kCopperId,             // ID index of Copper Medium
    kKaptonId,             // ID index of Kapton Medium
    kSteelId               // ID index of Steel medium
  };  

  TObjArray*         fBad;                  //! debugging - bad hits 

private:  
  /** Copy constructor 
      @param other Object to copy from */
  AliFMD(const AliFMD& other);
  /** Assignment operator 
      @param other Object to assign from
      @return Reference to this object  */
  AliFMD& operator=(const AliFMD& other);

  ClassDef(AliFMD,11)     // Base class FMD entry point
};

#endif
//____________________________________________________________________
//
// Local Variables:
//   mode: C++
// End:
//
// EOF
//
 AliFMD.h:1
 AliFMD.h:2
 AliFMD.h:3
 AliFMD.h:4
 AliFMD.h:5
 AliFMD.h:6
 AliFMD.h:7
 AliFMD.h:8
 AliFMD.h:9
 AliFMD.h:10
 AliFMD.h:11
 AliFMD.h:12
 AliFMD.h:13
 AliFMD.h:14
 AliFMD.h:15
 AliFMD.h:16
 AliFMD.h:17
 AliFMD.h:18
 AliFMD.h:19
 AliFMD.h:20
 AliFMD.h:21
 AliFMD.h:22
 AliFMD.h:23
 AliFMD.h:24
 AliFMD.h:25
 AliFMD.h:26
 AliFMD.h:27
 AliFMD.h:28
 AliFMD.h:29
 AliFMD.h:30
 AliFMD.h:31
 AliFMD.h:32
 AliFMD.h:33
 AliFMD.h:34
 AliFMD.h:35
 AliFMD.h:36
 AliFMD.h:37
 AliFMD.h:38
 AliFMD.h:39
 AliFMD.h:40
 AliFMD.h:41
 AliFMD.h:42
 AliFMD.h:43
 AliFMD.h:44
 AliFMD.h:45
 AliFMD.h:46
 AliFMD.h:47
 AliFMD.h:48
 AliFMD.h:49
 AliFMD.h:50
 AliFMD.h:51
 AliFMD.h:52
 AliFMD.h:53
 AliFMD.h:54
 AliFMD.h:55
 AliFMD.h:56
 AliFMD.h:57
 AliFMD.h:58
 AliFMD.h:59
 AliFMD.h:60
 AliFMD.h:61
 AliFMD.h:62
 AliFMD.h:63
 AliFMD.h:64
 AliFMD.h:65
 AliFMD.h:66
 AliFMD.h:67
 AliFMD.h:68
 AliFMD.h:69
 AliFMD.h:70
 AliFMD.h:71
 AliFMD.h:72
 AliFMD.h:73
 AliFMD.h:74
 AliFMD.h:75
 AliFMD.h:76
 AliFMD.h:77
 AliFMD.h:78
 AliFMD.h:79
 AliFMD.h:80
 AliFMD.h:81
 AliFMD.h:82
 AliFMD.h:83
 AliFMD.h:84
 AliFMD.h:85
 AliFMD.h:86
 AliFMD.h:87
 AliFMD.h:88
 AliFMD.h:89
 AliFMD.h:90
 AliFMD.h:91
 AliFMD.h:92
 AliFMD.h:93
 AliFMD.h:94
 AliFMD.h:95
 AliFMD.h:96
 AliFMD.h:97
 AliFMD.h:98
 AliFMD.h:99
 AliFMD.h:100
 AliFMD.h:101
 AliFMD.h:102
 AliFMD.h:103
 AliFMD.h:104
 AliFMD.h:105
 AliFMD.h:106
 AliFMD.h:107
 AliFMD.h:108
 AliFMD.h:109
 AliFMD.h:110
 AliFMD.h:111
 AliFMD.h:112
 AliFMD.h:113
 AliFMD.h:114
 AliFMD.h:115
 AliFMD.h:116
 AliFMD.h:117
 AliFMD.h:118
 AliFMD.h:119
 AliFMD.h:120
 AliFMD.h:121
 AliFMD.h:122
 AliFMD.h:123
 AliFMD.h:124
 AliFMD.h:125
 AliFMD.h:126
 AliFMD.h:127
 AliFMD.h:128
 AliFMD.h:129
 AliFMD.h:130
 AliFMD.h:131
 AliFMD.h:132
 AliFMD.h:133
 AliFMD.h:134
 AliFMD.h:135
 AliFMD.h:136
 AliFMD.h:137
 AliFMD.h:138
 AliFMD.h:139
 AliFMD.h:140
 AliFMD.h:141
 AliFMD.h:142
 AliFMD.h:143
 AliFMD.h:144
 AliFMD.h:145
 AliFMD.h:146
 AliFMD.h:147
 AliFMD.h:148
 AliFMD.h:149
 AliFMD.h:150
 AliFMD.h:151
 AliFMD.h:152
 AliFMD.h:153
 AliFMD.h:154
 AliFMD.h:155
 AliFMD.h:156
 AliFMD.h:157
 AliFMD.h:158
 AliFMD.h:159
 AliFMD.h:160
 AliFMD.h:161
 AliFMD.h:162
 AliFMD.h:163
 AliFMD.h:164
 AliFMD.h:165
 AliFMD.h:166
 AliFMD.h:167
 AliFMD.h:168
 AliFMD.h:169
 AliFMD.h:170
 AliFMD.h:171
 AliFMD.h:172
 AliFMD.h:173
 AliFMD.h:174
 AliFMD.h:175
 AliFMD.h:176
 AliFMD.h:177
 AliFMD.h:178
 AliFMD.h:179
 AliFMD.h:180
 AliFMD.h:181
 AliFMD.h:182
 AliFMD.h:183
 AliFMD.h:184
 AliFMD.h:185
 AliFMD.h:186
 AliFMD.h:187
 AliFMD.h:188
 AliFMD.h:189
 AliFMD.h:190
 AliFMD.h:191
 AliFMD.h:192
 AliFMD.h:193
 AliFMD.h:194
 AliFMD.h:195
 AliFMD.h:196
 AliFMD.h:197
 AliFMD.h:198
 AliFMD.h:199
 AliFMD.h:200
 AliFMD.h:201
 AliFMD.h:202
 AliFMD.h:203
 AliFMD.h:204
 AliFMD.h:205
 AliFMD.h:206
 AliFMD.h:207
 AliFMD.h:208
 AliFMD.h:209
 AliFMD.h:210
 AliFMD.h:211
 AliFMD.h:212
 AliFMD.h:213
 AliFMD.h:214
 AliFMD.h:215
 AliFMD.h:216
 AliFMD.h:217
 AliFMD.h:218
 AliFMD.h:219
 AliFMD.h:220
 AliFMD.h:221
 AliFMD.h:222
 AliFMD.h:223
 AliFMD.h:224
 AliFMD.h:225
 AliFMD.h:226
 AliFMD.h:227
 AliFMD.h:228
 AliFMD.h:229
 AliFMD.h:230
 AliFMD.h:231
 AliFMD.h:232
 AliFMD.h:233
 AliFMD.h:234
 AliFMD.h:235
 AliFMD.h:236
 AliFMD.h:237
 AliFMD.h:238
 AliFMD.h:239
 AliFMD.h:240
 AliFMD.h:241
 AliFMD.h:242
 AliFMD.h:243
 AliFMD.h:244
 AliFMD.h:245
 AliFMD.h:246
 AliFMD.h:247
 AliFMD.h:248
 AliFMD.h:249
 AliFMD.h:250
 AliFMD.h:251
 AliFMD.h:252
 AliFMD.h:253
 AliFMD.h:254
 AliFMD.h:255
 AliFMD.h:256
 AliFMD.h:257
 AliFMD.h:258
 AliFMD.h:259
 AliFMD.h:260
 AliFMD.h:261
 AliFMD.h:262
 AliFMD.h:263
 AliFMD.h:264
 AliFMD.h:265
 AliFMD.h:266
 AliFMD.h:267
 AliFMD.h:268
 AliFMD.h:269
 AliFMD.h:270
 AliFMD.h:271
 AliFMD.h:272
 AliFMD.h:273
 AliFMD.h:274
 AliFMD.h:275
 AliFMD.h:276
 AliFMD.h:277
 AliFMD.h:278
 AliFMD.h:279
 AliFMD.h:280
 AliFMD.h:281
 AliFMD.h:282
 AliFMD.h:283
 AliFMD.h:284
 AliFMD.h:285
 AliFMD.h:286
 AliFMD.h:287
 AliFMD.h:288
 AliFMD.h:289
 AliFMD.h:290
 AliFMD.h:291
 AliFMD.h:292
 AliFMD.h:293
 AliFMD.h:294
 AliFMD.h:295
 AliFMD.h:296
 AliFMD.h:297
 AliFMD.h:298
 AliFMD.h:299
 AliFMD.h:300
 AliFMD.h:301
 AliFMD.h:302
 AliFMD.h:303
 AliFMD.h:304
 AliFMD.h:305
 AliFMD.h:306
 AliFMD.h:307
 AliFMD.h:308
 AliFMD.h:309
 AliFMD.h:310
 AliFMD.h:311
 AliFMD.h:312
 AliFMD.h:313
 AliFMD.h:314
 AliFMD.h:315
 AliFMD.h:316
 AliFMD.h:317
 AliFMD.h:318
 AliFMD.h:319
 AliFMD.h:320
 AliFMD.h:321
 AliFMD.h:322
 AliFMD.h:323
 AliFMD.h:324
 AliFMD.h:325
 AliFMD.h:326
 AliFMD.h:327
 AliFMD.h:328
 AliFMD.h:329
 AliFMD.h:330
 AliFMD.h:331
 AliFMD.h:332
 AliFMD.h:333
 AliFMD.h:334
 AliFMD.h:335
 AliFMD.h:336
 AliFMD.h:337
 AliFMD.h:338
 AliFMD.h:339
 AliFMD.h:340
 AliFMD.h:341
 AliFMD.h:342
 AliFMD.h:343
 AliFMD.h:344
 AliFMD.h:345
 AliFMD.h:346
 AliFMD.h:347
 AliFMD.h:348
 AliFMD.h:349
 AliFMD.h:350
 AliFMD.h:351
 AliFMD.h:352
 AliFMD.h:353
 AliFMD.h:354
 AliFMD.h:355
 AliFMD.h:356
 AliFMD.h:357
 AliFMD.h:358
 AliFMD.h:359
 AliFMD.h:360
 AliFMD.h:361
 AliFMD.h:362
 AliFMD.h:363
 AliFMD.h:364
 AliFMD.h:365
 AliFMD.h:366
 AliFMD.h:367
 AliFMD.h:368
 AliFMD.h:369
 AliFMD.h:370
 AliFMD.h:371
 AliFMD.h:372
 AliFMD.h:373
 AliFMD.h:374
 AliFMD.h:375
 AliFMD.h:376
 AliFMD.h:377
 AliFMD.h:378
 AliFMD.h:379
 AliFMD.h:380
 AliFMD.h:381
 AliFMD.h:382
 AliFMD.h:383
 AliFMD.h:384
 AliFMD.h:385
 AliFMD.h:386
 AliFMD.h:387
 AliFMD.h:388
 AliFMD.h:389
 AliFMD.h:390
 AliFMD.h:391
 AliFMD.h:392
 AliFMD.h:393
 AliFMD.h:394
 AliFMD.h:395
 AliFMD.h:396
 AliFMD.h:397
 AliFMD.h:398
 AliFMD.h:399
 AliFMD.h:400
 AliFMD.h:401
 AliFMD.h:402
 AliFMD.h:403
 AliFMD.h:404
 AliFMD.h:405
 AliFMD.h:406
 AliFMD.h:407
 AliFMD.h:408
 AliFMD.h:409
 AliFMD.h:410
 AliFMD.h:411
 AliFMD.h:412
 AliFMD.h:413
 AliFMD.h:414
 AliFMD.h:415
 AliFMD.h:416
 AliFMD.h:417
 AliFMD.h:418
 AliFMD.h:419
 AliFMD.h:420
 AliFMD.h:421
 AliFMD.h:422
 AliFMD.h:423
 AliFMD.h:424
 AliFMD.h:425
 AliFMD.h:426
 AliFMD.h:427
 AliFMD.h:428
 AliFMD.h:429
 AliFMD.h:430
 AliFMD.h:431
 AliFMD.h:432
 AliFMD.h:433
 AliFMD.h:434
 AliFMD.h:435
 AliFMD.h:436
 AliFMD.h:437
 AliFMD.h:438
 AliFMD.h:439
 AliFMD.h:440
 AliFMD.h:441
 AliFMD.h:442
 AliFMD.h:443
 AliFMD.h:444
 AliFMD.h:445
 AliFMD.h:446
 AliFMD.h:447
 AliFMD.h:448
 AliFMD.h:449
 AliFMD.h:450
 AliFMD.h:451
 AliFMD.h:452
 AliFMD.h:453
 AliFMD.h:454
 AliFMD.h:455
 AliFMD.h:456
 AliFMD.h:457
 AliFMD.h:458
 AliFMD.h:459
 AliFMD.h:460
 AliFMD.h:461
 AliFMD.h:462
 AliFMD.h:463
 AliFMD.h:464
 AliFMD.h:465
 AliFMD.h:466
 AliFMD.h:467
 AliFMD.h:468
 AliFMD.h:469
 AliFMD.h:470
 AliFMD.h:471
 AliFMD.h:472
 AliFMD.h:473
 AliFMD.h:474
 AliFMD.h:475
 AliFMD.h:476
 AliFMD.h:477
 AliFMD.h:478
 AliFMD.h:479
 AliFMD.h:480
 AliFMD.h:481
 AliFMD.h:482
 AliFMD.h:483
 AliFMD.h:484
 AliFMD.h:485
 AliFMD.h:486
 AliFMD.h:487
 AliFMD.h:488
 AliFMD.h:489
 AliFMD.h:490
 AliFMD.h:491
 AliFMD.h:492
 AliFMD.h:493
 AliFMD.h:494
 AliFMD.h:495
 AliFMD.h:496
 AliFMD.h:497
 AliFMD.h:498
 AliFMD.h:499
 AliFMD.h:500
 AliFMD.h:501
 AliFMD.h:502
 AliFMD.h:503
 AliFMD.h:504
 AliFMD.h:505
 AliFMD.h:506
 AliFMD.h:507
 AliFMD.h:508
 AliFMD.h:509
 AliFMD.h:510
 AliFMD.h:511
 AliFMD.h:512
 AliFMD.h:513
 AliFMD.h:514
 AliFMD.h:515
 AliFMD.h:516
 AliFMD.h:517
 AliFMD.h:518
 AliFMD.h:519
 AliFMD.h:520
 AliFMD.h:521
 AliFMD.h:522
 AliFMD.h:523
 AliFMD.h:524
 AliFMD.h:525
 AliFMD.h:526
 AliFMD.h:527
 AliFMD.h:528
 AliFMD.h:529
 AliFMD.h:530
 AliFMD.h:531
 AliFMD.h:532
 AliFMD.h:533
 AliFMD.h:534
 AliFMD.h:535
 AliFMD.h:536
 AliFMD.h:537
 AliFMD.h:538
 AliFMD.h:539
 AliFMD.h:540
 AliFMD.h:541
 AliFMD.h:542
 AliFMD.h:543
 AliFMD.h:544
 AliFMD.h:545
 AliFMD.h:546
 AliFMD.h:547
 AliFMD.h:548
 AliFMD.h:549
 AliFMD.h:550
 AliFMD.h:551
 AliFMD.h:552
 AliFMD.h:553
 AliFMD.h:554
 AliFMD.h:555
 AliFMD.h:556
 AliFMD.h:557
 AliFMD.h:558
 AliFMD.h:559
 AliFMD.h:560
 AliFMD.h:561
 AliFMD.h:562
 AliFMD.h:563
 AliFMD.h:564
 AliFMD.h:565
 AliFMD.h:566
 AliFMD.h:567
 AliFMD.h:568
 AliFMD.h:569
 AliFMD.h:570
 AliFMD.h:571
 AliFMD.h:572
 AliFMD.h:573
 AliFMD.h:574
 AliFMD.h:575
 AliFMD.h:576
 AliFMD.h:577
 AliFMD.h:578
 AliFMD.h:579
 AliFMD.h:580
 AliFMD.h:581
 AliFMD.h:582
 AliFMD.h:583
 AliFMD.h:584
 AliFMD.h:585
 AliFMD.h:586
 AliFMD.h:587
 AliFMD.h:588
 AliFMD.h:589
 AliFMD.h:590
 AliFMD.h:591
 AliFMD.h:592
 AliFMD.h:593
 AliFMD.h:594
 AliFMD.h:595
 AliFMD.h:596
 AliFMD.h:597
 AliFMD.h:598
 AliFMD.h:599
 AliFMD.h:600
 AliFMD.h:601
 AliFMD.h:602
 AliFMD.h:603
 AliFMD.h:604
 AliFMD.h:605
 AliFMD.h:606
 AliFMD.h:607
 AliFMD.h:608
 AliFMD.h:609
 AliFMD.h:610
 AliFMD.h:611
 AliFMD.h:612
 AliFMD.h:613
 AliFMD.h:614
 AliFMD.h:615
 AliFMD.h:616
 AliFMD.h:617
 AliFMD.h:618
 AliFMD.h:619
 AliFMD.h:620
 AliFMD.h:621
 AliFMD.h:622
 AliFMD.h:623
 AliFMD.h:624
 AliFMD.h:625
 AliFMD.h:626
 AliFMD.h:627
 AliFMD.h:628
 AliFMD.h:629
 AliFMD.h:630
 AliFMD.h:631
 AliFMD.h:632
 AliFMD.h:633
 AliFMD.h:634
 AliFMD.h:635