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