ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitAna/DataTree/interface/Track.h
Revision: 1.16
Committed: Fri Aug 29 01:51:01 2008 UTC (16 years, 8 months ago) by loizides
Content type: text/plain
Branch: MAIN
Changes since 1.15: +15 -17 lines
Log Message:
Cosmetics.

File Contents

# Content
1 //--------------------------------------------------------------------------------------------------
2 // $Id: Track.h,v 1.15 2008/07/31 13:28:42 bendavid Exp $
3 //
4 // Track
5 //
6 // We store the CMSSW track parameterization
7 // Parameters associated to the 5D curvilinear covariance matrix:
8 // (qoverp, lambda, phi, dxy, dsz)
9 // defined as:
10 // qoverp = q / abs(p) = signed inverse of momentum [1/GeV]
11 // lambda = pi/2 - polar angle at the given point
12 // phi = azimuth angle at the given point
13 // dxy = -vx*sin(phi) + vy*cos(phi) [cm]
14 // dsz = vz*cos(lambda) - (vx*cos(phi)+vy*sin(phi))*sin(lambda) [cm]
15 //
16 // Format for fHits: (We do not use anything resembling reco::HitPattern from CMSSW because that
17 // data format requires 800 bits per track!)
18 // There is a one to one mapping between bits and tracker layers, where layers are enumerated
19 // seperately in the PXB, PXF, TIB, TID, TOB, TEC and r-phi and stereo modules are treated as
20 // seperate layers in those detectors which have them
21 // (TIB L1,L2, TID L1,L2, TOB L1,L2, TEC L1,L2,L5).
22 //
23 // A bit value of 1 indicates a hit in the corresponding layer, and 0 indicates no hit.
24 //
25 // Note that currently this only stores information about hits in the Tracker,
26 // but muon chamber information will likely be added as well.
27 //
28 // Bit-Layer assignments (starting from bit 0):
29 // Bit 0: PXB L1
30 // Bit 1: PXB L2
31 // Bit 2: PXB L3
32 // Bit 3: PXF L1
33 // Bit 4: PXF L2
34 // Bit 5: TIB L1 r-phi
35 // Bit 6: TIB L1 stereo
36 // Bit 7: TIB L2 r-phi
37 // Bit 8: TIB L2 stereo
38 // Bit 9: TIB L3 r-phi
39 // Bit 10: TIB L4 r-phi
40 // Bit 11: TID L1 r-phi
41 // Bit 12: TID L1 stereo
42 // Bit 13: TID L2 r-phi
43 // Bit 14: TID L2 stereo
44 // Bit 15: TID L3 r-phi
45 // Bit 16: TOB L1 r-phi
46 // Bit 17: TOB L1 stereo
47 // Bit 18: TOB L2 r-phi
48 // Bit 19: TOB L2 stereo
49 // Bit 20: TOB L3 r-phi
50 // Bit 21: TOB L4 r-phi
51 // Bit 22: TOB L5 r-phi
52 // Bit 23: TOB L6 r-phi
53 // Bit 24: TEC L1 r-phi
54 // Bit 25: TEC L1 stereo
55 // Bit 26: TEC L2 r-phi
56 // Bit 27: TEC L2 stereo
57 // Bit 28: TEC L3 r-phi
58 // Bit 29: TEC L4 r-phi
59 // Bit 30: TEC L5 r-phi
60 // Bit 31: TEC L5 stereo
61 // Bit 32: TEC L6 r-phi
62 // Bit 33: TEC L7 r-phi
63 // Bit 34: TEC L8 r-phi
64 // Bit 35: TEC L9 r-phi
65 //
66 // Authors: C.Loizides, J.Bendavid, C.Paus
67 //--------------------------------------------------------------------------------------------------
68
69 #ifndef DATATREE_TRACK_H
70 #define DATATREE_TRACK_H
71
72 #include "MitAna/DataTree/interface/DataObject.h"
73 #include "MitAna/DataTree/interface/MCParticle.h"
74 #include "MitAna/DataTree/interface/BitMask32.h"
75 #include "MitAna/DataTree/interface/BitMask64.h"
76 #include "MitAna/DataTree/interface/Types.h"
77
78 namespace mithep
79 {
80 class Track : public DataObject
81 {
82 public:
83 enum HitLayer { PXB1,
84 PXB2,
85 PXB3,
86 PXF1,
87 PXF2,
88 TIB1,
89 TIB1S,
90 TIB2,
91 TIB2S,
92 TIB3,
93 TIB4,
94 TID1,
95 TID1S,
96 TID2,
97 TID2S,
98 TID3,
99 TOB1,
100 TOB1S,
101 TOB2,
102 TOB2S,
103 TOB3,
104 TOB4,
105 TOB5,
106 TOB6,
107 TEC1,
108 TEC1S,
109 TEC2,
110 TEC2S,
111 TEC3,
112 TEC4,
113 TEC5,
114 TEC5S,
115 TEC6,
116 TEC7,
117 TEC8,
118 TEC9 };
119
120 Track() : fQOverP(0), fQOverPErr(0), fLambda(0), fLambdaErr(0),
121 fPhi0(0), fPhi0Err(0), fDxy(0), fDxyErr(0), fDsz(0), fDszErr(0),
122 fChi2(0), fNdof(0) {}
123 Track(Double_t qOverP, Double_t lambda, Double_t phi0, Double_t dxy, Double_t dsz) :
124 fQOverP(qOverP), fQOverPErr(0), fLambda(lambda), fLambdaErr(0),
125 fPhi0(phi0), fPhi0Err(0), fDxy(dxy), fDxyErr(0), fDsz(dsz), fDszErr(0),
126 fChi2(0), fNdof(0) {}
127 ~Track() {}
128
129 Double_t QOverP() const { return fQOverP; }
130 Double_t QOverPErr() const { return fQOverPErr; }
131 Double_t Lambda() const { return fLambda; }
132 Double_t LambdaErr() const { return fLambdaErr; }
133 Double_t Phi0() const { return fPhi0; }
134 Double_t Phi0Err() const { return fPhi0Err; }
135 Double_t Dxy() const { return fDxy; }
136 Double_t DxyErr() const { return fDxyErr; }
137 Double_t Dsz() const { return fDsz; }
138 Double_t DszErr() const { return fDszErr; }
139
140
141
142 Int_t Charge() const { return (fQOverP>0) ? 1 : -1; }
143 Double_t Chi2() const { return fChi2; }
144 void ClearHit(HitLayer l) { fHits.ClearBit(l); }
145 Double_t D0() const { return -fDxy; }
146 Double_t D0Err() const { return fDxyErr; }
147 Bool_t Hit(HitLayer l) const { return fHits.TestBit(l); }
148 BitMask64 &Hits() { return fHits; }
149 const BitMask64 &Hits() const { return fHits; }
150 ULong64_t HitMask() const { return fHits.Bits(); }
151 ThreeVector Mom() const { return ThreeVector(Px(),Py(),Pz()); }
152 UInt_t Ndof() const { return fNdof; }
153 Double_t P2() const { return P()*P(); }
154 Double_t P() const { return TMath::Abs(1./fQOverP); }
155 Double_t Px() const { return Pt()*TMath::Cos(fPhi0); }
156 Double_t Py() const { return Pt()*TMath::Sin(fPhi0); }
157 Double_t Pz() const { return P()*TMath::Sin(fLambda); }
158 Double_t Phi() const { return fPhi0; }
159 Double_t Pt() const { return TMath::Abs(TMath::Cos(fLambda)/fQOverP); }
160 //Double_t PtErr() const { return fPtErr; }
161 void SetChi2(Double_t chi2) { fChi2 = chi2; }
162 void SetHit(HitLayer l) { fHits.SetBit(l); }
163 void SetHits(BitMask64 hits) { fHits = hits; }
164 void SetHits(ULong64_t hitMask) { fHits.SetBits(hitMask); }
165 void SetNdof(UInt_t dof) { fNdof = dof; }
166 void SetStat(BitMask32 stat) { fStat = stat; }
167 void SetStat(UInt_t statBits) { fStat.SetBits(statBits); }
168 BitMask32 &Stat() { return fStat; }
169 const BitMask32 &Stat() const { return fStat; }
170 UInt_t StatBits() const { return fStat.Bits(); }
171 Double_t Theta() const { return (TMath::PiOver2() - fLambda); }
172 Double_t Z0() const { return fDsz/TMath::Cos(fLambda); }
173 //Double_t Z0Err() const { return fZ0Err; }
174
175 FourVector Mom4(Double_t m) const { return FourVector(Px(),Py(),Pz(),E(m)); }
176 Double_t E2(Double_t m) const { return P2()+m*m; }
177 Double_t E(Double_t m) const { return TMath::Sqrt(E2(m)); }
178 UInt_t NHits() const { return fHits.NBitsSet(); }
179
180 void SetHelix (Double_t qOverP, Double_t lambda, Double_t phi0,
181 Double_t dXy, Double_t dSz);
182 void SetErrors(Double_t qOverPErr, Double_t lambdaErr, Double_t phi0Err,
183 Double_t dXyErr, Double_t dSzErr);
184
185 const MCParticle *MCPart() const;
186 void SetMCPart(MCParticle *p) { fMCParticleRef = p; }
187
188 protected:
189 // Constant which is store in the file
190 BitMask64 fHits; // Mostly Hit informations
191 BitMask32 fStat; // Storage for various interesting things
192 Double_t fQOverP, fQOverPErr;
193 Double_t fLambda, fLambdaErr;
194 Double_t fPhi0,fPhi0Err; // Follow track parameters/uncertainties
195 Double_t fDxy, fDxyErr;
196 Double_t fDsz, fDszErr;
197
198 Double_t fChi2; //chi squared of track fit
199 UInt_t fNdof; //number of dof of track fit
200
201 TRef fMCParticleRef; //reference to sim particle (for monte carlo)
202
203 ClassDef(Track, 1) // Track class
204 };
205 }
206
207 //--------------------------------------------------------------------------------------------------
208 inline
209 void mithep::Track::SetHelix(Double_t qOverP, Double_t lambda, Double_t phi0,
210 Double_t dxy, Double_t dsz)
211 {
212 // Set helix parameters.
213
214 fQOverP = qOverP;
215 fLambda = lambda;
216 fPhi0 = phi0;
217 fDxy = dxy;
218 fDsz = dsz;
219 }
220
221 //--------------------------------------------------------------------------------------------------
222 inline
223 void mithep::Track::SetErrors(Double_t qOverPErr, Double_t lambdaErr, Double_t phi0Err,
224 Double_t dxyErr, Double_t dszErr)
225 {
226 // Set helix errors.
227
228 fQOverPErr = qOverPErr;
229 fLambdaErr = lambdaErr;
230 fPhi0Err = phi0Err;
231 fDxyErr = dxyErr;
232 fDszErr = dszErr;
233 }
234
235 //--------------------------------------------------------------------------------------------------
236 inline
237 const mithep::MCParticle *mithep::Track::MCPart() const
238 {
239 // Get reference to simulated particle.
240
241 return static_cast<const MCParticle*>(fMCParticleRef.GetObject());
242 }
243 #endif