]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICH.h
Including RVersion.h
[u/mrichter/AliRoot.git] / RICH / AliRICH.h
1 #ifndef AliRICH_h
2 #define AliRICH_h
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 #include <TObjArray.h>
8 #include <TClonesArray.h>
9 #include <TLorentzVector.h>
10 #include <AliDetector.h>
11 #include <AliHit.h>
12 #include <AliDigit.h>
13 #include "AliRICHConst.h"
14 #include "AliRICHChamber.h"
15
16 class AliRICHRawCluster;
17 class AliRICHRecHit1D;
18 class AliRICHRecHit3D;
19
20 //__________________AliRICHhit______________________________________________________________________
21 //__________________________________________________________________________________________________
22 //__________________________________________________________________________________________________
23 class AliRICHhit : public AliHit
24 {
25 public:
26   inline   AliRICHhit();
27   inline   AliRICHhit(Int_t fIshunt, Int_t track, Int_t *vol, Float_t *hits);
28   inline   AliRICHhit(Int_t track,Int_t iPID,Int_t iChamber,TLorentzVector x4,Float_t eloss);
29   virtual ~AliRICHhit()         {;}
30     
31   Int_t   Chamber()             {return fChamber;}
32   Int_t   Particle()            {return fParticle;}    
33   Float_t Theta()               {return fTheta;}
34   Float_t Phi()                 {return fPhi;}
35   Float_t Tlength()             {return fTlength;}
36   Float_t Eloss()               {return fEloss;}
37   Float_t Loss()                {return fLoss;}
38   Float_t PHfirst()             {return fPHfirst;}
39   Float_t PHlast()              {return fPHlast;}
40   Float_t MomX()                {return fMomX;}
41   Float_t MomY()                {return fMomY;}
42   Float_t MomZ()                {return fMomZ;}
43   Float_t CerenkovAngle()       {return fMomX;}
44   Float_t MomFreoX()            {return fMomX;}
45   Float_t MomFreoY()            {return fMomY;}
46   Float_t MomFreoZ()            {return fMomZ;}
47 protected:
48   Int_t     fChamber;                      //chamber number
49   Int_t     fParticle;                     //particle code
50   Float_t   fTheta,fPhi ;                  //incident theta phi angles in degrees      
51   Float_t   fTlength;                      //track length inside the chamber
52   Float_t   fEloss;                        //ionisation energy loss in gas   
53   Float_t   fPHfirst;                      //first padhit
54   Float_t   fPHlast;                       //last padhit
55   Float_t   fLoss;                         // did it hit the freon?
56   Float_t   fMomX,fMomY,fMomZ;             //momentum at photochatode entry point
57   Float_t   fNPads;                        // Pads hit
58   Float_t   fCerenkovAngle;                // Dummy cerenkov angle
59   Float_t   fMomFreoX,fMomFreoY,fMomFreoZ; //momentum at freon entry point
60   ClassDef(AliRICHhit,1)                   //RICH hit class
61 };//class AliRICHhit
62
63   //__________________________________________________________________________________________________
64 AliRICHhit::AliRICHhit()
65            :AliHit() 
66 {//default ctor  
67   fChamber=fParticle=kBad;
68   fTheta=fPhi=fTlength=fEloss=fPHfirst=fPHlast=fLoss=kBad;
69   fMomX=fMomY=fMomZ=fNPads=fCerenkovAngle=fMomFreoX=fMomFreoY=fMomFreoZ=kBad;
70 }//AliRICHhit::default ctor
71 //__________________________________________________________________________________________________
72 AliRICHhit::AliRICHhit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hit):
73             AliHit(shunt, track)
74 {//ctor
75   fChamber=vol[0];
76   fParticle=(Int_t)hit[0];
77   fX=hit[1];fY=hit[2];fZ=hit[3];
78   fTheta=hit[4];fPhi=hit[5];
79   fTlength=hit[6];
80   fEloss=hit[7];
81   fPHfirst=(Int_t)hit[8];
82   fPHlast=(Int_t)hit[9];
83   fLoss=hit[13];
84   fMomX=hit[14];fMomY=hit[15];fMomZ=hit[16];
85   fNPads=hit[17];
86   fCerenkovAngle=hit[18];
87   fMomFreoX=hit[19];fMomFreoY=hit[20];fMomFreoZ=hit[21];
88 }//AliRICHhit::ctor
89 //__________________________________________________________________________________________________
90 AliRICHhit::AliRICHhit(Int_t track,Int_t iPID,Int_t iChamber,TLorentzVector x4,Float_t eloss):
91             AliHit(0, track)
92 {//ctor
93   fChamber=iChamber;
94   fParticle=iPID;
95   fX=x4.X();fY=x4.Y();fZ=x4.Z();
96   fEloss=eloss;
97 }//AliRICHhit::ctor
98
99 //__________________AliRICHCerenkov_________________________________________________________________
100 //__________________________________________________________________________________________________
101 //__________________________________________________________________________________________________
102 class AliRICHCerenkov: public AliHit 
103 {
104 public:
105   inline   AliRICHCerenkov();
106   inline   AliRICHCerenkov(Int_t fIshunt, Int_t track, Int_t *vol, Float_t *Cerenkovs);
107   virtual ~AliRICHCerenkov() {;}
108 public:
109   Int_t     fChamber;          //chamber number
110   Float_t   fTheta,fPhi;       //incident theta phi angles in degrees      
111   Float_t   fTlength;          //track length inside the chamber
112   Float_t   fEloss;            //ionisation energy loss in gas
113   Int_t     fPHfirst;          //first padhit
114   Int_t     fPHlast;           //last padhit
115   Int_t     fCMother;          //index of mother particle
116   Float_t   fLoss;             //nature of particle loss
117   Float_t   fIndex;            //index of photon
118   Float_t   fProduction;       //point of production
119   Float_t   fMomX,fMomY,fMomZ; //local Momentum
120   Float_t   fNPads;           // Pads hit
121   Float_t   fCerenkovAngle;   // Cerenkov Angle
122     
123   ClassDef(AliRICHCerenkov,1)  //RICH cerenkov class
124 };//class AliRICHCerenkov
125
126 //__________________________________________________________________________________________________
127 AliRICHCerenkov::AliRICHCerenkov()
128 {//ctor
129     fChamber=kBad;
130     fX=fY=fZ=fTheta=fPhi=fTlength=fEloss=kBad;
131     fPHfirst=fPHlast=fCMother=kBad;
132     fLoss=fIndex=fProduction=fMomX=fMomY=fMomZ=fNPads=fCerenkovAngle=kBad;
133 }//AliRICHCerenkov::ctor
134 //__________________________________________________________________________________________________
135 AliRICHCerenkov::AliRICHCerenkov(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits)
136                 :AliHit(shunt, track)
137 {//ctor
138     fChamber=vol[0];
139     fX=hits[1];fY=hits[2];fZ=hits[3];
140     fTheta=hits[4];fPhi=hits[5];
141     fTlength=hits[6];
142     fEloss=hits[7];
143     fPHfirst=(Int_t)hits[8];fPHlast=(Int_t)hits[9];
144     fCMother=Int_t(hits[10]);
145     fIndex = hits[11];
146     fProduction = hits[12];  
147     fLoss=hits[13];
148     fMomX=hits[14];fMomY=hits[15];fMomZ=hits[16];
149     fNPads=hits[17];
150     fCerenkovAngle=hits[18];
151 }//AliRICHCerenkov::ctor
152
153 //__________________AliRICHdigit____________________________________________________________________
154 //__________________________________________________________________________________________________
155 //__________________________________________________________________________________________________
156 class AliRICHdigit :public AliDigit
157 {
158 public:
159            AliRICHdigit() {fPadX=fPadY=fChamber=fAdc=fTracks[0]=fTracks[1]=fTracks[2]=kBad;}
160   inline   AliRICHdigit(Int_t iC,Int_t iX,Int_t iY,Int_t iAdc,Int_t iT1,Int_t iT2,Int_t iT3);
161   virtual ~AliRICHdigit() {;}  
162   Int_t C()   const{return fChamber;}
163   Int_t X()   const{return fPadX;}
164   Int_t Y()   const{return fPadY;}
165   Int_t Pad() const{return fPad;}
166   Int_t Adc() const{return fAdc;}
167 protected:
168   Int_t fChamber;  //module number 
169   Int_t fPadX;     //pad number along X
170   Int_t fPadY;     //pad number along Y
171   Int_t fPad;      //pad number 1000*X+Y
172   Int_t fAdc;      //ADC value
173   ClassDef(AliRICHdigit,1) //RICH digit class       
174 };//class AliRICHdigit
175 //__________________________________________________________________________________________________
176 AliRICHdigit::AliRICHdigit(Int_t iC,Int_t iX,Int_t iY,Int_t iAdc,Int_t iT0,Int_t iT1,Int_t iT2)
177 {
178   fChamber=iC;fPadX=iX;fPadY=iY;fPad=1000*fPadX+fPadY;fAdc=iAdc;
179   fTracks[0]=iT0;fTracks[1]=iT1;fTracks[2]=iT2;
180 }//AliRICHdigit::ctor  
181
182 //__________________AliRICH_________________________________________________________________________
183 //__________________________________________________________________________________________________
184 //__________________________________________________________________________________________________
185 class AliRICHParam;
186 class AliRICHSDigit;
187
188 class AliRICH : public AliDetector 
189 {
190 public:
191             AliRICH();                                            
192             AliRICH(const char *name, const char *title);         
193             AliRICH(const AliRICH& RICH):AliDetector(RICH) {;}   
194   virtual  ~AliRICH();                                            
195           
196   AliRICH&  operator=(const AliRICH&)                 {return *this;}
197   virtual Int_t  IsVersion()const =0;            
198   
199   inline  void    AddHit(Int_t track, Int_t *vol, Float_t *hits);//virtual
200   inline  void    AddHit(Int_t track,Int_t iPID,Int_t iChamber,TLorentzVector x4,Float_t eloss);
201   inline  void    AddCerenkov(Int_t track, Int_t *vol, Float_t *cerenkovs);
202   inline  void    AddSDigit(Int_t iC,Int_t iX,Int_t iY,Int_t iAdc,Int_t iT0,Int_t iT1=kBad,Int_t iT2=kBad);
203   inline  void    AddDigit(Int_t iC,Int_t iX,Int_t iY,Int_t iAdc,Int_t iT0,Int_t iT1=kBad,Int_t iT2=kBad);
204   
205   inline  void    ResetHits();    //virtual
206   inline  void    ResetSDigits(); 
207           void    ResetDigits();  //virtual
208   
209   virtual void    Hits2SDigits();   
210   virtual void    SDigits2Digits();
211   virtual void    Digits2Reco();   
212           
213   AliRICHChamber* C(Int_t i)                       const{return (AliRICHChamber*)fChambers->At(i-1);}
214   AliRICHParam*   Param()                          const{return fpParam;}
215   TClonesArray*   SDigits()                        const{return fSDigits;}
216   TClonesArray*   Cerenkovs()                      const{return fCerenkovs;}
217   
218           void    CreateChambers();         
219   virtual void    CreateMaterials(); //GEANT materials definition
220           Float_t AbsoCH4(Float_t x);
221           Float_t Fresnel(Float_t ene,Float_t pdoti, Bool_t pola);
222   virtual void    BuildGeometry();   //TNode ROOT variant for event display
223   virtual void    CreateGeometry();  //GEANT volumes tree for simulation  
224   
225   virtual void    StepManager()=0;
226           void    GenerateFeedbacks(Float_t eloss);
227             
228   AliRICHChamber& Chamber(Int_t id)      {return *((AliRICHChamber *) (*fChambers)[id]);}
229   
230           void  AddDigits(Int_t id, Int_t *tracks, Int_t *charges, Int_t *digits);
231           void  AddRawCluster(Int_t id, const AliRICHRawCluster& cluster);
232           void  AddRecHit1D(Int_t id, Float_t* rechit, Float_t* photons, Int_t* padsx, Int_t* padsy);
233           void  AddRecHit3D(Int_t id, Float_t* rechit, Float_t omega, Float_t theta, Float_t phi);
234           void  ResetRawClusters();
235           void  ResetRecHits1D();
236           void  ResetRecHits3D();
237   virtual void  FindClusters(Int_t nev);
238           Int_t DistancetoPrimitive(Int_t /*px*/, Int_t /*py*/)      {return 9999;}
239    
240   virtual void   MakeBranch(Option_t *opt=" ");
241   virtual void   MakeBranchInTreeD(TTree *treeD, const char *file=0);
242   virtual void   SetTreeAddress();
243               
244   
245   TObjArray     *Dchambers()                     {return fDchambers;}
246   TObjArray     *RecHits3D()                const{return fRecHits3D;}
247   TObjArray     *RecHits1D()                const{return fRecHits1D;}
248   Int_t         *Ndch()                          {return fNdch;}
249   Int_t         *Nrechits1D()                    {return fNrechits1D;} 
250   Int_t         *Nrechits3D()                    {return fNrechits3D;} 
251   TClonesArray  *DigitsAddress(Int_t id)         {return ((TClonesArray *) (*fDchambers)[id]);}
252   TClonesArray  *RecHitsAddress1D(Int_t id) const{return ((TClonesArray *) (*fRecHits1D)[id]);}
253   TClonesArray  *RecHitsAddress3D(Int_t id) const{return ((TClonesArray *) (*fRecHits3D)[id]);}
254   TClonesArray  *RawClustAddress(Int_t id)  const{return ((TClonesArray *) (*fRawClusters)[id]);}    
255   AliRICHSDigit* FirstPad(AliRICHhit *hit, TClonesArray *clusters);
256   AliRICHSDigit* NextPad(TClonesArray *clusters);
257
258  
259   virtual void Print(Option_t *option)const;
260     
261 protected:
262   AliRICHParam         *fpParam;             //main RICH parametrization     
263   TObjArray            *fChambers;           //list of RICH chambers
264   Int_t                 fNsdigits;           //Current number of sdigits
265   Int_t                 fNcerenkovs;         //Current number of cerenkovs
266   TClonesArray         *fSDigits;            //! List of sdigits
267   TObjArray            *fDchambers;          //! Array of lists of digits
268   TClonesArray         *fCerenkovs;          //! List of cerenkovs
269   Int_t                 fNdch[kNCH];         //Array of current numbers of digits
270   TObjArray            *fRawClusters;        //!List of raw clusters
271   TObjArray            *fRecHits1D;          //!List of rec. hits
272   TObjArray            *fRecHits3D;          //!List of rec. hits
273   Int_t                 fNrawch[kNCH];       //Array of current numbers of raw clusters
274   Int_t                 fNrechits1D[kNCH];   //Array of current numbers of rec hits 1D
275   Int_t                 fNrechits3D[kNCH];   //Array of current numbers of rec hits 3D 
276
277   Int_t fCkovNumber;                         // Number of Cerenkov photons
278   Int_t fFreonProd;                          // Cerenkovs produced in freon
279   Int_t fFeedbacks;                          // Number of feedback photons
280     
281   ClassDef(AliRICH,2)                        //Main RICH class 
282 };//class AliRICH
283   
284 //__________________________________________________________________________________________________
285 void AliRICH::AddHit(Int_t track, Int_t *vol, Float_t *hits)
286 {//Adds the current hit to the RICH hits list
287   TClonesArray &tmp=*fHits;
288   new(tmp[fNhits++])AliRICHhit(fIshunt,track,vol,hits);
289 }
290 //__________________________________________________________________________________________________
291 void AliRICH::AddHit(Int_t track,Int_t iPID,Int_t iChamber,TLorentzVector x4,Float_t eloss)
292 {//Adds the current hit to the RICH hits list
293   TClonesArray &tmp=*fHits;
294   new(tmp[fNhits++])AliRICHhit(track,iPID,iChamber,x4,eloss);
295 }
296 //__________________________________________________________________________________________________
297 void AliRICH::AddCerenkov(Int_t track, Int_t *vol, Float_t *cerenkovs)
298 {//Adds the current RICH cerenkov hit to the Cerenkovs list   
299   TClonesArray &tmp=*fCerenkovs;
300   new(tmp[fNcerenkovs++]) AliRICHCerenkov(fIshunt,track,vol,cerenkovs);
301 }
302 //__________________________________________________________________________________________________
303 void AliRICH::ResetHits()
304 {//Resets hits and cerenkovs
305   AliDetector::ResetHits();
306   fNcerenkovs=0;
307   if(fCerenkovs)fCerenkovs->Clear();
308 }
309 //__________________________________________________________________________________________________
310 void AliRICH::AddSDigit(Int_t iC,Int_t iX,Int_t iY,Int_t iAdc,Int_t iT0,Int_t iT1,Int_t iT2)
311 {//Adds the current Sdigit to the RICH list of Sdigits   
312   TClonesArray &tmp=*fSDigits;
313   new(tmp[fNsdigits++])AliRICHdigit(iC,iX,iY,iAdc,iT0,iT1,iT2);
314
315 //__________________________________________________________________________________________________
316 void AliRICH::ResetSDigits()
317 {//Resets sdigits
318   fNsdigits=0;
319   if(fSDigits)fSDigits->Clear();
320 }
321 //__________________________________________________________________________________________________
322 void AliRICH::AddDigit(Int_t iC,Int_t iX,Int_t iY,Int_t iAdc,Int_t iT0,Int_t iT1,Int_t iT2)
323 {//Adds the current Digit to the RICH list of digits
324
325    TClonesArray &tmp=*fDigits;
326    new(tmp[fNdigits++]) AliRICHdigit(iC,iX,iY,iAdc,iT0,iT1,iT2);
327 }
328 #endif