]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/ESD/AliESDMuonTrack.h
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / STEER / ESD / AliESDMuonTrack.h
1 #ifndef ALIESDMUONTRACK_H
2 #define ALIESDMUONTRACK_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 /// \class AliESDMuonTrack
10 /// \brief Class to describe the MUON tracks in the Event Summary Data class
11 //  Author: G.Martinez
12
13
14 #include <TMath.h>
15 #include <TMatrixD.h>
16 #include <TDatabasePDG.h>
17 #include <TArrayI.h>
18
19 #include "AliVParticle.h"
20
21 class AliESDEvent;
22 class TClonesArray;
23 class TLorentzVector;
24
25 class AliESDMuonTrack : public AliVParticle {
26 public:
27   AliESDMuonTrack(); //Constructor
28   virtual ~AliESDMuonTrack(); // Destructor
29   AliESDMuonTrack(const AliESDMuonTrack& esdm);
30   AliESDMuonTrack& operator=(const AliESDMuonTrack& esdm);
31   virtual void Copy(TObject &obj) const;
32
33   virtual void Clear(Option_t* opt = "");
34   
35   void Reset();
36   
37   // Return kTRUE if the track contain tracker data
38   Bool_t ContainTrackerData() const {return (fMuonClusterMap>0) ? kTRUE : kFALSE;}
39   // Return kTRUE if the track contain trigger data
40   Bool_t ContainTriggerData() const {return (LoCircuit()>0) ? kTRUE : kFALSE;}
41   
42   // Get and Set methods for data at vertex
43   Double_t GetInverseBendingMomentum(void) const {return fInverseBendingMomentum;}
44   void     SetInverseBendingMomentum(Double_t InverseBendingMomentum) 
45                 {fInverseBendingMomentum = InverseBendingMomentum;}
46   Double_t GetThetaX(void) const {return fThetaX;}
47   void     SetThetaX(Double_t ThetaX) {fThetaX = ThetaX;}
48   Double_t GetThetaY(void) const {return fThetaY;}
49   void     SetThetaY(Double_t ThetaY) {fThetaY = ThetaY;}
50   Double_t GetZ(void) const {return fZ;}
51   void     SetZ(Double_t Z) {fZ = Z;}
52   Double_t GetBendingCoor(void) const {return fBendingCoor;}
53   void     SetBendingCoor(Double_t BendingCoor) {fBendingCoor = BendingCoor;}
54   Double_t GetNonBendingCoor(void) const {return fNonBendingCoor;}
55   void     SetNonBendingCoor(Double_t NonBendingCoor) {fNonBendingCoor = NonBendingCoor;}
56   
57   // Get and Set methods for data at Distance of Closest Approach in the vertex plane
58   Double_t GetInverseBendingMomentumAtDCA(void) const {return fInverseBendingMomentumAtDCA;}
59   void     SetInverseBendingMomentumAtDCA(Double_t InverseBendingMomentum) 
60                 {fInverseBendingMomentumAtDCA = InverseBendingMomentum;}
61   Double_t GetThetaXAtDCA(void) const {return fThetaXAtDCA;}
62   void     SetThetaXAtDCA(Double_t ThetaX) {fThetaXAtDCA = ThetaX;}
63   Double_t GetThetaYAtDCA(void) const {return fThetaYAtDCA;}
64   void     SetThetaYAtDCA(Double_t ThetaY) {fThetaYAtDCA = ThetaY;}
65   Double_t GetBendingCoorAtDCA(void) const {return fBendingCoorAtDCA;}
66   void     SetBendingCoorAtDCA(Double_t BendingCoor) {fBendingCoorAtDCA = BendingCoor;}
67   Double_t GetNonBendingCoorAtDCA(void) const {return fNonBendingCoorAtDCA;}
68   void     SetNonBendingCoorAtDCA(Double_t NonBendingCoor) {fNonBendingCoorAtDCA = NonBendingCoor;}
69   Double_t GetDCA(void) const {return TMath::Sqrt(fNonBendingCoorAtDCA*fNonBendingCoorAtDCA +
70                                                   fBendingCoorAtDCA*fBendingCoorAtDCA);}
71   
72   // Get and Set methods for data at first station
73   Double_t GetInverseBendingMomentumUncorrected(void) const {return fInverseBendingMomentumUncorrected;}
74   void     SetInverseBendingMomentumUncorrected(Double_t InverseBendingMomentum) 
75                 {fInverseBendingMomentumUncorrected = InverseBendingMomentum;}
76   Double_t GetThetaXUncorrected(void) const {return fThetaXUncorrected;}
77   void     SetThetaXUncorrected(Double_t ThetaX) {fThetaXUncorrected = ThetaX;}
78   Double_t GetThetaYUncorrected(void) const {return fThetaYUncorrected;}
79   void     SetThetaYUncorrected(Double_t ThetaY) {fThetaYUncorrected = ThetaY;}
80   Double_t GetZUncorrected(void) const {return fZUncorrected;}
81   void     SetZUncorrected(Double_t Z) {fZUncorrected = Z;}
82   Double_t GetBendingCoorUncorrected(void) const {return fBendingCoorUncorrected;}
83   void     SetBendingCoorUncorrected(Double_t BendingCoor) {fBendingCoorUncorrected = BendingCoor;}
84   Double_t GetNonBendingCoorUncorrected(void) const {return fNonBendingCoorUncorrected;}
85   void     SetNonBendingCoorUncorrected(Double_t NonBendingCoor) {fNonBendingCoorUncorrected = NonBendingCoor;}
86   
87   // Get and Set methods for covariance matrix of data at first station
88   void     GetCovariances(TMatrixD& cov) const;
89   void     SetCovariances(const TMatrixD& cov);
90   void     GetCovarianceXYZPxPyPz(Double_t cov[21]) const;
91   
92   // Get and Set methods for the transverse position r of the track at the end of the absorber
93   Double_t GetRAtAbsorberEnd() const { return fRAtAbsorberEnd; }
94   void     SetRAtAbsorberEnd(Double_t r) { fRAtAbsorberEnd = r; }
95   
96   // Get and Set methods for global tracking info
97   Double_t GetChi2(void) const {return fChi2;}
98   void     SetChi2(Double_t Chi2) {fChi2 = Chi2;}
99   UChar_t  GetNHit(void) const {return fNHit;}
100   Int_t    GetNDF() const;
101   Double_t GetNormalizedChi2() const;
102   
103   // Get and Set methods for trigger matching
104   Int_t    GetMatchTrigger() const;
105   Bool_t   MatchTriggerDigits(Bool_t fromTrack) const;
106   Double_t GetChi2MatchTrigger() const {return fChi2MatchTrigger;}
107   void     SetChi2MatchTrigger(Double_t Chi2MatchTrigger) {fChi2MatchTrigger = Chi2MatchTrigger;}
108   UShort_t GetHitsPatternInTrigCh() const {return fHitsPatternInTrigCh;}
109   void     SetHitsPatternInTrigCh(UShort_t hitsPatternInTrigCh) {fHitsPatternInTrigCh = hitsPatternInTrigCh;}
110   UInt_t   GetHitsPatternInTrigChTrk() const {return fHitsPatternInTrigChTrk;}
111   void     SetHitsPatternInTrigChTrk(UInt_t hitsPatternInTrigChTrk) {fHitsPatternInTrigChTrk = hitsPatternInTrigChTrk;}
112   void     SetLocalTrigger(Int_t locTrig) { fLocalTrigger = locTrig; }
113   Int_t    LoCircuit(void) const { return fLocalTrigger & 0xFF;       }
114   Int_t    LoStripX(void) const  { return fLocalTrigger >>  8 & 0x1F; }
115   Int_t    LoStripY(void) const  { return fLocalTrigger >> 13 & 0x0F; }
116   Int_t    LoDev(void)    const  { return fLocalTrigger >> 17 & 0x1F; }
117   Int_t    LoLpt(void)    const  { return fLocalTrigger >> 22 & 0x03; }
118   Int_t    LoHpt(void)    const  { return fLocalTrigger >> 24 & 0x03; }
119   Int_t    GetTriggerWithoutChamber(void) const { return fLocalTrigger >> 26 & 0xF; }
120   Bool_t TriggerFiredWithoutChamber(Int_t ich) const { return GetTriggerWithoutChamber() >> (3 - ich) & 0x1; }
121
122   // Get and Set methods for the hit strips pattern in the trigger chambers
123   UShort_t GetTriggerX1Pattern() const { return fX1Pattern; }
124   UShort_t GetTriggerY1Pattern() const { return fY1Pattern; }
125   UShort_t GetTriggerX2Pattern() const { return fX2Pattern; }
126   UShort_t GetTriggerY2Pattern() const { return fY2Pattern; }
127   UShort_t GetTriggerX3Pattern() const { return fX3Pattern; }
128   UShort_t GetTriggerY3Pattern() const { return fY3Pattern; }
129   UShort_t GetTriggerX4Pattern() const { return fX4Pattern; }
130   UShort_t GetTriggerY4Pattern() const { return fY4Pattern; }
131   void     SetTriggerX1Pattern(UShort_t pat) { fX1Pattern = pat; }
132   void     SetTriggerY1Pattern(UShort_t pat) { fY1Pattern = pat; }
133   void     SetTriggerX2Pattern(UShort_t pat) { fX2Pattern = pat; }
134   void     SetTriggerY2Pattern(UShort_t pat) { fY2Pattern = pat; }
135   void     SetTriggerX3Pattern(UShort_t pat) { fX3Pattern = pat; }
136   void     SetTriggerY3Pattern(UShort_t pat) { fY3Pattern = pat; }
137   void     SetTriggerX4Pattern(UShort_t pat) { fX4Pattern = pat; }
138   void     SetTriggerY4Pattern(UShort_t pat) { fY4Pattern = pat; }
139
140   // Get and Set methods for muon cluster map
141   UInt_t   GetMuonClusterMap() const {return fMuonClusterMap;}
142   void     SetMuonClusterMap(UInt_t muonClusterMap) {fMuonClusterMap = muonClusterMap;}
143   void     AddInMuonClusterMap(Int_t chamber) {fMuonClusterMap |= BIT(chamber);}
144   Bool_t   IsInMuonClusterMap(Int_t chamber) const {return (Bool_t) ((fMuonClusterMap & BIT(chamber)) != 0);}
145   
146   // Identify the tracks sharing cluster(s) with another (use the last bit of fMuonClusterMap)
147   void     Connected(Bool_t flag = kTRUE) {flag ? SETBIT(fMuonClusterMap,31) : CLRBIT(fMuonClusterMap,31);}
148   Bool_t   IsConnected() const {return TESTBIT(fMuonClusterMap,31);}
149   
150   // Methods to fill and get the Id of associated clusters
151   void     AddClusterId(UInt_t clusterId);
152   Int_t    GetNClusters() const {return static_cast<Int_t>(fNHit);}
153   UInt_t   GetClusterId(Int_t i) const {return (fClustersId && i >= 0 && i < fNHit) ? static_cast<UInt_t>(fClustersId->At(i)) : 0;}
154   
155   // Method to transfer clusters to the new ESD structure
156   Bool_t   IsOldTrack() {return (fClusters);}
157   void     MoveClustersToESD(AliESDEvent &esd);
158   
159   // Methods to compute track momentum
160   Double_t Px() const;
161   Double_t Py() const;
162   Double_t Pz() const;
163   Double_t P() const;
164   Bool_t   PxPyPz(Double_t p[3]) const { p[0] = Px(); p[1] = Py(); p[2] = Pz(); return kTRUE; }
165   void     LorentzP(TLorentzVector& vP) const;
166   Double_t PxAtDCA() const;
167   Double_t PyAtDCA() const;
168   Double_t PzAtDCA() const;
169   Double_t PAtDCA() const;
170   Bool_t   PxPyPzAtDCA(Double_t p[3]) const { p[0] = Px(); p[1] = Py(); p[2] = Pz(); return kTRUE; }
171   void     LorentzPAtDCA(TLorentzVector& vP) const;
172   Double_t PxUncorrected() const;
173   Double_t PyUncorrected() const;
174   Double_t PzUncorrected() const;
175   Double_t PUncorrected() const;
176   Bool_t   PxPyPzUncorrected(Double_t p[3]) const { p[0] = Px(); p[1] = Py(); p[2] = Pz(); return kTRUE; }
177   void     LorentzPUncorrected(TLorentzVector& vP) const;
178   
179   // additional methods to comply with AliVParticle
180   Double_t Xv() const {return -999.;} // put reasonable values here
181   Double_t Yv() const {return -999.;} //
182   Double_t Zv() const {return -999.;} //
183   Bool_t   XvYvZv(Double_t x[3]) const { x[0] = Xv(); x[1] = Yv(); x[2] = Zv(); return kTRUE; }  
184   Double_t Pt() const { return TMath::Sqrt(Px()*Px() + Py()*Py()); }
185   Double_t OneOverPt() const { return (Pt() != 0.) ? 1./Pt() : FLT_MAX; }
186   Double_t Phi() const { return TMath::Pi()+TMath::ATan2(-Py(), -Px()); }
187   Double_t Theta() const { return TMath::ATan2(Pt(), Pz()); }
188   Double_t E() const { return TMath::Sqrt(M()*M() + P()*P()); }
189   Double_t M() const { return TDatabasePDG::Instance()->GetParticle("mu-")->Mass(); }
190   Double_t Eta() const { return -TMath::Log(TMath::Tan(0.5 * Theta()));}
191   Double_t Y() const { return (Pz()/E() != 1.) ? TMath::ATanH(Pz()/E()) : FLT_MAX; }
192   Short_t  Charge() const { return (Short_t)TMath::Sign(1., GetInverseBendingMomentum()); }
193
194
195   // Dummy
196   const Double_t *PID() const { return (Double_t*)0x0; }
197   Int_t    PdgCode() const {return 0;}
198   
199   /// Set the corresponding MC track number
200   void  SetLabel(Int_t label) {fLabel = label;}
201   /// Return the corresponding MC track number
202   Int_t GetLabel() const {return fLabel;}
203
204   /// Additional methods to decode hit pattern
205   /// The hit pattern is a UShort_t with:
206   /// <pre>
207   ///   0    |   1 0 0 0 1  |  1 1  |   1  1  0  1  |   1  1  0  1   
208   ///        |              |       |               | 
209   /// unused |  RPC (0-17)  |  flag |  Bend plane   | Non-bend plane
210   ///        |      or      |       | Match chamber | Match chamber
211   ///        | further info |       |  11 12 13 14  |  11 12 13 14
212   ///        |    (20-24)   |       |               |
213   /// </pre>
214   enum EAliTriggerChPatternFlag {
215     kNoEff,    ///< Track is not good for chamber efficiency evaluation
216     kChEff,    ///< Track crosses different RPCs
217     kSlatEff,  ///< Track crosses the same RPC in all planes
218     kBoardEff  ///< Track crosses the same board in all planes
219   };
220   enum EAliTriggerChPatternInfo {
221     kCrossDifferentSlats  = 20, ///< The RPC cannot be univoquely determined
222     kTrackMatchesManyPads = 21, ///< Track not good for effciency calculation since it matches many pads
223     kTrackMatchesFewPads  = 22, ///< Track not good for effciency calculation since it matches pads in less than 3/4 chambers
224     kTrackOutsideGeometry = 23, ///< Problems in pattern determination since track extrapolation is outside trigger chambers
225     kTrackerTrackPattern  = 24, ///< The pattern was calculated from a tracker track not matching trigger track
226     kTrackMatchesMasks    = 25  ///< Track not good for effciency calculation since it matches masked pads
227   };
228   /// Set hits pattern
229   static void SetFiredChamber(UInt_t& pattern, Int_t cathode, Int_t chamber);
230   /// Add efficiency flag and crossed RPC or info on rejected track
231   static void AddEffInfo(UInt_t& pattern, Int_t slatOrInfo, Int_t board = 0, EAliTriggerChPatternFlag effType = kNoEff);
232   /// Chamber was hit
233   static Bool_t IsChamberHit(UInt_t pattern, Int_t cathode, Int_t chamber);
234   /// Get Efficiency flag
235   static Int_t GetEffFlag(UInt_t pattern);
236   /// Getting crossed slat or info
237   static Int_t GetSlatOrInfo(UInt_t pattern);
238   /// Getting crossed board
239   static Int_t GetCrossedBoard(UInt_t pattern);
240
241   AliESDEvent* GetESDEvent() const {return fESDEvent;}
242   void         SetESDEvent(AliESDEvent* evt) {fESDEvent = evt;}  
243   
244 protected:
245   // parameters at vertex
246   Double32_t fInverseBendingMomentum; ///< Inverse bending momentum (GeV/c ** -1) times the charge 
247   Double32_t fThetaX;                 ///< Angle of track at vertex in X direction (rad)
248   Double32_t fThetaY;                 ///< Angle of track at vertex in Y direction (rad)
249   Double32_t fZ;                      ///< Z coordinate (cm)
250   Double32_t fBendingCoor;            ///< bending coordinate (cm)
251   Double32_t fNonBendingCoor;         ///< non bending coordinate (cm)
252   
253   // parameters at Distance of Closest Approach in the vertex plane
254   Double32_t fInverseBendingMomentumAtDCA; ///< Inverse bending momentum (GeV/c ** -1) times the charge 
255   Double32_t fThetaXAtDCA;                 ///< Angle of track at vertex in X direction (rad)
256   Double32_t fThetaYAtDCA;                 ///< Angle of track at vertex in Y direction (rad)
257   Double32_t fBendingCoorAtDCA;            ///< bending coordinate (cm)
258   Double32_t fNonBendingCoorAtDCA;         ///< non bending coordinate (cm)
259   
260   // parameters at first tracking station
261   Double32_t fInverseBendingMomentumUncorrected; ///< Inverse bending momentum (GeV/c ** -1) times the charge 
262   Double32_t fThetaXUncorrected;                 ///< Angle of track at vertex in X direction (rad)
263   Double32_t fThetaYUncorrected;                 ///< Angle of track at vertex in Y direction (rad)
264   Double32_t fZUncorrected;                      ///< Z coordinate (cm)
265   Double32_t fBendingCoorUncorrected;            ///< bending coordinate (cm)
266   Double32_t fNonBendingCoorUncorrected;         ///< non bending coordinate (cm)
267   
268   /// reduced covariance matrix of UNCORRECTED track parameters, ordered as follow:      <pre>
269   /// [0] =  <X,X>
270   /// [1] =<X,ThetaX>  [2] =<ThetaX,ThetaX>
271   /// [3] =  <X,Y>     [4] =  <Y,ThetaX>     [5] =  <Y,Y>
272   /// [6] =<X,ThetaY>  [7] =<ThetaX,ThetaY>  [8] =<Y,ThetaY>  [9] =<ThetaY,ThetaY>
273   /// [10]=<X,InvP_yz> [11]=<ThetaX,InvP_yz> [12]=<Y,InvP_yz> [13]=<ThetaY,InvP_yz> [14]=<InvP_yz,InvP_yz>  </pre>
274   Double32_t fCovariances[15]; ///< \brief reduced covariance matrix of parameters AT FIRST CHAMBER
275   
276   Double32_t fRAtAbsorberEnd; ///< transverse position r of the track at the end of the absorber
277   
278   // global tracking info
279   Double32_t fChi2;                ///< chi2 in the MUON track fit
280   Double32_t fChi2MatchTrigger;    ///< chi2 of trigger/track matching
281   Int_t      fLocalTrigger;        ///< packed local trigger information
282
283   // hit strips pattern in the trigger chambers
284   UShort_t fX1Pattern;             ///< x-strips pattern in st6/ch1
285   UShort_t fY1Pattern;             ///< y-strips pattern in st6/ch1
286   UShort_t fX2Pattern;             ///< x-strips pattern in st6/ch2
287   UShort_t fY2Pattern;             ///< y-strips pattern in st6/ch2
288   UShort_t fX3Pattern;             ///< x-strips pattern in st7/ch1
289   UShort_t fY3Pattern;             ///< y-strips pattern in st7/ch1
290   UShort_t fX4Pattern;             ///< x-strips pattern in st7/ch2
291   UShort_t fY4Pattern;             ///< y-strips pattern in st7/ch2
292   
293   UInt_t   fMuonClusterMap;        ///< Map of clusters in tracking chambers
294   UShort_t fHitsPatternInTrigCh;   ///< Word containing info on the hits left in trigger chambers
295   UInt_t  fHitsPatternInTrigChTrk; ///< Trigger hit map from tracker track extrapolation
296   UChar_t  fNHit;                  ///< number of clusters attached to the track
297   
298   mutable TClonesArray* fClusters; ///< Array of clusters attached to the track -- deprecated
299   
300   TArrayI* fClustersId;            ///< Array of clusters'Id attached to the track
301   
302   Int_t fLabel;                    ///< point to the corresponding MC track
303
304   AliESDEvent*   fESDEvent; //!Pointer back to event to which the track belongs
305   
306   ClassDef(AliESDMuonTrack,15) // MUON ESD track class 
307 };
308
309 #endif