]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONEventReconstructor.h
New ITS detailed geometry to be used for the PPR. Only sensitive volumes added for...
[u/mrichter/AliRoot.git] / MUON / AliMUONEventReconstructor.h
1 #ifndef ALIMUONEVENTRECONSTRUCTOR_H
2 #define ALIMUONEVENTRECONSTRUCTOR_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 #include <TROOT.h>
9
10 class AliMUONHit;
11 class AliMUONHitForRec;
12 class AliMUONSegment;
13 class TClonesArray;
14 class TFile;
15 class TTree;
16 class AliMUONRecoEvent;
17
18 // Constants which should be elsewhere ????
19 const Int_t kMaxMuonTrackingChambers = 10;
20 const Int_t kMaxMuonTrackingStations = kMaxMuonTrackingChambers / 2;
21
22 class AliMUONEventReconstructor : public TObject {
23
24  public:
25   AliMUONEventReconstructor(void); // Constructor
26   virtual ~AliMUONEventReconstructor(void); // Destructor
27   AliMUONEventReconstructor (const AliMUONEventReconstructor& Reconstructor); // copy constructor
28   AliMUONEventReconstructor& operator=(const AliMUONEventReconstructor& Reconstructor); // assignment operator
29
30   // Parameters for event reconstruction: public methods
31   // Get and Set, Set to defaults
32   Double_t GetMinBendingMomentum(void) {return fMinBendingMomentum;}
33   void SetMinBendingMomentum(Double_t MinBendingMomentum) {fMinBendingMomentum = MinBendingMomentum;}
34   Double_t GetMaxSigma2Distance(void) {return fMaxSigma2Distance;}
35   void SetMaxSigma2Distance(Double_t MaxSigma2Distance) {fMaxSigma2Distance = MaxSigma2Distance;}
36   Double_t GetBendingResolution(void) {return fBendingResolution;}
37   void SetBendingResolution(Double_t BendingResolution) {fBendingResolution = BendingResolution;}
38   Double_t GetNonBendingResolution(void) {return fNonBendingResolution;}
39   void SetNonBendingResolution(Double_t NonBendingResolution) {fNonBendingResolution = NonBendingResolution;}
40   Double_t GetChamberThicknessInX0(void) {return fChamberThicknessInX0;}
41   void SetChamberThicknessInX0(Double_t ChamberThicknessInX0) {fChamberThicknessInX0 = ChamberThicknessInX0;}
42   Double_t GetSimpleBValue(void) {return fSimpleBValue;}
43   void SetSimpleBValue(Double_t SimpleBValue) {fSimpleBValue = SimpleBValue;}
44   Double_t GetSimpleBLength(void) {return fSimpleBLength;}
45   void SetSimpleBLength(Double_t SimpleBLength) {fSimpleBLength = SimpleBLength;}
46   Double_t GetSimpleBPosition(void) {return fSimpleBPosition;}
47   void SetSimpleBPosition(Double_t SimpleBPosition) {fSimpleBPosition = SimpleBPosition;}
48   Int_t GetRecGeantHits(void) {return fRecGeantHits;}
49   void SetRecGeantHits(Int_t RecGeantHits) {fRecGeantHits = RecGeantHits;}
50   Double_t GetEfficiency(void) {return fEfficiency;}
51   void SetEfficiency(Double_t Efficiency) {fEfficiency = Efficiency;}
52   Int_t GetPrintLevel(void) {return fPrintLevel;}
53   void SetPrintLevel(Int_t PrintLevel) {fPrintLevel = PrintLevel;}
54   void SetReconstructionParametersToDefaults(void);
55
56   // Parameters for GEANT background events
57   TFile* GetBkgGeantFile(void) {return fBkgGeantFile;}
58   void SetBkgGeantFile(Text_t *BkgGeantFileName); // set background file for GEANT hits
59   void NextBkgGeantEvent(void); // next event in background file for GEANT hits
60
61   // Hits for reconstruction
62   Int_t GetNHitsForRec(void) {return fNHitsForRec;} // Number
63
64   // Reconstructed tracks
65   Int_t GetNRecTracks() {return fNRecTracks;} // Number
66   void SetNRecTracks(Int_t NRecTracks) {fNRecTracks = NRecTracks;}
67   TClonesArray* GetRecTracksPtr(void) {return fRecTracksPtr;} // Array
68
69   // Hits on reconstructed tracks
70   Int_t GetNRecTrackHits() {return fNRecTrackHits;} // Number
71   void SetNRecTrackHits(Int_t NRecTrackHits) {fNRecTrackHits = NRecTrackHits;}
72   TClonesArray* GetRecTrackHitsPtr(void) {return fRecTrackHitsPtr;} // Array
73
74   // Functions
75   Double_t GetImpactParamFromBendingMomentum(Double_t BendingMomentum);
76   Double_t GetBendingMomentumFromImpactParam(Double_t ImpactParam);
77   void EventReconstruct(void);
78   void EventDump(void);  // dump reconstructed event
79   void FillEvent();      // fill and write tree of reconstructed events
80
81  protected:
82
83  private:
84
85   // Parameters for event reconstruction
86   Double_t fMinBendingMomentum; // minimum value (GeV/c) of momentum in bending plane
87   Double_t fMaxSigma2Distance; // maximum square distance in units of the variance (maximum chi2)
88   Double_t fRMin[kMaxMuonTrackingChambers]; // minimum radius (cm)
89   Double_t fRMax[kMaxMuonTrackingChambers]; // maximum radius (cm)
90   Double_t fSegmentMaxDistBending[kMaxMuonTrackingStations]; // maximum distance (cm) for segments in bending plane
91   Double_t fSegmentMaxDistNonBending[kMaxMuonTrackingStations]; // maximum distance (cm) for segments in non bending plane
92   Double_t fBendingResolution; // chamber resolution (cm) in bending plane
93   Double_t fNonBendingResolution; // chamber resolution (cm) in non bending plane
94   Double_t fChamberThicknessInX0; // chamber thickness in number of radiation lengths
95                                   // how to take it from simulation ????
96   Double_t fSimpleBValue; // simple magnetic field: value (kG)
97   Double_t fSimpleBLength; // simple magnetic field: length (cm)
98   Double_t fSimpleBPosition; // simple magnetic field: Z central position (cm)
99   Int_t fRecGeantHits; // reconstruction from raw clusters (0) or from GEANT hits (1)
100   Double_t fEfficiency; // chamber efficiency (used for GEANT hits only)
101   Int_t fPrintLevel; // print level
102
103   // Parameters for GEANT background events
104   // should be in AliMUON class ????
105   TFile *fBkgGeantFile; // pointer to file
106   TTree *fBkgGeantTK; // pointer to tree TK
107   TClonesArray *fBkgGeantParticles;   // pointer to list of particles in tree TK
108   TTree *fBkgGeantTH; // pointer to tree TH
109   TClonesArray *fBkgGeantHits;   // pointer to list of hits in tree TH
110   Int_t fBkgGeantEventNumber;   // event number
111   
112   // Hits for reconstruction (should be in AliMUON ????)
113   TClonesArray *fHitsForRecPtr; // pointer to the array of hits for reconstruction
114   Int_t fNHitsForRec; // number of hits for reconstruction
115   // Information per chamber (should be in AliMUONChamber ????)
116   Int_t fNHitsForRecPerChamber[kMaxMuonTrackingChambers]; // number of HitsForRec
117   Int_t fIndexOfFirstHitForRecPerChamber[kMaxMuonTrackingChambers]; // index (0...) of first HitForRec
118
119   // Segments inside a station
120   TClonesArray *fSegmentsPtr[kMaxMuonTrackingStations]; // array of pointers to the segments for each station
121   Int_t fNSegments[kMaxMuonTrackingStations]; // number of segments for each station
122
123   // Reconstructed tracks
124   TClonesArray *fRecTracksPtr; // pointer to array of reconstructed tracks
125   Int_t fNRecTracks; // number of reconstructed tracks
126
127   // Track hits on reconstructed tracks
128   TClonesArray *fRecTrackHitsPtr; // pointer to array of hits on reconstructed tracks
129   Int_t fNRecTrackHits; // number of hits on reconstructed tracks
130
131   // Objects needed for tree writing
132   AliMUONRecoEvent *fRecoEvent; // the reconstructed event
133   TTree            *fEventTree; // tree of reconstructed events
134   TFile            *fTreeFile;  // file where the tree is outputed
135
136   // Functions
137   void ResetHitsForRec(void);
138   void MakeEventToBeReconstructed(void);
139   void AddHitsForRecFromGEANT(TTree *TH);
140   void AddHitsForRecFromBkgGEANT(TTree *TH, TClonesArray *Hits);
141   AliMUONHitForRec* NewHitForRecFromGEANT(AliMUONHit* Hit, Int_t TrackNumber, Int_t HitNumber, Int_t Signal);
142 /*   void AddHitsForRecFromCathodeCorrelations(TTree* TC); */
143   void AddHitsForRecFromRawClusters(TTree* TR);
144   void SortHitsForRecWithIncreasingChamber();
145   void MakeSegments(void);
146   void ResetSegments(void);
147   void MakeSegmentsPerStation(Int_t Station);
148   void MakeTracks(void);
149   void ResetTrackHits(void);
150   void ResetTracks(void);
151   Int_t MakeTrackCandidatesWithTwoSegments(AliMUONSegment *BegSegment);
152   Int_t MakeTrackCandidatesWithOneSegmentAndOnePoint(AliMUONSegment *BegSegment);
153   void MakeTrackCandidates(void);
154   void FollowTracks(void);
155   void RemoveDoubleTracks(void);
156
157   ClassDef(AliMUONEventReconstructor, 0) // MUON event reconstructor in ALICE
158     };
159         
160 #endif