]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICH.h
call to Init() added to deafult ctor
[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 /* $Id$ */
8
9
10 ////////////////////////////////////////////////
11 //  Manager and hits classes for set:RICH     //
12 ////////////////////////////////////////////////
13 #include "AliDetector.h"
14 #include "AliRICHConst.h"
15 #include "AliRICHChamber.h"
16 static const int kNCH=7;
17
18 class TFile;
19
20 class AliRICHHit;
21 class AliRICHSDigit;
22 class AliRICHRawCluster;
23 class AliRICHRecHit1D;
24 class AliRICHRecHit3D;
25 class AliRICHClusterFinder;
26 class AliRICHDetect;
27 class AliRICHChamber;
28 class AliRICHCerenkov;
29 class AliSegmentation;
30 class AliRICHResponse;
31 class AliRICHGeometry;
32 class AliRICHMerger;
33
34 class AliRICH : public AliDetector 
35 {
36    
37 enum EDebugBits {kDebugStart=BIT(0),kDebugParam=BIT(1),kDebugHit=BIT(2),kDebugDigit=BIT(3),kDebugReco=BIT(4)};
38    
39 public:
40 // ctor dtor staff      
41    AliRICH(); // default ctor
42    AliRICH(const char *name, const char *title); // named ctor
43    AliRICH(const AliRICH& RICH);   
44    virtual       ~AliRICH();
45    
46     virtual void   AddHit(Int_t track, Int_t *vol, Float_t *hits);
47     virtual void   AddCerenkov(Int_t track, Int_t *vol, Float_t *cerenkovs);
48     virtual void   AddSDigit(Int_t *clhits);
49     virtual void   AddDigits(Int_t id, Int_t *tracks, Int_t *charges, Int_t *digits);
50     virtual void   AddRawCluster(Int_t id, const AliRICHRawCluster& cluster);
51     virtual void   AddRecHit1D(Int_t id, Float_t* rechit, Float_t* photons, Int_t* padsx, Int_t* padsy);
52     virtual void   AddRecHit3D(Int_t id, Float_t* rechit);
53
54
55     virtual void   BuildGeometry();   // TNode ROOT variant for event display
56     virtual void   CreateGeometry();  // GEANT volumes tree fro simulation
57     virtual void   CreateMaterials(); // GEANT materials definition
58     virtual Float_t AbsoCH4(Float_t x);
59     virtual Float_t Fresnel(Float_t ene,Float_t pdoti, Bool_t pola);
60     virtual void   StepManager();
61     Int_t          DistancetoPrimitive(Int_t px, Int_t py);
62     virtual Int_t  IsVersion() const =0;
63 //
64     TClonesArray  *SDigits() {return fSDigits;}
65     TClonesArray  *Cerenkovs() {return fCerenkovs;}
66     virtual void   MakeBranch(Option_t *opt=" ", const char *file=0);
67     virtual void   MakeBranchInTreeD(TTree *treeD, const char *file=0);
68     void           SetTreeAddress();
69     virtual void   ResetHits();
70     virtual void   ResetDigits();
71     virtual void   ResetRawClusters();
72     virtual void   ResetRecHits1D();
73     virtual void   ResetRecHits3D();
74     virtual void   FindClusters(Int_t nev,Int_t lastEntry);
75 // Converters    
76    virtual void   Hits2SDigits();
77    virtual void   SDigits2Digits();
78    virtual void   SDigits2Digits(Int_t nev, Int_t flag);
79    virtual void   Digits2Reco();
80
81 // Models for chambers
82    virtual void     SetGeometryModel(Int_t iChamberN, AliRICHGeometry *pRICHGeo)
83                                                             {((AliRICHChamber*)fChambers->At(iChamberN))->GeometryModel(pRICHGeo);}
84    AliRICHGeometry* GetGeometryModel(Int_t iChamberN=0)const{return ((AliRICHChamber*)fChambers->At(iChamberN))->GetGeometryModel();}
85     
86    virtual void     SetSegmentationModel(Int_t iChamberN, AliSegmentation *pAliSeg)
87                                                                        {((AliRICHChamber*)fChambers->At(iChamberN))->SetSegmentationModel(pAliSeg);}
88    AliSegmentation* GetSegmentationModel(Int_t iChamberN=0)const{return ((AliRICHChamber*)fChambers->At(iChamberN))->GetSegmentationModel();}
89          
90    virtual void     SetResponseModel(Int_t iChamberN, AliRICHResponse *pRICHRes)
91                                                                    {((AliRICHChamber*)fChambers->At(iChamberN))->ResponseModel(pRICHRes);}
92    AliRICHResponse* GetResponseModel(Int_t iChamberN)  const{return ((AliRICHChamber*)fChambers->At(iChamberN))->GetResponseModel();}
93
94    virtual void     SetReconstructionModel(Int_t id, AliRICHClusterFinder *pRICHReco)
95                                                            {    ((AliRICHChamber*)fChambers->At(id))->SetReconstructionModel(pRICHReco);}
96 // Debug staff
97    void   SetDebugLevel(Int_t level) {fDebugLevel=level;}
98    Int_t  GetDebugLevel()       const{return fDebugLevel;}
99    
100    void   SetDebugStart()     {fDebugLevel+=kDebugStart;}        // Controls debug message at the entring point of methods
101    void ResetDebugStart()     {fDebugLevel-=kDebugStart;}        // Controls debug message at the entring point of methods
102    Bool_t  IsDebugStart()const{return fDebugLevel&kDebugStart;}  // Controls debug message at the entring point of methods
103    
104    void   SetDebugParam()     {fDebugLevel+=kDebugParam;}        // Controls debug printout for the parameters
105    void ResetDebugParam()     {fDebugLevel-=kDebugParam;}        // Controls debug printout for the parameters
106    Bool_t  IsDebugParam()const{return fDebugLevel&kDebugParam;}  // Controls debug printout for the parameters
107    
108    void   SetDebugHit()       {fDebugLevel+=kDebugHit;}          // Controls debug printout for hits
109    void ResetDebugHit()       {fDebugLevel-=kDebugHit;}          // Controls debug printout for hits
110    Bool_t  IsDebugHit()  const{return fDebugLevel&kDebugHit;}    // Controls debug printout for hits
111    
112    void   SetDebugDigit()     {fDebugLevel+=kDebugDigit;}        // Controls debug printout for digits
113    void ResetDebugDigit()     {fDebugLevel-=kDebugDigit;}        // Controls debug printout for digits
114    Bool_t  IsDebugDigit()const{return fDebugLevel&kDebugDigit;}  // Controls debug printout for digits
115    
116    void   SetDebugReco()      {fDebugLevel+=kDebugReco;}         // Controls debug printout for reco
117    void ResetDebugReco()      {fDebugLevel-=kDebugReco;}         // Controls debug printout for reco
118    Bool_t  IsDebugReco() const{return fDebugLevel&kDebugReco;}   // Controls debug printout for reco
119    
120 // Set Merger
121     virtual void   SetMerger(AliRICHMerger* thisMerger) {fMerger=thisMerger;}  
122 // Response Simulation
123     virtual Int_t  Hits2SDigits(Float_t xhit,Float_t yhit,Float_t eloss,Int_t id, ResponseType res);
124 // Return reference to Chamber #id
125     virtual        AliRICHChamber& Chamber(Int_t id) {return *((AliRICHChamber *) (*fChambers)[id]);}
126 // Retrieve pad hits for a given Hit
127     virtual        AliRICHSDigit* FirstPad(AliRICHHit *hit, TClonesArray *clusters);
128     virtual        AliRICHSDigit* NextPad(TClonesArray *clusters);
129 // Return pointers to digits 
130     TObjArray     *Dchambers() {return fDchambers;}
131     Int_t         *Ndch() {return fNdch;}
132     virtual TClonesArray *DigitsAddress(Int_t id) {return ((TClonesArray *) (*fDchambers)[id]);}
133 // Return pointers to rec. hits
134     TObjArray     *RecHits1D() {return fRecHits1D;}
135     Int_t         *Nrechits1D() {return fNrechits1D;}
136     virtual TClonesArray *RecHitsAddress1D(Int_t id) {return ((TClonesArray *) (*fRecHits1D)[id]);}
137     TObjArray     *RecHits3D() {return fRecHits3D;}
138     Int_t         *Nrechits3D() {return fNrechits3D;}
139     virtual TClonesArray *RecHitsAddress3D(Int_t id) {return ((TClonesArray *) (*fRecHits3D)[id]);}
140     
141 // Return pointers to reconstructed clusters
142     virtual TClonesArray *RawClustAddress(Int_t id) {return ((TClonesArray *) (*fRawClusters)[id]);}    
143 // Assignment operator
144     AliRICH& operator=(const AliRICH& rhs);
145
146 // Analysis routines
147     // Full events
148     virtual void DiagnosticsFE(Int_t evNumber1=0,Int_t evNumber2=0);
149     // Single events
150     virtual void DiagnosticsSE(Int_t diaglevel,Int_t evNumber1=0,Int_t evNumber2=0);
151     
152    virtual void Print(Option_t *option)const;
153     
154 protected:
155     TObjArray            *fChambers;           // List of Tracking Chambers
156     Int_t                 fNSDigits;           // Number of clusters
157     Int_t                 fNcerenkovs;         // Number of cerenkovs
158     TClonesArray         *fSDigits;            // List of clusters
159     TObjArray            *fDchambers;          // List of digits
160     TClonesArray         *fCerenkovs;          // List of cerenkovs
161     Int_t                 fNdch[kNCH];         // Number of digits
162     TObjArray            *fRawClusters;        // List of raw clusters
163     TObjArray            *fRecHits1D;          // List of rec. hits
164     TObjArray            *fRecHits3D;          // List of rec. hits
165     Int_t                 fNrawch[kNCH];       // Number of raw clusters
166     Int_t                 fNrechits1D[kNCH];   // Number of rec hits 
167     Int_t                 fNrechits3D[kNCH];   // Number of rec hits 
168     Int_t                 fDebugLevel;         // Source debugging level
169
170     Int_t fCkovNumber;                         // Number of Cerenkov photons
171     Int_t fCkovQuarz;                          // Cerenkovs crossing quartz
172     Int_t fCkovGap;                            // Cerenkovs crossing gap
173     Int_t fCkovCsi;                            // Cerenkovs crossing csi
174     Int_t fLostRfreo;                          // Cerenkovs reflected in freon
175     Int_t fLostRquar;                          // Cerenkovs reflected in quartz
176     Int_t fLostAfreo;                          // Cerenkovs absorbed in freon 
177     Int_t fLostAquarz;                         // Cerenkovs absorbed in quartz
178     Int_t fLostAmeta;                          // Cerenkovs absorbed in methane
179     Int_t fLostCsi;                            // Cerenkovs below csi quantum efficiency 
180     Int_t fLostWires;                          // Cerenkovs lost in wires
181     Int_t fFreonProd;                          // Cerenkovs produced in freon
182     Float_t fMipx;                             // x coord. of MIP
183     Float_t fMipy;                             // y coord. of MIP
184     Int_t fFeedbacks;                          // Number of feedback photons
185     Int_t fLostFresnel;                        // Cerenkovs lost by Fresnel reflection
186
187
188 // Background event for event mixing
189     Text_t *fFileName;                         //! File with background hits
190     AliRICHMerger *fMerger;                    //! pointer to merger
191     
192     ClassDef(AliRICH,1)                        //Hits manager for set:RICH
193 };
194     
195 inline void AliRICH::Print(Option_t *option)const
196 {
197    if(IsDebugParam()){
198       GetGeometryModel(0)->Print();
199       GetSegmentationModel(0)->Print();
200       GetResponseModel(0)->Print();
201    }
202 }//inline void AliRICH::Print(Option_t *option)const
203
204 #endif