]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUON.h
Minor syntax for the Alpha OSF
[u/mrichter/AliRoot.git] / MUON / AliMUON.h
1 #ifndef MUON_H
2 #define MUON_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 /* $Id$ */
7
8 ////////////////////////////////////////////////
9 //  Manager and hits classes for set:MUON     //
10 ////////////////////////////////////////////////
11 #include "AliDetector.h"
12 #include "AliHit.h"
13 #include "AliMUONConst.h"
14 #include "AliDigit.h"
15 #include "AliMUONchamber.h"
16 #include "AliMUONSegRes.h"
17 #include <TVector.h>
18 #include <TObjArray.h>
19 #include <TArrayF.h>
20 #include <TFile.h>
21 #include <TTree.h>
22 typedef enum {simple, medium, big} Cluster_t;
23
24 static const int NCH=14;
25
26 class AliMUONcluster;
27 class AliMUONRawCluster;
28 class AliMUONClusterFinder;
29 class AliMUONcorrelation;
30
31
32 //----------------------------------------------
33
34  
35 class AliMUONcluster : public TObject {
36 public:
37
38    Int_t     fHitNumber;    // Hit number
39    Int_t     fCathode;      // Cathode number
40    Int_t     fQ  ;          // Total charge      
41    Int_t     fPadX  ;       // Pad number along X
42    Int_t     fPadY  ;       // Pad number along Y
43    Int_t     fQpad  ;       // Charge per pad
44    Int_t     fRSec  ;       // R -sector of pad
45  
46 public:
47    AliMUONcluster() {
48       fHitNumber=fQ=fPadX=fPadY=fQpad=fRSec=0;   
49 }
50    AliMUONcluster(Int_t *clhits);
51    virtual ~AliMUONcluster() {;}
52  
53    ClassDef(AliMUONcluster,1)  //Cluster object for set:MUON
54 };
55
56  
57 class AliMUONreccluster : public TObject {
58 public:
59
60    Int_t     fTracks[3];      //labels of overlapped tracks
61
62    Int_t       fQ  ;          // Q of cluster (in ADC counts)     
63    Float_t     fX  ;          // X of cluster
64    Float_t     fY  ;          // Y of cluster
65
66 public:
67    AliMUONreccluster() {
68        fTracks[0]=fTracks[1]=fTracks[2]=-1; 
69        fQ=0; fX=fY=0; 
70    }
71    virtual ~AliMUONreccluster() {;}
72  
73    ClassDef(AliMUONreccluster,1)  //Cluster object for set:MUON
74 };
75
76 //_____________________________________________________________________________
77
78 class AliMUONdigit : public TObject {
79  public:
80     Int_t     fPadX;        // Pad number along x
81     Int_t     fPadY ;       // Pad number along y
82     Int_t     fSignal;      // Signal amplitude
83     Int_t     fTcharges[10];   // charge per track making this digit (up to 10)
84     Int_t     fTracks[10];     // primary tracks making this digit (up to 10)
85     Int_t     fPhysics;        // physics contribution to signal 
86     Int_t     fHit;            // hit number - temporary solution
87
88
89  
90  public:
91     AliMUONdigit() {}
92     AliMUONdigit(Int_t *digits);
93     AliMUONdigit(Int_t *tracks, Int_t *charges, Int_t *digits);
94     virtual ~AliMUONdigit();
95     
96     ClassDef(AliMUONdigit,1)  //Digits for set:MUON
97 };
98 //_____________________________________________________________________________
99
100 class AliMUONlist : public AliMUONdigit {
101  public:
102     Int_t          fChamber;       // chamber number of pad
103     TObjArray     *fTrackList; 
104  public:
105     AliMUONlist() {fTrackList=0;}
106     AliMUONlist(Int_t rpad, Int_t *digits);
107     virtual ~AliMUONlist() {delete fTrackList;}
108     TObjArray  *TrackList()   {return fTrackList;}
109     ClassDef(AliMUONlist,1)  //Digits for set:MUON
110 };
111 //___________________________________________
112
113
114 //___________________________________________
115  
116 class AliMUONhit : public AliHit {
117  public:
118     Int_t     fChamber;       // Chamber number
119     Float_t   fParticle;      // Geant3 particle type
120     Float_t   fTheta ;        // Incident theta angle in degrees      
121     Float_t   fPhi   ;        // Incident phi angle in degrees
122     Float_t   fTlength;       // Track length inside the chamber
123     Float_t   fEloss;         // ionisation energy loss in gas   
124     Int_t     fPHfirst;       // first padhit
125     Int_t     fPHlast;        // last padhit
126
127 // modifs perso
128     Float_t   fPTot;          // hit momentum P
129     Float_t   fCxHit;            // Px/P
130     Float_t   fCyHit;            // Py/P
131     Float_t   fCzHit;            // Pz/P
132
133  public:
134     AliMUONhit() {}
135     AliMUONhit(Int_t fIshunt, Int_t track, Int_t *vol, Float_t *hits);
136     virtual ~AliMUONhit() {}
137     
138     ClassDef(AliMUONhit,1)  //Hits object for set:MUON
139 };
140
141 class AliMUON : public  AliDetector {
142  public:
143     AliMUON();
144     AliMUON(const char *name, const char *title);
145     virtual       ~AliMUON();
146     virtual void   AddHit(Int_t, Int_t*, Float_t*);
147     virtual void   AddCluster(Int_t*);
148     virtual void   AddDigits(Int_t, Int_t*, Int_t*, Int_t*);
149     virtual void   AddRawCluster(Int_t, const AliMUONRawCluster&);
150     virtual void   AddCathCorrel(Int_t, Int_t*, Float_t*, Float_t*);
151     virtual void   BuildGeometry();
152     virtual void   CreateGeometry() {}
153     virtual void   CreateMaterials() {}
154     virtual void   StepManager();
155     Int_t          DistancetoPrimitive(Int_t px, Int_t py);
156     virtual Int_t  IsVersion() const =0;
157 //
158     TClonesArray  *Clusters() {return fClusters;}
159     virtual  void  MakeTreeC(Option_t *option="C");
160     void           GetTreeC(Int_t);
161     virtual void   MakeBranch(Option_t *opt=" ");
162     void           SetTreeAddress();
163     virtual void   ResetHits();
164     virtual void   ResetDigits();
165     virtual void   ResetRawClusters();
166     virtual void   ResetCorrelation();
167     virtual void   FindClusters(Int_t,Int_t);
168     virtual void   Digitise(Int_t,Int_t,Option_t *opt1=" ",Option_t *opt2=" ",Text_t *name=" ");
169     virtual void   CathodeCorrelation(Int_t);
170     virtual void   SortTracks(Int_t *,Int_t *,Int_t);
171 //
172 // modifs perso
173
174     void     InitTracking(Double_t &, Double_t &, Double_t &);
175     void     Reconst(Int_t &,Int_t &,Int_t,Int_t &,Int_t&,Int_t&, Option_t *option,Text_t *filename);
176     void     FinishEvent();
177     void     CloseTracking();
178     void     SetCutPxz(Double_t p) {fSPxzCut=p;}
179     void     SetSigmaCut(Double_t p) {fSSigmaCut=p;}
180     void     SetXPrec(Double_t p) {fSXPrec=p;}
181     void     SetYPrec(Double_t p) {fSYPrec=p;}
182     Double_t GetCutPxz() {return fSPxzCut;}
183     Double_t GetSigmaCut() {return fSSigmaCut;}
184     Double_t GetXPrec() {return fSXPrec;}
185     Double_t GetYPrec() {return fSYPrec;}
186 // fin modifs perso 
187     
188 // Configuration Methods (per station id)
189 //
190 // Set Chamber Segmentation Parameters
191 // id refers to the station and isec to the cathode plane   
192     virtual void   SetPADSIZ(Int_t id, Int_t isec, Float_t p1, Float_t p2);
193
194 // Set Signal Generation Parameters
195     virtual void   SetSigmaIntegration(Int_t id, Float_t p1);
196     virtual void   SetChargeSlope(Int_t id, Float_t p1);
197     virtual void   SetChargeSpread(Int_t id, Float_t p1, Float_t p2);
198     virtual void   SetMaxAdc(Int_t id, Float_t p1);
199 // Set Segmentation and Response Model
200     virtual void   SetSegmentationModel(Int_t id, Int_t isec, AliMUONsegmentation *segmentation);
201     virtual void   SetResponseModel(Int_t id, AliMUONresponse *response);
202     virtual void   SetNsec(Int_t id, Int_t nsec);
203 // Set Reconstruction Model
204     virtual void   SetReconstructionModel(Int_t id, AliMUONClusterFinder *reconstruction);
205 // Set Stepping Parameters
206     virtual void   SetMaxStepGas(Float_t p1);
207     virtual void   SetMaxStepAlu(Float_t p1);
208     virtual void   SetMaxDestepGas(Float_t p1);
209     virtual void   SetMaxDestepAlu(Float_t p1);
210     virtual void   SetMuonAcc(Bool_t acc=0, Float_t angmin=2, Float_t angmax=9);
211 // Response Simulation
212     virtual void   MakePadHits(Float_t xhit,Float_t yhit,Float_t eloss,Int_t id);
213 // Return reference to Chamber #id
214     virtual AliMUONchamber& Chamber(Int_t id) {return *((AliMUONchamber *) (*fChambers)[id]);}
215 // Retrieve pad hits for a given Hit
216     virtual AliMUONcluster* FirstPad(AliMUONhit *, TClonesArray *);
217     virtual AliMUONcluster* NextPad(TClonesArray *);
218 // Return pointers to digits 
219     TObjArray            *Dchambers() {return fDchambers;}
220     Int_t                *Ndch() {return fNdch;}
221     virtual TClonesArray *DigitsAddress(Int_t id) {return ((TClonesArray *) (*fDchambers)[id]);}
222 // Return pointers to reconstructed clusters
223     TObjArray            *RawClusters() {return fRawClusters;}
224     Int_t                *Nrawch() {return fNrawch;}
225     virtual TClonesArray *RawClustAddress(Int_t id) {return ((TClonesArray *) (*fRawClusters)[id]);}
226
227 // modifs perso
228     AliMUONRawCluster *RawCluster(Int_t ichamber, Int_t icathod, Int_t icluster);
229     
230     
231     // Return pointers to list of correlated clusters
232     TObjArray            *CathCorrel() {return fCathCorrel;}
233     Int_t                *Ncorch() {return fNcorch;}
234     virtual TClonesArray *CathCorrelAddress(Int_t id)
235         {return ((TClonesArray *) (*fCathCorrel)[id]);}
236
237 // Return pointer to TreeC
238     TTree      *TreeC() {return fTreeC;} 
239  protected:
240     TObjArray            *fChambers;           // List of Tracking Chambers
241     Int_t                fNclusters;           // Number of clusters
242     TClonesArray         *fClusters;           // List of clusters
243     TObjArray            *fDchambers;          // List of digits
244     Int_t                *fNdch;               // Number of digits
245
246
247     TObjArray            *fRawClusters;            // List of raw clusters
248     Int_t                *fNrawch;                 // Number of raw clusters
249     TObjArray            *fCathCorrel;             // List of correlated clusters
250     Int_t                *fNcorch;                 // Number of correl clusters
251     TTree                *fTreeC;                  // Cathode correl index tree
252
253 //
254     Bool_t   fAccCut;          //Transport acceptance cut
255     Float_t  fAccMin;          //Minimum acceptance cut used during transport
256     Float_t  fAccMax;          //Minimum acceptance cut used during transport
257 //  
258
259 //  Stepping Parameters
260    Float_t fMaxStepGas;      // Maximum step size inside the chamber gas
261    Float_t fMaxStepAlu;      // Maximum step size inside the chamber aluminum
262    Float_t fMaxDestepGas;    // Maximum relative energy loss in gas
263    Float_t fMaxDestepAlu;    // Maximum relative energy loss in aluminum
264 //
265 // modifs perso
266 //  Parameters for reconstruction program
267    Double_t fSPxzCut;        // Pxz cut  (GeV/c) to begin the track finding
268    Double_t fSSigmaCut;      // Number of sig. delimiting the searching areas
269    Double_t fSXPrec;         // Chamber precision in X (cm) 
270    Double_t fSYPrec;         // Chamber precision in Y (cm)
271
272    Text_t *fFileName;
273    
274  protected:
275
276    ClassDef(AliMUON,1)  //Hits manager for set:MUON
277 };
278 //___________________________________________
279
280 class AliMUONRawCluster : public TObject {
281 public:
282
283    Int_t     fTracks[3];      //labels of overlapped tracks
284    Int_t       fQ  ;          // Q of cluster (in ADC counts)     
285    Float_t     fX  ;          // X of cluster
286    Float_t     fY  ;          // Y of cluster
287    Int_t       fPeakSignal;
288    Int_t       fIndexMap[50];  //indeces of digits
289    Int_t       fOffsetMap[50]; //Emmanuel special
290    Float_t     fContMap[50];   //Contribution from digit
291    Int_t       fPhysicsMap[50];
292    Int_t       fMultiplicity;  //cluster multiplicity
293    Int_t       fNcluster[2];
294    Int_t       fClusterType;   
295  public:
296    AliMUONRawCluster() {
297        fTracks[0]=fTracks[1]=fTracks[2]=-1; 
298        fQ=0; fX=fY=0; fMultiplicity=0;
299        for (int k=0;k<50;k++) {
300            fIndexMap[k]=-1;
301            fOffsetMap[k]=0;
302            fContMap[k]=0;
303            fPhysicsMap[k]=-1;
304        }
305        fNcluster[0]=fNcluster[1]=-1;
306    }
307    virtual ~AliMUONRawCluster() {}
308
309    Float_t GetRadius() {return TMath::Sqrt(fX*fX+fY*fY);}
310
311    Bool_t IsSortable() const {return kTRUE;}
312    Int_t  Compare(TObject *obj);
313    Int_t PhysicsContribution();
314    static Int_t BinarySearch(Float_t r, TArrayF, Int_t from, Int_t upto);
315    static void  SortMin(Int_t *,Float_t *,Float_t *,Float_t *,Float_t *,Int_t);
316  
317    ClassDef(AliMUONRawCluster,1)  //Cluster object for set:MUON
318 };
319
320 //___________________________________________
321 class AliMUONcorrelation : public TObject {
322 public:
323
324   // correlation starts from the 1-st cathode  
325   // last number in arrays corresponds to cluster on 1-st cathode
326
327    Int_t       fCorrelIndex[4];  // entry number in TreeR for the associated 
328                                  // cluster candidates on the 2-nd cathode
329    Float_t     fX[4]  ;          // X of clusters on the 2-nd cathode  
330    Float_t     fY[4]  ;          // Y of clusters
331
332 public:
333    AliMUONcorrelation() {
334        fCorrelIndex[0]=fCorrelIndex[1]=fCorrelIndex[2]=fCorrelIndex[3]=0;
335        fX[0]=fX[1]=fX[2]=fX[3]=0; fY[0]=fY[1]=fY[2]=fY[3]=0; 
336    }
337    AliMUONcorrelation(Int_t *idx, Float_t *x, Float_t *y);
338    virtual ~AliMUONcorrelation() {}
339    ClassDef(AliMUONcorrelation,1)  //Cathode correlation object for set:MUON
340 };
341
342 #endif
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357