]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CPV/AliCPV.h
This commit was generated by cvs2svn to compensate for changes in r275,
[u/mrichter/AliRoot.git] / CPV / AliCPV.h
1 #ifndef CPV_H
2 #define CPV_H
3 ////////////////////////////////////////////////
4 //  Manager and hits classes for set:CPV      //
5 //                                            //
6 //  Author: Yuri Kharlov, IHEP, Protvino      //
7 //  e-mail: Yuri.Kharlov@cern.ch              //
8 //  Last modified: 17 September 1999          //
9 ////////////////////////////////////////////////
10  
11 // --- ROOT system ---
12 #include <TArray.h> 
13 #include <TRandom.h> 
14 #include <TH2.h>
15 #include <TVector3.h>
16 #include <TLorentzVector.h>
17
18 // --- galice header files ---
19 #include "AliDetector.h"
20 #include "AliHit.h"
21 #include "AliRun.h"
22
23 //==============================================================================
24 //                              AliCPVExactHit
25 //==============================================================================
26
27 class AliCPVExactHit : public TObject {
28   
29 public:
30   TLorentzVector fMomentum;   // 4-momentum of the particle
31   Float_t        fXYhit[2];   // Hit's X,Y-coordinates
32   Int_t          fIpart;      // Hit's particle type
33   
34 public:
35   virtual ~AliCPVExactHit() {}
36            AliCPVExactHit() {}
37            AliCPVExactHit(TLorentzVector p, Float_t *xy, Int_t ipart);
38
39   TLorentzVector GetMomentum()  {return  fMomentum; }
40   Float_t        GetXY(Int_t i) {return  fXYhit[i]; }
41   Int_t          GetIpart()     {return  fIpart;    }
42   void           Print();
43
44   ClassDef(AliCPVExactHit,1)  // Hits object for set:CPV
45 };
46  
47 //==============================================================================
48 //                              AliCPVHit
49 //==============================================================================
50
51 class AliCPVHit : public TObject {
52   
53 public:
54   Float_t        fXYhit[2];   // Hit's X,Y-coordinates
55   
56 public:
57   virtual ~AliCPVHit() {}
58            AliCPVHit() {}
59            AliCPVHit(Float_t *xy);
60
61   Float_t  GetXY(Int_t i) {return  fXYhit[i]; }
62   void     Print();
63
64   ClassDef(AliCPVHit,1)  // Hits object for set:CPV
65 };
66  
67 //==============================================================================
68 //                              AliCPVCradle
69 //==============================================================================
70
71 class AliCPVCradle : public TObject {
72
73 public:
74
75   virtual   ~AliCPVCradle(void);
76              AliCPVCradle(void);
77              AliCPVCradle(Int_t   Geometry  ,
78                           Float_t PadZSize  ,
79                           Float_t PadPhiSize,
80                           Float_t Radius    ,
81                           Float_t Thickness ,
82                           Float_t Angle     ,
83                           Int_t   Nz        ,
84                           Int_t   Nphi     );
85
86   Float_t  GetPadZSize    (void) const {return fPadZSize;  }
87   Float_t  GetPadPhiSize  (void) const {return fPadPhiSize;}
88   Float_t  GetRadius      (void) const {return fRadius;    }
89   Float_t  GetThickness   (void) const {return fThickness; }
90   Int_t    GetNz          (void) const {return fNz;        }
91   Int_t    GetNphi        (void) const {return fNphi;      }
92   Float_t  GetAngle       (void) const {return fAngle;     }
93
94   void     AddHit(TLorentzVector p, Float_t *xy, Int_t ipart);
95   void     Clear(Option_t *opt="");
96   void     Print(Option_t *opt="");
97   
98   void     Reconstruction(Float_t min_distance, Float_t min_signal);
99   
100   TClonesArray *GetHitExact         (void) {return fHitExact;        }
101   TClonesArray *GetHitReconstructed (void) {return fHitReconstructed;}
102
103 private:
104
105   Float_t       fPadZSize;          // Pad size in beam direction in cm
106   Float_t       fPadPhiSize;        // Pad size in phi direction in cm
107   Float_t       fRadius;            // Distance from IP to CPV in cm
108   Float_t       fThickness;         // CPV thickness in cm
109   Float_t       fAngle;             // Position of CRADLE center in degrees
110   Int_t         fNz;                // Cells amount in beam direction
111   Int_t         fNphi;              // Cells amount in phi direction
112   
113   TClonesArray *fHitExact;          // List of exact hits in the cradle
114   TClonesArray *fHitReconstructed;  // List of reconstructed hits inthe cradle
115   Int_t         fNhitsExact;        // Number of exact hits
116   Int_t         fNhitsReconstructed;// Number of reconstructed hits
117   
118   ClassDef(AliCPVCradle,1)          // CPV cradle
119 };
120
121 //==============================================================================
122 //                            AliCPV
123 //==============================================================================
124
125 class AliCPV : public AliDetector {
126
127 public:
128
129                 AliCPV();
130                 AliCPV(const char *name, const char *title);
131   virtual      ~AliCPV();
132
133   virtual void  Init           ();
134   virtual void  SetGeometry    (Int_t ncradles, Int_t nz, Int_t nphi, Float_t angle);
135   virtual void  BuildGeometry  ();
136   virtual void  CreateGeometry () {}
137   virtual void  CreateMaterials();
138   void          FinishEvent    (void);
139   void          FinishRun      (void);
140
141   void          ResetDigits    (void);
142   void          Print          (Option_t *opt="");
143
144   AliCPVCradle &GetCradle(int n) {return *(AliCPVCradle*)fCradles->operator[](n);}
145   void          Reconstruction(Float_t min_distance, Float_t min_signal);
146
147   virtual void  StepManager()=0;
148   virtual void  AddCPVCradles();
149
150   virtual Int_t GetCPVIdtmed (void){Int_t *idtmed = fIdtmed->GetArray()-1999;
151                                     return idtmed[2000];}
152
153   Float_t  GetPadZSize          (void) const {return fPadZSize;  }
154   Float_t  GetPadPhiSize        (void) const {return fPadPhiSize;}
155   Float_t  GetRadius            (void) const {return fRadius;    }
156   Float_t  GetThickness         (void) const {return fThickness; }
157   Int_t    GetNz                (void) const {return fNz;        }
158   Int_t    GetNphi              (void) const {return fNphi;      }
159   Int_t    GetNofCradles        (void) const {return fNCradles;  }
160   Float_t  GetAngle             (void) const {return fAngle;     }
161
162   Int_t                 fDebugLevel;
163
164 private:
165
166   Float_t               fPadZSize;           // Pad size along beam       [cm]
167   Float_t               fPadPhiSize;         // Pad size across beam      [cm]
168   Float_t               fRadius;             // Distance of CPV from IP   [cm]
169   Float_t               fThickness;          // CPV thickness             [cm]
170   Float_t               fAngle;              // Angle between CPV cradles [deg]
171   Int_t                 fNz;                 // Number of pads along beam
172   Int_t                 fNphi;               // Number of pads across beam
173
174   Int_t                 fNCradles;           // Number of cradles
175   TClonesArray         *fCradles;            // Array of CPV cradles
176
177   ClassDef(AliCPV,1)  // Detector CPV
178 };
179  
180 #endif