StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StPicoTrack.h
1 
8 #ifndef StPicoTrack_h
9 #define StPicoTrack_h
10 
11 // C++ headers
12 #include <cmath>
13 
14 // ROOT headers
15 #include "TObject.h"
16 #include "TVector3.h"
17 
18 // PicoDst headers
19 #include "StPicoHelix.h"
20 #include "StPicoPhysicalHelix.h"
21 
22 #if defined (_VANILLA_ROOT_)
23 #include "SystemOfUnits.h"
24 #include "PhysicalConstants.h"
25 #else
26 #include "StarClassLibrary/SystemOfUnits.h"
27 #include "StarClassLibrary/PhysicalConstants.h"
28 #endif
29 
30 #if defined (__TFG__VERSION__)
31 #include "StPicoTrackCovMatrix.h"
32 #endif /* __TFG__VERSION__ */
33 
34 //_________________
35 class StPicoTrack : public TObject {
36 
37  public:
39  StPicoTrack();
43  virtual ~StPicoTrack();
45  virtual void Print(const Char_t *option = "") const;
46 
47 #if defined (__TFG__VERSION__)
48  enum { eTopologyMap = 3 };
49 #else /* ! __TFG__VERSION__ */
50  enum { eTopologyMap = 2 };
51 #endif /* __TFG__VERSION__ */
52 
53  //
54  // Getters
55  //
56 
58  Int_t id() const { return mId; }
60  Float_t chi2() const { return mChi2 / 1000.f; }
62  TVector3 pMom() const { return TVector3(mPMomentumX, mPMomentumY, mPMomentumZ); }
64  TVector3 gMom() const { return TVector3(mGMomentumX, mGMomentumY, mGMomentumZ); }
66  TVector3 origin() const { return TVector3(mOriginX, mOriginY, mOriginZ); }
68  Float_t gPt() const { return gMom().Perp(); }
70  Float_t gPtot() const { return gMom().Mag(); }
72  Float_t pPt() const { return isPrimary() ? pMom().Perp() : 0.; }
74  Float_t pPtot() const { return isPrimary() ? pMom().Mag() : 0.; }
76  TVector3 gMom(TVector3 pVtx, Float_t B) const;
78  StPicoPhysicalHelix helix(Float_t const B) const;
79 
80  // Next functions return DCA (or its components) of the global track
81  // to the point with coordinates (pVtxX, pVtxY, pVtxZ)
82 
84  Float_t gDCAx(Float_t pVtxX) const { return (mOriginX - pVtxX); }
86  Float_t gDCAy(Float_t pVtxY) const { return (mOriginY - pVtxY); }
88  Float_t gDCAz(Float_t pVtxZ) const { return (mOriginZ - pVtxZ); }
90  Float_t gDCAxy(Float_t pVtxX, Float_t pVtxY) const;
92  Float_t gDCA(Float_t pVtxX, Float_t pVtxY, Float_t pVtxZ) const;
94  TVector3 gDCA(TVector3 pVtx) const;
100  Float_t gDCAs(TVector3 pVtx) const;
102  Short_t charge() const { return (mNHitsFit > 0) ? 1 : -1; }
104  Int_t nHits() const { return (mNHitsFit > 0) ? (Int_t)mNHitsFit : (Int_t)(-1 * mNHitsFit); }
106  Int_t nHitsFit() const { return (mNHitsFit > 0) ? (Int_t)mNHitsFit : (Int_t)(-1 * mNHitsFit); }
108  Int_t nHitsMax() const { return (Int_t)mNHitsMax; }
110  Int_t nHitsPoss() const { return nHitsMax(); }
112  Int_t nHitsDedx() const { return (Int_t)mNHitsDedx; }
114  UInt_t hftHitsMap() const { return topologyMap(0) >> 1 & 0x7F; }
116  Float_t dEdx() const { return mDedx; }
118  Float_t dEdxError() const { return mDedxError; }
119 
120 #if defined (__TFG__VERSION__)
121  Float_t dEdxPull(Float_t mass, UChar_t fit = 1, Int_t charge = 1) const;
122  Float_t dEdxPullToF(Float_t mass, UChar_t fit = 1, Int_t charge = 1) const;
123  Float_t dEdxPullPion() const { return dEdxPull(0.13956995,1); }
124  Float_t dEdxPullKaon() const { return dEdxPull(0.493677,1); }
125  Float_t dEdxPullProton() const { return dEdxPull(0.93827231,1); }
126  Float_t dEdxPullElectron() const { return dEdxPull(0.51099907e-3,1); }
127 #endif
128  Float_t dNdx() const { return mDnDx; }
131  Float_t dNdxError() const { return mDnDxError; }
133  Char_t status() const { return mStatus; }
134 
136  Float_t nSigmaPion() const { return (Float_t)mNSigmaPion / 1000.f; }
138  Float_t nSigmaKaon() const { return (Float_t)mNSigmaKaon / 1000.f; }
140  Float_t nSigmaProton() const { return (Float_t)mNSigmaProton / 1000.f; }
142  Float_t nSigmaElectron() const { return (Float_t)mNSigmaElectron / 1000.f; }
143 
145  UInt_t topologyMap(UInt_t idx) const { return (idx>1) ? 0 : mTopologyMap[idx]; }
146 #if !defined (__TFG__VERSION__)
147  ULong64_t iTpcTopologyMap() const { return mTopoMap_iTpc; }
149 #endif
150 
151 
153  Bool_t hasPxl1Hit() const { return hftHitsMap() >> 0 & 0x1; }
155  Bool_t hasPxl2Hit() const { return hftHitsMap() >> 1 & 0x3; }
157  Bool_t hasIstHit() const { return hftHitsMap() >> 3 & 0x3; }
159  Bool_t hasSstHit() const { return hftHitsMap() >> 5 & 0x3; }
161  Bool_t isHft() const { return hasPxl1Hit() && hasPxl2Hit() && (hasIstHit() || hasSstHit()); }
163  Bool_t isHFTTrack() const { return isHft(); }
165  Bool_t hasHft4Layers() const { return hasPxl1Hit() && hasPxl2Hit() && hasIstHit() && hasSstHit(); }
167  Bool_t isTofTrack() const { return (mBTofPidTraitsIndex<0) ? false : true; }
169  Bool_t isBemcTrack() const { return (mBEmcPidTraitsIndex<0) ? false : true; }
171  Bool_t isMtdTrack() const { return (mMtdPidTraitsIndex<0) ? false : true; }
173  Bool_t isETofTrack() const { return (mETofPidTraitsIndex<0) ? false : true; }
175  Bool_t isBemcMatchedTrack() const { return mBEmcMatchedTowerIndex != 0; }
176 
178  Bool_t isPrimary() const { return ( pMom().Mag()>0 ); }
179 
181  Int_t bemcPidTraitsIndex() const { return mBEmcPidTraitsIndex; }
183  Int_t bTofPidTraitsIndex() const { return mBTofPidTraitsIndex; }
185  Int_t mtdPidTraitsIndex() const { return mMtdPidTraitsIndex; }
187  Int_t eTofPidTraitsIndex() const { return mETofPidTraitsIndex; }
188 
198  Int_t bemcTowerIndex() const { return abs(mBEmcMatchedTowerIndex) - 1; }
201  Bool_t isBemcMatchedExact() const { return mBEmcMatchedTowerIndex > 0;}
203  Int_t idTruth() const { return mIdTruth; }
205  Int_t qaTruth() const { return mQATruth; }
207  Int_t vertexIndex() const { return (Int_t)mVertexIndex; }
208 
209  //
210  // Setters
211  //
212 
214  void setId(Int_t id) { mId = (UShort_t)id; }
216  void setChi2(Float_t chi2);
218  void setPrimaryMomentum(Double_t px, Double_t py, Double_t pz)
219  { mPMomentumX = (Float_t)px; mPMomentumY = (Float_t)py; mPMomentumZ = (Float_t)pz; }
221  void setPrimaryMomentum(Float_t px, Float_t py, Float_t pz)
222  { mPMomentumX = px; mPMomentumY = py; mPMomentumZ = pz; }
224  void setPrimaryMomentum(TVector3 mom)
225  { mPMomentumX = (Float_t)mom.X(); mPMomentumY = (Float_t)mom.Y(); mPMomentumZ = (Float_t)mom.Z(); }
227  void setGlobalMomentum(Double_t px, Double_t py, Double_t pz)
228  { mGMomentumX = (Float_t)px; mGMomentumY = (Float_t)py; mGMomentumZ = (Float_t)pz; }
230  void setGlobalMomentum(Float_t px, Float_t py, Float_t pz)
231  { mGMomentumX = px; mGMomentumY = py; mGMomentumZ = pz; }
233  void setGlobalMomentum(TVector3 mom)
234  { mGMomentumX = (Float_t)mom.X(); mGMomentumY = (Float_t)mom.Y(); mGMomentumZ = (Float_t)mom.Z(); }
236  void setOrigin(Double_t x, Double_t y, Double_t z)
237  { mOriginX = (Float_t)x; mOriginY = (Float_t)y; mOriginZ = (Float_t)z; }
239  void setOrigin(Float_t x, Float_t y, Float_t z)
240  { mOriginX = x; mOriginY = y; mOriginZ = z; }
242  void setOrigin(TVector3 orig)
243  { mOriginX = (Float_t)orig.X(); mOriginY = (Float_t)orig.Y(); mOriginZ = (Float_t)orig.Z(); }
244 
246  void setDedx(Float_t dEdx);
249 
251  void setDndx(Float_t dNdx) { mDnDx = dNdx;}
255  void setStatus(Char_t k = 0) { mStatus = k; }
256 
258  void setNHitsFit(Int_t nhits) { mNHitsFit = (Char_t)nhits; }
260  void setNHitsPossible(Int_t nhits);
262  void setNHitsMax(Int_t nhits);
264  void setNHitsDedx(Int_t nhits);
266  void setNSigmaPion(Float_t ns);
268  void setNSigmaKaon(Float_t ns);
270  void setNSigmaProton(Float_t ns);
272  void setNSigmaElectron(Float_t ns);
274  void setTopologyMap(Int_t id, UInt_t word);
275 #if !defined (__TFG__VERSION__)
276  void setiTpcTopologyMap(ULong64_t map) { mTopoMap_iTpc = map; }
278 #endif
279 
281  void setBEmcPidTraitsIndex(Int_t index) { mBEmcPidTraitsIndex = (Short_t)index; }
283  void setBTofPidTraitsIndex(Int_t index) { mBTofPidTraitsIndex = (Short_t)index; }
285  void setMtdPidTraitsIndex(Int_t index) { mMtdPidTraitsIndex = (Short_t)index; }
287  void setETofPidTraitsIndex(Int_t index) { mETofPidTraitsIndex = (Short_t)index; }
289  void setBEmcMatchedTowerIndex(Int_t index) { mBEmcMatchedTowerIndex = (Short_t)index; }
291  void setMcTruth(Int_t index, Int_t qa) { mIdTruth = (UShort_t)index; mQATruth = (UShort_t)qa; }
293  void setVertexIndex(Int_t index);
294 
295  protected:
296 
298  UShort_t mId;
300  UShort_t mChi2;
302  Float_t mPMomentumX;
304  Float_t mPMomentumY;
306  Float_t mPMomentumZ;
308  Float_t mGMomentumX;
310  Float_t mGMomentumY;
312  Float_t mGMomentumZ;
314  Float_t mOriginX;
316  Float_t mOriginY;
318  Float_t mOriginZ;
319 
321  Float16_t mDedx;
323  Float16_t mDedxError;
325  Float_t mDnDx;
327  Float_t mDnDxError;
328 
330  Char_t mNHitsFit;
332  UChar_t mNHitsMax;
334  UChar_t mNHitsDedx;
336  Short_t mNSigmaPion;
338  Short_t mNSigmaKaon;
340  Short_t mNSigmaProton;
344  UInt_t mTopologyMap[eTopologyMap];
345 
359 
360 #if !defined (__TFG__VERSION__)
361  ULong64_t mTopoMap_iTpc;
363 #endif
364 
368  Char_t mStatus;
369 
371  UShort_t mIdTruth;
373  UShort_t mQATruth;
375  Char_t mVertexIndex;
376 
377 #if !defined (__TFG__VERSION__)
378  ClassDef(StPicoTrack, 9)
379 #else
380  ClassDef(StPicoTrack, 10)
381 #endif
382 };
383 
384 #endif
Float_t gPtot() const
Return global total momentum (GeV/c)
Definition: StPicoTrack.h:70
Bool_t isPrimary() const
Return if track is primary.
Definition: StPicoTrack.h:178
Float_t dEdxError() const
Return dE/dx error of the track (in GeV/cm)
Definition: StPicoTrack.h:118
UChar_t mNHitsMax
Possible number of hits (in TPC)
Definition: StPicoTrack.h:332
ULong64_t mTopoMap_iTpc
Topology map for the iTPC.
Definition: StPicoTrack.h:362
void setNSigmaKaon(Float_t ns)
Set nSigma(kaon)
Short_t mBEmcPidTraitsIndex
Index of the BEMC pidTratis in the event.
Definition: StPicoTrack.h:347
UInt_t topologyMap(UInt_t idx) const
Return track topology map (return 0 in case when requested index is &gt;1)
Definition: StPicoTrack.h:145
Float_t mOriginX
Track origin x (DCAx to the primary vertex) in cm.
Definition: StPicoTrack.h:314
void setStatus(Char_t k=0)
Set status of the track (0 - not fitted to any vertex, 1 - fitted to a vertex)
Definition: StPicoTrack.h:255
StPicoTrack()
Default constructor.
Definition: StPicoTrack.cxx:27
Short_t mNSigmaPion
nSigma(pion) (encoding = nsigma * 1000)
Definition: StPicoTrack.h:336
Short_t mNSigmaElectron
nSigma(electron) (encoding = nsigma * 1000)
Definition: StPicoTrack.h:342
Float_t mDnDx
Fitted dN/dx.
Definition: StPicoTrack.h:325
Float_t mPMomentumY
Py momentum (GeV/c) of the primary track ( 0 if not primary )
Definition: StPicoTrack.h:304
Float_t nSigmaKaon() const
Return nSigma(kaon)
Definition: StPicoTrack.h:138
Short_t mNSigmaProton
nSigma(proton) (encoding = nsigma * 1000)
Definition: StPicoTrack.h:340
UShort_t mChi2
Chi2 of the track (encoding = chi2*1000)
Definition: StPicoTrack.h:300
Float_t pPt() const
Return transverse momentum (GeV/c) of the primary track.
Definition: StPicoTrack.h:72
Float_t nSigmaPion() const
Return nSigma(pion)
Definition: StPicoTrack.h:136
void setTopologyMap(Int_t id, UInt_t word)
Set topology map (2 words)
Bool_t hasSstHit() const
Return if the track has an inner SST hit.
Definition: StPicoTrack.h:159
void setDedx(Float_t dEdx)
Set dE/dx of the track.
Bool_t hasHft4Layers() const
Return if track has hits in all 4 layers of silicon.
Definition: StPicoTrack.h:165
Float_t dNdxError() const
Return dN/dx error of the track.
Definition: StPicoTrack.h:131
void setNSigmaProton(Float_t ns)
Set nSigma(proton)
Short_t mMtdPidTraitsIndex
Index of the MTD pidTratis in the event.
Definition: StPicoTrack.h:351
Bool_t isMtdTrack() const
Return if track has MTD hit.
Definition: StPicoTrack.h:171
UChar_t mNHitsDedx
Number of hits used for dE/dx estimation (in TPC)
Definition: StPicoTrack.h:334
void setiTpcTopologyMap(ULong64_t map)
Set iTPC topology map.
Definition: StPicoTrack.h:277
void setGlobalMomentum(TVector3 mom)
Set momentum of the global track.
Definition: StPicoTrack.h:233
Int_t nHitsPoss() const
Return number of hits possible.
Definition: StPicoTrack.h:110
Float_t mOriginY
Track origin y (DCAy to the primary vertex) in cm.
Definition: StPicoTrack.h:316
void setOrigin(Float_t x, Float_t y, Float_t z)
Set origin of the track (DCA point to the primary vertex)
Definition: StPicoTrack.h:239
Float16_t mDedxError
dE/dx error (in GeV/cm)
Definition: StPicoTrack.h:323
Float_t mDnDxError
Fitted dN/dx error.
Definition: StPicoTrack.h:327
void setMcTruth(Int_t index, Int_t qa)
Set index of the corresonding MC track.
Definition: StPicoTrack.h:291
void setNSigmaPion(Float_t ns)
Set nSigma(pion)
Bool_t isHft() const
Return if trach has hit in HFT.
Definition: StPicoTrack.h:161
Int_t eTofPidTraitsIndex() const
Return index to the corresponding ETOF PID trait.
Definition: StPicoTrack.h:187
void setMtdPidTraitsIndex(Int_t index)
Set index to MTD PID traits.
Definition: StPicoTrack.h:285
Int_t bemcTowerIndex() const
Definition: StPicoTrack.h:198
TVector3 gMom() const
Return momentum (GeV/c) of the global tracks at the point of DCA to the primary vertex.
Definition: StPicoTrack.h:64
Float_t dNdx() const
Return dN/dx of the track.
Definition: StPicoTrack.h:129
Float_t gDCAs(TVector3 pVtx) const
void setOrigin(Double_t x, Double_t y, Double_t z)
Set origin of the track (DCA point to the primary vertex)
Definition: StPicoTrack.h:236
void setChi2(Float_t chi2)
Set chi2 of the track.
Int_t bemcPidTraitsIndex() const
Return index to the corresponding BEMC PID trait.
Definition: StPicoTrack.h:181
Bool_t isBemcMatchedExact() const
Definition: StPicoTrack.h:201
Float_t gDCAy(Float_t pVtxY) const
Return signed distance in y direction between the value and y position of the DCA point (gDCAy - y) ...
Definition: StPicoTrack.h:86
Holds information about track parameters.
Definition: StPicoTrack.h:35
Short_t mNSigmaKaon
nSigma(kaon) (encoding = nsigma * 1000)
Definition: StPicoTrack.h:338
UInt_t hftHitsMap() const
Return a map of hits in HFT.
Definition: StPicoTrack.h:114
Bool_t hasPxl1Hit() const
Return if the track has an inner PXL hit.
Definition: StPicoTrack.h:153
Bool_t hasPxl2Hit() const
Return if the track has an outer PXL hit.
Definition: StPicoTrack.h:155
void setNHitsPossible(Int_t nhits)
Set nHitsPoss.
Float_t gDCAx(Float_t pVtxX) const
Return signed distance in x direction (cm) between the value and x position of the DCA point (gDCAx -...
Definition: StPicoTrack.h:84
void setNSigmaElectron(Float_t ns)
Set nSigma(electron)
TVector3 origin() const
Return space coordinate (cm) of DCA to the primary vertex.
Definition: StPicoTrack.h:66
Float_t mGMomentumX
Px component of the momentum (GeV/c) of the global track at DCA to primary vertex.
Definition: StPicoTrack.h:308
UShort_t mIdTruth
MC track id.
Definition: StPicoTrack.h:371
Float_t chi2() const
Return chi2 of the track.
Definition: StPicoTrack.h:60
Short_t mBTofPidTraitsIndex
Index of the BTOF pidTratis in the event.
Definition: StPicoTrack.h:349
UShort_t mQATruth
MC track quality (percentage of hits coming from corresponding MC track)
Definition: StPicoTrack.h:373
Float_t gDCAxy(Float_t pVtxX, Float_t pVtxY) const
Return distance in xy direction (cm) between the (x,y) point and the DCA point to primary vertex...
void setPrimaryMomentum(Float_t px, Float_t py, Float_t pz)
Set momentum of the primary track.
Definition: StPicoTrack.h:221
void setDedxError(Float_t dEdxError)
Set dE/dx error of the track.
Definition: StPicoTrack.h:248
Float_t nSigmaElectron() const
Return nSigma(electron)
Definition: StPicoTrack.h:142
void setNHitsFit(Int_t nhits)
Set nHitsFit ( charge * nHitsFit )
Definition: StPicoTrack.h:258
Bool_t isHFTTrack() const
Return if trach has hit in HFT.
Definition: StPicoTrack.h:163
void setETofPidTraitsIndex(Int_t index)
Set index to ETOF PID traits.
Definition: StPicoTrack.h:287
void setNHitsMax(Int_t nhits)
Set nHitsPoss.
void setBTofPidTraitsIndex(Int_t index)
Set index to BTOF PID traits.
Definition: StPicoTrack.h:283
Char_t mStatus
Definition: StPicoTrack.h:368
Int_t nHitsDedx() const
Return number of hits used for dE/dx measurement.
Definition: StPicoTrack.h:112
Char_t status() const
Return if the track was fitted to any vertex (0 - not fitted, 1 - fitted)
Definition: StPicoTrack.h:133
Int_t vertexIndex() const
Return parent vertex index (-2 if not fitted to any vertex)
Definition: StPicoTrack.h:207
Float_t gPt() const
Return global tranverse momentum (GeV/c)
Definition: StPicoTrack.h:68
void setNHitsDedx(Int_t nhits)
Set nHitsDedx.
void setBEmcPidTraitsIndex(Int_t index)
Set index to BEMC PID traits.
Definition: StPicoTrack.h:281
void setDndx(Float_t dNdx)
Set dN/dx of the track.
Definition: StPicoTrack.h:251
Float_t mPMomentumZ
Pz momentum (GeV/c) of the primary track ( 0 if not primary )
Definition: StPicoTrack.h:306
TVector3 pMom() const
Return momentum (GeV/c) of the primary track. Return (0,0,0) if not primary track.
Definition: StPicoTrack.h:62
Float_t dEdx() const
Return dE/dx (in keV/cm) of the track.
Definition: StPicoTrack.h:116
Helis parametrization for the particle.
void setVertexIndex(Int_t index)
Set vertex index to which the track was fitted.
Float_t pPtot() const
Return total momentum (GeV/c) of the primary track.
Definition: StPicoTrack.h:74
Bool_t isBemcTrack() const
Return if track has BEMC hit.
Definition: StPicoTrack.h:169
Float_t mGMomentumZ
Pz component of the momentum (GeV/c) of the global track at DCA to primary vertex.
Definition: StPicoTrack.h:312
void setGlobalMomentum(Double_t px, Double_t py, Double_t pz)
Set momentum of the global track.
Definition: StPicoTrack.h:227
Float_t mGMomentumY
Py component of the momentum (GeV/c) of the global track at DCA to primary vertex.
Definition: StPicoTrack.h:310
Float_t nSigmaProton() const
Return nSigma(proton)
Definition: StPicoTrack.h:140
Bool_t isBemcMatchedTrack() const
Checks if track matches any BEMC tower.
Definition: StPicoTrack.h:175
void setId(Int_t id)
Set track ID.
Definition: StPicoTrack.h:214
Bool_t isETofTrack() const
Return if track has ETOF hit.
Definition: StPicoTrack.h:173
Float16_t mDedx
dE/dx in KeV/cm (dE/dx * 1e6)
Definition: StPicoTrack.h:321
void setPrimaryMomentum(Double_t px, Double_t py, Double_t pz)
Set momentum of the primary track.
Definition: StPicoTrack.h:218
Int_t nHitsFit() const
Return number of hits fit.
Definition: StPicoTrack.h:106
Float_t gDCA(Float_t pVtxX, Float_t pVtxY, Float_t pVtxZ) const
Return distance in xyz direction (cm) between the (x,y,z) point and the DCA point to primary vertex...
Int_t qaTruth() const
Qualtiy of the MC track.
Definition: StPicoTrack.h:205
Float_t mOriginZ
Track origin z (DCAy to the primary vertex) in cm.
Definition: StPicoTrack.h:318
Float_t gDCAz(Float_t pVtxZ) const
Return signed distance in z direction (cm) between the value and z position of the DCA point (gDCAz -...
Definition: StPicoTrack.h:88
virtual ~StPicoTrack()
Destructor.
Definition: StPicoTrack.cxx:95
Bool_t hasIstHit() const
Return if the track has an IST hit.
Definition: StPicoTrack.h:157
Short_t mETofPidTraitsIndex
Index of the ETOF pidTratis in the event.
Definition: StPicoTrack.h:353
Int_t idTruth() const
Index of the corresponding MC track.
Definition: StPicoTrack.h:203
void setGlobalMomentum(Float_t px, Float_t py, Float_t pz)
Set momentum of the global track.
Definition: StPicoTrack.h:230
Int_t nHitsMax() const
Return number of hits possible.
Definition: StPicoTrack.h:108
Bool_t isTofTrack() const
Return if track has TOF hit.
Definition: StPicoTrack.h:167
UInt_t mTopologyMap[eTopologyMap]
Toplogy Map data0 and data1. See StEvent/StTrackTopologyMap.cxx.
Definition: StPicoTrack.h:344
void setBEmcMatchedTowerIndex(Int_t index)
Set index of the BEMC tower that matches the track.
Definition: StPicoTrack.h:289
Int_t id() const
Return unique Id of the track.
Definition: StPicoTrack.h:58
Char_t mNHitsFit
Charge * nHitsFit.
Definition: StPicoTrack.h:330
void setOrigin(TVector3 orig)
Set origin of the track (DCA point to the primary vertex)
Definition: StPicoTrack.h:242
Int_t mtdPidTraitsIndex() const
Return index to the corresponding MTD PID trait.
Definition: StPicoTrack.h:185
virtual void Print(const Char_t *option="") const
Print track parameters.
UShort_t mId
Unique track ID.
Definition: StPicoTrack.h:298
Int_t bTofPidTraitsIndex() const
Return index to the corresponding BTOF PID trait.
Definition: StPicoTrack.h:183
Int_t nHits() const
Return number of hits.
Definition: StPicoTrack.h:104
Short_t mBEmcMatchedTowerIndex
Definition: StPicoTrack.h:358
Float_t mPMomentumX
Px momentum (GeV/c) of the primary track ( 0 if not primary )
Definition: StPicoTrack.h:302
Char_t mVertexIndex
Parent vertex index. -2 if no vertex.
Definition: StPicoTrack.h:375
void setPrimaryMomentum(TVector3 mom)
Set momentum of the primary track.
Definition: StPicoTrack.h:224
StPicoPhysicalHelix helix(Float_t const B) const
Helix at point of DCA to StPicoEvent::mPrimaryVertex.
void setDndxError(Float_t dNdxError)
Set dN/dx error of the track.
Definition: StPicoTrack.h:253
Short_t charge() const
Return charge of the track (encoded in nHitsFit as: nHitsFit * charge)
Definition: StPicoTrack.h:102
ULong64_t iTpcTopologyMap() const
Return topology map for iTPC.
Definition: StPicoTrack.h:148