StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StPicoTrack.cxx
1 //
2 // StPicoTrack holds information about the reconstructed tracks
3 //
4 
5 // C+++ headers
6 #include <limits>
7 
8 // ROOT headers
9 #include "TMath.h"
10 
11 // PicoDst headers
12 #include "StPicoMessMgr.h"
13 #include "StPicoTrack.h"
14 
15 #if defined (__TFG__VERSION__)
16 #include "TF1.h"
17 #include "St_base/StMessMgr.h"
18 #include "StMuDSTMaker/COMMON/StMuTrack.h"
19 #include "StBichsel/Bichsel.h"
20 #include "StBichsel/StdEdxModel.h"
21 #include "StBichsel/StdEdxPull.h"
22 #endif /* __TFG__VERSION__ */
23 
24 ClassImp(StPicoTrack)
25 
26 //_________________
27 StPicoTrack::StPicoTrack() : TObject(),
28  mId(0),
29  mChi2(std::numeric_limits<unsigned short>::max()),
30  mPMomentumX(0), mPMomentumY(0), mPMomentumZ(0),
31  mGMomentumX(0), mGMomentumY(0), mGMomentumZ(0),
32  mOriginX(0),mOriginY(0), mOriginZ(0),
33  mDedx(0), mDedxError(0),
34  mDnDx(0),mDnDxError(0),
35  mNHitsFit(0), mNHitsMax(0), mNHitsDedx(0),
36  mNSigmaPion( std::numeric_limits<short>::min() ),
37  mNSigmaKaon( std::numeric_limits<short>::min() ),
38  mNSigmaProton( std::numeric_limits<short>::min() ),
39  mNSigmaElectron( std::numeric_limits<short>::min() ),
40  mTopologyMap{}, mBEmcPidTraitsIndex(-1), mBTofPidTraitsIndex(-1),
41  mMtdPidTraitsIndex(-1), mETofPidTraitsIndex(-1),
42  mBEmcMatchedTowerIndex(-1)
43 #if !defined (__TFG__VERSION__)
44  , mTopoMap_iTpc(0),
45 #endif
46  mStatus(0),
47  mIdTruth(0), mQATruth(0), mVertexIndex(-1)
48  {
49  // Default constructor
50  /* empty */
51 }
52 
53 //_________________
55 
56  mId = track.mId;
57  mChi2 = track.mChi2;
58  mPMomentumX = track.mPMomentumX;
59  mPMomentumY = track.mPMomentumY;
60  mPMomentumZ = track.mPMomentumZ;
61  mGMomentumX = track.mGMomentumX;
62  mGMomentumY = track.mGMomentumY;
63  mGMomentumZ = track.mGMomentumZ;
64  mOriginX = track.mOriginX;
65  mOriginY = track.mOriginY;
66  mOriginZ = track.mOriginZ;
67  mDedx = track.mDedx;
68  mDedxError = track.mDedxError;
69  mDnDx = track.mDnDx;
70  mDnDxError = track.mDnDxError;
71  mNHitsFit = track.mNHitsFit;
72  mNHitsMax = track.mNHitsMax;
73  mNHitsDedx = track.mNHitsDedx;
74  mNSigmaPion = track.mNSigmaPion;
75  mNSigmaKaon = track.mNSigmaKaon;
78  for(int iIter=0; iIter<eTopologyMap; iIter++) {
79  mTopologyMap[iIter] = track.mTopologyMap[iIter];
80  }
85 #if !defined (__TFG__VERSION__)
87 #endif
88  mStatus = track.mStatus;
89  mIdTruth = track.mIdTruth;
90  mQATruth = track.mQATruth;
91  mVertexIndex = track.mVertexIndex;
92 }
93 
94 //_________________
96  /* emtpy */
97 }
98 
99 //_________________
100 void StPicoTrack::Print(const Char_t* option __attribute__((unused))) const {
101  LOG_INFO << "id: " << id() << " chi2: " << chi2() << "\n"
102  << "pMom: " << pMom().X() << " " << pMom().Y() << " " << pMom().Z() << "\n"
103  << "gMom: " << gMom().X() << " " << gMom().Y() << " " << gMom().Z() << "\n"
104  << "origin: " << origin().X() << " " << origin().Y() << " " << origin().Z() << "\n"
105  << "nHitsFit: " << nHitsFit()
106  << " nHitsdEdx: " << nHitsDedx() << "\n"
107  << "nSigma pi/K/p/e: " << nSigmaPion() << "/" << nSigmaKaon() << "/"
108  << nSigmaProton() << "/" << nSigmaElectron() << "\n"
109  << "Hit index in BEMC/BTof/MTD/ETof: " << mBEmcPidTraitsIndex << "/"
110  << mBTofPidTraitsIndex << "/" << mMtdPidTraitsIndex << "/" << mETofPidTraitsIndex << "\n"
111  << "idTruth: " << idTruth() << " qaTruth: " << qaTruth() << "\n"
112  << endm;
113 }
114 
115 //_________________
116 Float_t StPicoTrack::gDCAxy(Float_t x, Float_t y) const {
117  return TMath::Sqrt( (mOriginX-x)*(mOriginX-x) + (mOriginY-y)*(mOriginY-y) );
118 }
119 
120 //_________________
121 Float_t StPicoTrack::gDCA(Float_t x, Float_t y, Float_t z) const {
122  return TMath::Sqrt( (mOriginX-x)*(mOriginX-x) +
123  (mOriginY-y)*(mOriginY-y) +
124  (mOriginZ-z)*(mOriginZ-z) );
125 }
126 
127 //_________________
128 TVector3 StPicoTrack::gDCA(TVector3 pVtx) const {
129  return (origin() - pVtx);
130 }
131 
132 //_________________
133 void StPicoTrack::setChi2(Float_t chi2) {
134  mChi2 = ( (chi2 * 1000.) > std::numeric_limits<unsigned short>::max() ?
135  std::numeric_limits<unsigned short>::max() :
136  (UShort_t)( TMath::Nint( chi2 * 1000. ) ) );
137 }
138 
139 //_________________
140 void StPicoTrack::setDedx(Float_t dEdx) {
141  // In keV/cm
142  mDedx = dEdx * 1.e6;
143 }
144 
145 //_________________
146 void StPicoTrack::setNHitsMax(Int_t nhits) {
147  mNHitsMax = (UChar_t)nhits;
148 }
149 
150 //_________________
151 void StPicoTrack::setNHitsPossible(Int_t nhits) {
152  // For those who wants to have standard terminology
153  setNHitsMax(nhits);
154 }
155 
156 //_________________
157 void StPicoTrack::setNHitsDedx(Int_t nhits) {
158  mNHitsDedx = (UChar_t)nhits;
159 }
160 
161 //_________________
162 void StPicoTrack::setTopologyMap(Int_t id, UInt_t word) {
163  if(id==0 || id==1 || id==2) {
164  mTopologyMap[id] = word;
165  }
166  else {
167  // Shouldn't here be a protection?
168  }
169 }
170 
171 //_________________
172 void StPicoTrack::setNSigmaPion(Float_t ns) {
173  mNSigmaPion = ( fabs(ns * 1000.) > std::numeric_limits<short>::max() ?
174  ( (ns > 0) ? std::numeric_limits<short>::max() :
175  std::numeric_limits<short>::min() ) :
176  (Short_t)( TMath::Nint( ns * 1000. ) ) );
177 }
178 
179 //_________________
180 void StPicoTrack::setNSigmaKaon(Float_t ns) {
181  mNSigmaKaon = ( fabs(ns * 1000.) > std::numeric_limits<short>::max() ?
182  ( (ns > 0) ? std::numeric_limits<short>::max() :
183  std::numeric_limits<short>::min() ) :
184  (Short_t)( TMath::Nint( ns * 1000. ) ) );
185 }
186 
187 //_________________
189  mNSigmaProton = ( fabs(ns * 1000.) > std::numeric_limits<short>::max() ?
190  ( (ns > 0) ? std::numeric_limits<short>::max() :
191  std::numeric_limits<short>::min() ) :
192  (Short_t)( TMath::Nint( ns * 1000. ) ) );
193 }
194 
195 //_________________
197  mNSigmaElectron = ( fabs(ns * 1000.) > std::numeric_limits<short>::max() ?
198  ( (ns > 0) ? std::numeric_limits<short>::max() :
199  std::numeric_limits<short>::min() ) :
200  (Short_t)( TMath::Nint( ns * 1000. ) ) );
201 }
202 
203 //_________________
204 TVector3 StPicoTrack::gMom(TVector3 pVtx, Float_t const B) const {
205  StPicoPhysicalHelix gHelix = helix(B);
206  return gHelix.momentumAt( gHelix.pathLength( pVtx ), B * kilogauss );
207 }
208 
209 //_________________
210 StPicoPhysicalHelix StPicoTrack::helix(Float_t const B) const {
211  return StPicoPhysicalHelix( gMom(), origin(), B * kilogauss,
212  static_cast<float>( charge() ) );
213 }
214 
215 #if defined (__TFG__VERSION__)
216 //_________________
217 Float_t StPicoTrack::dEdxPull(Float_t mass, UChar_t fit, Int_t charge) const {
218  Float_t z = -999.;
219  Float_t momentum = gMom().Mag();
220  Float_t betagamma = momentum * TMath::Abs(charge) / mass;
221  Float_t dedx_measured, dedx_resolution = -1;
222  if (! fit) { // I70
223  dedx_measured = 1e-6*dEdx();
224  dedx_resolution = dEdxError();
225  }
226  else if ( fit == 1) { // Ifit
227  dedx_measured = 1e-6*dEdx();
228  dedx_resolution = dEdxError();
229  }
230  else { // dNdx
231  dedx_measured = dNdx();
232  dedx_resolution = dNdxError();
233  }
234  if (dedx_resolution <= 0) return z;
235  z = StdEdxPull::Eval(dedx_measured,dedx_resolution,betagamma,fit,charge);
236  return z;
237 }
238 
239 //_________________
240 Float_t StPicoTrack::dEdxPullToF(Float_t mass, UChar_t fit, Int_t charge) const {
241  Float_t z = -999.;
242  Float_t momentum = gMom().Mag();
243  Float_t betagamma = momentum * TMath::Abs(charge) / mass;
244  Float_t dedx_measured, dedx_resolution = -1;
245  Float_t dedx_predicted = -1;
246  static TF1 *ToFCor[2][3] = {0};
247 
248 #if 1
249  if (! ToFCor[0][0]) {
250  /* Maksym, 11/22/17
251  dE/dx correction Pion:
252  p0 = -0.0583149 +/- 0.000381198
253  p1 = 0.418013 +/- 0.00692832
254  p2 = -2.13934 +/- 0.0494599
255  p3 = 8.50604 +/- 0.182424
256  p4 = -19.3807 +/- 0.384702
257  p5 = 24.777 +/- 0.480388
258  p6 = -17.7556 +/- 0.350305
259  p7 = 6.67151 +/- 0.137513
260  p8 = -1.02484 +/- 0.0224013
261 
262  ToFCor[0][0] = new TF1("dEdxCorPion","pol8",-0.05,1.7);
263  ToFCor[0][0]->SetParameters(-0.0583149, 0.418013, -2.13934, 8.50604, -19.3807,
264  24.777, -17.7556, 6.67151, -1.02484);
265  */
266 
267  ToFCor[0][0] = new TF1("dEdxCorPion","pol8",0,1.4);
268  ToFCor[0][0]->SetParameters(-0.0686941, 1.00512, -7.01202, 26.5077, -55.6186,
269  67.2198, -46.5966, 17.194, -2.61537);
270  /*
271  Kaon:
272  p0 = 0.0159834 +/- 5.0974e-05
273  p1 = -0.115402 +/- 0.00053059
274  p2 = 0.429126 +/- 0.00363743
275  p3 = 0.569058 +/- 0.021112
276  p4 = -3.00964 +/- 0.0443932
277  p5 = 1.96151 +/- 0.223594
278  p6 = 4.58199 +/- 0.54901
279  p7 = -7.40423 +/- 0.54138
280  p8 = 3.00975 +/- 0.190495
281 
282  ToFCor[0][1] = new TF1("dEdxCorKaon","pol8",-0.30,1.1);
283  ToFCor[0][1]->SetParameters(0.0159834,-0.115402, 0.429126, 0.569058, -3.00964,
284  1.96151, 4.58199, -7.40423, 3.00975);
285  */
286 
287  ToFCor[0][1] = new TF1("dEdxCorKaon","pol8",-0.30,1);
288  ToFCor[0][1]->SetParameters(0.0103183, -0.113391, 0.347118, 0.939172, -2.85298,
289  -2.16397, 12.1903, -12.2671, 3.9317);
290  /*
291  Proton:
292 
293  p0 = 0.0224256 +/- 3.17809e-05
294  p1 = -0.111028 +/- 0.000293992
295  p2 = 0.389417 +/- 0.00153563
296  p3 = 0.714018 +/- 0.00741038
297  p4 = -3.01122 +/- 0.0184805
298  p5 = -0.463794 +/- 0.0523714
299  p6 = 6.12037 +/- 0.082
300  p7 = -1.32531 +/- 0.100664
301  p8 = -2.51178 +/- 0.128703
302 
303  ToFCor[0][2] = new TF1("dEdxCorProt","pol8",-0.65,0.9);
304  ToFCor[0][2]->SetParameters(0.0224256,-0.111028, 0.389417, 0.714018, -3.01122,
305  -0.463794, 6.12037, -1.32531, -2.51178);
306  */
307 
308  ToFCor[0][2] = new TF1("dEdxCorProt","pol8",-0.5,0.75);
309  ToFCor[0][2]->SetParameters(0.0241608, -0.122085, 0.409931, 0.772181, -2.82003,
310  -1.15634, 5.60995, 0.425809, -3.1966);
311  /*
312  ================================================================================
313  dNdX correction: Pion:
314  p0 = -0.127892 +/- 0.000352695
315  p1 = 1.13166 +/- 0.00649171
316  p2 = -6.28713 +/- 0.0471639
317  p3 = 22.8583 +/- 0.177002
318  p4 = -48.3574 +/- 0.379105
319  p5 = 59.3737 +/- 0.479759
320  p6 = -41.8384 +/- 0.353823
321  p7 = 15.6944 +/- 0.140222
322  p8 = -2.42872 +/- 0.0230262
323 
324  ToFCor[1][0] = new TF1("dNdxCorPion","pol8",-0.05,1.7);
325  ToFCor[1][0]->SetParameters(-0.127892, 1.13166, -6.28713, 22.8583, -48.3574, 59.3737,
326  -41.8384, 15.6944,-2.42872);
327  */
328 
329  ToFCor[1][0] = new TF1("dNdxCorPion","pol8",0,1.4);
330  ToFCor[1][0]->SetParameters(-0.0790884, 0.960891, -7.18579, 29.6238, -65.671,
331  82.2409, -58.5028, 22.0458, -3.41535);
332  /*
333  Kaon:
334 
335  p0 = -0.0324453 +/- 6.52723e-05
336  p1 = 0.104644 +/- 0.000668736
337  p2 = 0.232821 +/- 0.00461983
338  p3 = -0.205617 +/- 0.0263313
339  p4 = -0.942393 +/- 0.0551658
340  p5 = 2.28523 +/- 0.279332
341  p6 = -1.43567 +/- 0.678682
342  p7 = -0.558613 +/- 0.664332
343  p8 = 0.608146 +/- 0.232324
344  12/04/17
345  0.0116721,-0.106367, 0.224163, 1.12377,-0.176388, -13.855, 27.2605, -13.9522, -1.64127
346 
347  ToFCor[1][1] = new TF1("dNdxCorKaon","pol8",-0.30,1.1);
348  ToFCor[1][1]->SetParameters(-0.0324453, 0.104644, 0.232821, -0.205617, -0.942393,
349  2.28523, -1.43567, -0.558613, 0.608146);
350  ToFCor[1][1] = new TF1("dNdxCorKaon","pol8",-0.30,1.1);
351  ToFCor[1][1]->SetParameters(0.0116721,-0.106367, 0.224163, 1.12377,-0.176388, -13.855,
352  27.2605, -13.9522, -1.64127);
353  */
354 
355 
356  ToFCor[1][1] = new TF1("dNdxCorKaon","pol8",-0.3,1);
357  ToFCor[1][1]->SetParameters(-0.0165074, 0.00880377, 0.181553, 0.723291, -1.69666,
358  -2.88222, 9.83642, -8.58367, 2.45608);
359  /*
360  Proton:
361  p0 = -0.0216107 +/- 4.06732e-05
362  p1 = 0.095628 +/- 0.000380259
363  p2 = 0.221496 +/- 0.00191754
364  p3 = -0.281845 +/- 0.00960846
365  p4 = -1.81734 +/- 0.0226373
366  p5 = 3.55535 +/- 0.0683598
367  p6 = 1.78527 +/- 0.100257
368  p7 = -8.3076 +/- 0.132161
369  p8 = 5.20165 +/- 0.161309
370 
371  ToFCor[1][2] = new TF1("dNdxCorProt","pol8",-0.65,0.9);
372  ToFCor[1][2]->SetParameters(-0.0216107, 0.095628, 0.221496, -0.281845, -1.81734,
373  3.55535, 1.78527, -8.3076, 5.20165);
374  */
375 
376  ToFCor[1][2] = new TF1("dNdxCorProt","pol8",-0.5,0.75);
377  ToFCor[1][2]->SetParameters(-0.0254615, 0.0610245, 0.364898, 0.0783227, -2.27553,
378  0.317491, 4.55173, -1.44629, -1.70712);
379  }
380 #endif
381 
382  Int_t hyp = -1;
383  if (mass < 1.0) hyp = 2;
384  if (mass < 0.5) hyp = 1;
385  if (mass < 0.2) hyp = 0;
386  if (mass < 0.1) hyp = -1; // no ocrrection for electron
387  Double_t bgL10 = TMath::Log10(betagamma);
388  if ( fit <= 1) { // I70
389  dedx_measured = 1e-6*dEdx();
390  dedx_resolution = dEdxError();
391  dedx_predicted = 1.e-6 * charge * charge * TMath::Exp(Bichsel::Instance()->GetMostProbableZ(bgL10));
392  if (ToFCor[0][hyp] &&hyp >= 0 && bgL10 > ToFCor[0][hyp]->GetXmin() &&
393  bgL10 < ToFCor[0][hyp]->GetXmax()) {
394  dedx_predicted *= TMath::Exp(ToFCor[0][hyp]->Eval(bgL10));
395  }
396  }
397  else { // dNdx
398  dedx_measured = dNdx();
399  dedx_resolution = dNdxError();
400  dedx_predicted = StdEdxModel::instance()->dNdx(betagamma,charge);
401  if (ToFCor[0][hyp] && hyp >= 0 && bgL10 > ToFCor[0][hyp]->GetXmin() &&
402  bgL10 < ToFCor[0][hyp]->GetXmax()) {
403  dedx_predicted *= TMath::Exp(ToFCor[1][hyp]->Eval(bgL10));
404  }
405  }
406  if (dedx_resolution <= 0) return z;
407  z = TMath::Log(dedx_measured/dedx_predicted) / dedx_resolution;
408  return z;
409 }
410 
411 #endif /* __TFG__VERSION__ */
412 
413 //_________________
414 Float_t StPicoTrack::gDCAs(TVector3 point) const {
415  // Signed DCA is defined for tracks with primary partners
416  // and with non-zero global track momentum
417  if ( (gMom().Mag() == 0 ) || ( !isPrimary() ) ) {
418  return -999;
419  }
420  // Momentum of the global track
421  TVector3 dir = gMom();
422  // Unit vector
423  dir = dir.Unit();
424  Float_t cosl = dir.Perp();
425  // Return DCA vector to the point (origin - point)
426  TVector3 dca = gDCA( point );
427  return -dir.Y()/cosl * dca.X() + dir.X()/cosl * dca.Y();
428 }
429 
430 //_________________
431 void StPicoTrack::setVertexIndex(Int_t index) {
432  if ( index<=-2 || index>std::numeric_limits<char>::max() ) {
433  mVertexIndex = -2;
434  }
435  else {
436  mVertexIndex = (Char_t)index;
437  }
438 }
TVector3 momentumAt(Double_t, Double_t) const
Return momemtum at S.
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
Float_t mOriginX
Track origin x (DCAx to the primary vertex) in cm.
Definition: StPicoTrack.h:314
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 nSigmaPion() const
Return nSigma(pion)
Definition: StPicoTrack.h:136
void setTopologyMap(Int_t id, UInt_t word)
Set topology map (2 words)
void setDedx(Float_t dEdx)
Set dE/dx of the track.
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
UChar_t mNHitsDedx
Number of hits used for dE/dx estimation (in TPC)
Definition: StPicoTrack.h:334
Float_t mOriginY
Track origin y (DCAy to the primary vertex) in cm.
Definition: StPicoTrack.h:316
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 setNSigmaPion(Float_t ns)
Set nSigma(pion)
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 setChi2(Float_t chi2)
Set chi2 of the track.
Holds information about track parameters.
Definition: StPicoTrack.h:35
Short_t mNSigmaKaon
nSigma(kaon) (encoding = nsigma * 1000)
Definition: StPicoTrack.h:338
void setNHitsPossible(Int_t nhits)
Set nHitsPoss.
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
pair< Double_t, Double_t > pathLength(Double_t r) const
path length at given r (cylindrical r)
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...
Float_t nSigmaElectron() const
Return nSigma(electron)
Definition: StPicoTrack.h:142
void setNHitsMax(Int_t nhits)
Set nHitsPoss.
Char_t mStatus
Definition: StPicoTrack.h:368
Int_t nHitsDedx() const
Return number of hits used for dE/dx measurement.
Definition: StPicoTrack.h:112
void setNHitsDedx(Int_t nhits)
Set nHitsDedx.
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 mGMomentumZ
Pz component of the momentum (GeV/c) of the global track at DCA to primary vertex.
Definition: StPicoTrack.h:312
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
Float16_t mDedx
dE/dx in KeV/cm (dE/dx * 1e6)
Definition: StPicoTrack.h:321
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
virtual ~StPicoTrack()
Destructor.
Definition: StPicoTrack.cxx:95
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
UInt_t mTopologyMap[eTopologyMap]
Toplogy Map data0 and data1. See StEvent/StTrackTopologyMap.cxx.
Definition: StPicoTrack.h:344
Int_t id() const
Return unique Id of the track.
Definition: StPicoTrack.h:58
Char_t mNHitsFit
Charge * nHitsFit.
Definition: StPicoTrack.h:330
virtual void Print(const Char_t *option="") const
Print track parameters.
UShort_t mId
Unique track ID.
Definition: StPicoTrack.h:298
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
StPicoPhysicalHelix helix(Float_t const B) const
Helix at point of DCA to StPicoEvent::mPrimaryVertex.
Short_t charge() const
Return charge of the track (encoded in nHitsFit as: nHitsFit * charge)
Definition: StPicoTrack.h:102