StFemtoTrack.h 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284
  1. /**
  2. * \class StFemtoTrack
  3. * \brief Holds track information
  4. *
  5. * The StFemtoTrack class describes track (primary or global). The class
  6. * is designed in the both standalone and specific environment modes.
  7. * For the former one needs _VANILLA_ROOT_ global variable
  8. * to be defined at the compilation and run times.
  9. *
  10. * \author Grigory Nigmatkulov; e-mail: nigmatkulov@gmail.com
  11. * \date June 28, 2018
  12. */
  13. #ifndef StFemtoTrack_h
  14. #define StFemtoTrack_h
  15. // ROOT headers
  16. #include "TObject.h"
  17. #include "TVector3.h"
  18. #include "TBufferFile.h"
  19. #include "TMath.h"
  20. // FemtoDst headers
  21. #include "StFemtoHelix.h"
  22. #include "StFemtoPhysicalHelix.h"
  23. #ifdef _VANILLA_ROOT_
  24. #include "SystemOfUnits.h"
  25. #include "PhysicalConstants.h"
  26. #else
  27. #include "StarClassLibrary/SystemOfUnits.h"
  28. #include "StarClassLibrary/PhysicalConstants.h"
  29. #endif
  30. //_________________
  31. class StFemtoTrack : public TObject {
  32. public:
  33. /// Default constructor
  34. StFemtoTrack();
  35. /// Copy constructor
  36. StFemtoTrack(const StFemtoTrack &copy);
  37. /// Destructor
  38. virtual ~StFemtoTrack() { /* empty */}
  39. /// Print track information
  40. void print();
  41. //
  42. // Setters
  43. //
  44. /// Set track unique ID
  45. void setId(const Short_t& id) { mId = (UShort_t)id; }
  46. /// Set number of hits with next encoding: nHits * charge
  47. void setNHits(const Char_t& nHits) { mNHits = nHits; }
  48. /// Set number of hits with next encoding: nHits * charge
  49. void setNHits(const Short_t& nHits) { mNHits = (Char_t)nHits; }
  50. /// Set number of hits with next encoding: nHits * charge
  51. void setNHits(const Int_t& nHits) { mNHits = (Char_t)nHits; }
  52. /// Set possible number of hits in TPC
  53. void setNHitsPoss(const UShort_t& nHits) { mNHitsPoss = (UChar_t)nHits; }
  54. /// Set number of hits that were used in dE/dx estimated
  55. void setNHitsDedx(const UShort_t& nHits) { mNHitsDedx = (UChar_t)mNHitsDedx; }
  56. /// Set chi2 of the track reconstruction
  57. void setChi2(const Float_t& chi2);
  58. /// Set nSigma(electron)
  59. void setNSigmaElectron(const Float_t& ns);
  60. /// Set nSigma(pion)
  61. void setNSigmaPion(const Float_t& ns);
  62. /// Set nSigma(kaon)
  63. void setNSigmaKaon(const Float_t& ns);
  64. /// Set nSigma(proton)
  65. void setNSigmaProton(const Float_t& ns);
  66. /// Set dE/dx
  67. void setDedx(const Float_t& dEdx);
  68. /// Set dE/dx
  69. void setDedx(const Double_t& dEdx);
  70. /// Set the first (out of two) word of topology map
  71. void setMap0(const Int_t& map0) { mMap[0] = (UInt_t)map0; }
  72. /// Set the second (out of two) word of topology map
  73. void setMap1(const Int_t& map1) { mMap[1] = (UInt_t)map1; }
  74. /// Set topology map
  75. void setTopologyMap(const Int_t& i, const Int_t& map) { if (i==0 || i==1) { mMap[i]=(UInt_t)map; } }
  76. /// Set X component of primary track momentum
  77. void setPx(const Float_t& px) { mPrimMomX = px; }
  78. /// Set Y component of primary track momentum
  79. void setPy(const Float_t& py) { mPrimMomY = py; }
  80. /// Set Z component of primary track momentum
  81. void setPz(const Float_t& pz) { mPrimMomZ = pz; }
  82. /// Set three-momentum of the primary track (px,py,pz)
  83. void setP(const Float_t& px, const Float_t& py, const Float_t& pz)
  84. { setPrimaryMomentum( px, py, pz ); }
  85. /// Set three-momentum of the primary track (px,py,pz)
  86. void setPrimaryMomentum(const Float_t& px, const Float_t& py, const Float_t& pz)
  87. {mPrimMomX = px; mPrimMomY = py; mPrimMomZ = pz; }
  88. /// Set position (x,y,z) of the track at DCA to primary vertex
  89. void setOriginXYZ(const Float_t& x, const Float_t& y, const Float_t& z)
  90. { mOriginX = x; mOriginY = y; mOriginZ = z; }
  91. /// Set X position of the track origin at DCA to primary vertex
  92. void setOriginX(const Float_t& x) { mOriginX = x; }
  93. /// Set Y position of the track origin at DCA to primary vertex
  94. void setOriginY(const Float_t& y) { mOriginY = y; }
  95. /// Set Z position of the track origin at DCA to primary vertex
  96. void setOriginZ(const Float_t& z) { mOriginZ = z; }
  97. /// Set momentum of the global track (px,py,pz) at DCA to primary vertex
  98. void setGlobalP(const Float_t& px, const Float_t& py, const Float_t& pz)
  99. { mGlobMomX = px; mGlobMomY = py; mGlobMomZ = pz; }
  100. /// Set TOF beta
  101. void setBeta(const Float_t& beta);
  102. //
  103. // Getters
  104. //
  105. /// Return track unique ID
  106. Short_t id() const { return (Short_t)mId; }
  107. /// Return track type: 0 - global; 1 - primary
  108. Short_t type() const { return ( isPrimary() ) ? 1 : 0; }
  109. /// Check if the track is primary
  110. Bool_t isPrimary() const { return (pMom().Mag()>0) ? true : false; }
  111. /// Return track flag. Not store it but 300 represents
  112. /// that track was reconstructed in TPC. Need it for historical reasons
  113. Short_t flag() const { return 300; }
  114. /// Return number of hits
  115. UShort_t nHits() const { return (UShort_t)TMath::Abs(mNHits); }
  116. /// Return number of hits used in a track fit
  117. UShort_t nHitsFit() const { return nHits(); }
  118. /// Return possible number of hits in TPC
  119. UShort_t nHitsPoss() const { return (UShort_t)mNHitsPoss; }
  120. /// Return number of hits used in dE/dx estimation
  121. UShort_t nHitsDedx() const { return (UShort_t)mNHitsDedx; }
  122. /// Return nSigma(electron)
  123. Float_t nSigmaElectron() const { return (Float_t)mNSigmaElectron / 1000.; }
  124. /// Return nSigma(pion)
  125. Float_t nSigmaPion() const { return (Float_t)mNSigmaPion / 1000.; }
  126. /// Return nSigma(kaon)
  127. Float_t nSigmaKaon() const { return (Float_t)mNSigmaKaon / 1000.; }
  128. /// Return nSigma(proton)
  129. Float_t nSigmaProton() const { return (Float_t)mNSigmaProton / 1000.; }
  130. /// Dummy method that returns probability to be an electron
  131. Float_t pidProbElectron() const { return 0.5; }
  132. /// Dummy method that returns probability to be an pion
  133. Float_t pidProbPion() const { return 0.5; }
  134. /// Dummy method that returns probability to be an kaon
  135. Float_t pidProbKaon() const { return 0.5; }
  136. /// Dummy method that returns probability to be an proton
  137. Float_t pidProbProton() const { return 0.5; }
  138. /// Return dE/dx of the track in [GeV * cm]
  139. Double_t dEdx() const { return (Double_t)mDedx * 1e-9; }
  140. /// Return the first word of topology map
  141. UInt_t map0() const { return (UInt_t)mMap[0]; }
  142. /// Return the second word of topology map
  143. UInt_t map1() const { return (UInt_t)mMap[1]; }
  144. /// Return the word of topology map
  145. UInt_t map(const Int_t& word) { return (word==0) ? mMap[0] : mMap[1]; }
  146. /// Return chi2 of the track
  147. Float_t chi2() const { return (Float_t)mChi2 / 10.; }
  148. /// Return charge of the track
  149. Short_t charge() const { return (mNHits>0) ? +1 : -1; }
  150. /// Return origin position (x,y,z) of the global track at DCA to the primary vertex
  151. TVector3 origin() const { return TVector3( mOriginX, mOriginY, mOriginZ); }
  152. /// Return transverse distance between (pointX,pointY) and (originX,originY)
  153. Float_t gDCAxy(const Float_t& pointX, const Float_t& pointY) const;
  154. /// Return longitudinal distance between (pointZ) and (originZ)
  155. Float_t gDCAz(const Float_t& pointZ) const;
  156. /// Return distance between point (pointX,pointY,pointZ) and origin (originX,originY,originZ)
  157. Float_t gDCA(const Float_t& pVtxX, const Float_t& pVtxY, const Float_t& pVtxZ) const;
  158. /// Return 3D distance between point (pointX,pointY,pointZ) and origin (originX,originY,originZ)
  159. TVector3 gDCA(TVector3 point) const;
  160. /// Return 3-momentum of the primary track. (0,0,0) if the track is global
  161. TVector3 pMom() const { return TVector3(mPrimMomX, mPrimMomY, mPrimMomZ); }
  162. /// Return magnitude of primary track 3-momentum. 0 if the track is global
  163. Float_t p() const { return ptot(); }
  164. /// Return magnitude of primary track 3-momentum. 0 if the track is global
  165. Float_t ptot() const { return ( isPrimary() ) ? pMom().Mag() : 0.; }
  166. /// Return squared magnitude of primary track 3-momentum. 0 if the track is global
  167. Float_t ptot2() const { return ( isPrimary() ) ? pMom().Mag2() : 0.; }
  168. /// Return transverse momentum of primary track. 0 if the track is global
  169. Float_t pt() const { return ( isPrimary() ) ? pMom().Perp() : 0.; }
  170. /// Return squared transverse momentum of primary track. 0 if the track is global
  171. Float_t pt2() const { return ( isPrimary() ) ? pMom().Perp2() : 0.; }
  172. /// Return pseudorapidity of primary track. 0 if the track is global
  173. Float_t eta() const { return ( isPrimary() ) ? pMom().PseudoRapidity() : 0.; }
  174. /// Return azimuthal angle of primary track. 0 if the track is global
  175. Float_t phi() const { return ( isPrimary() ) ? pMom().Phi() : 0.; }
  176. /// Return the 3-momentum of global track at the dca point to the
  177. /// point (usually it is the primary vertex). B - magnetic field
  178. /// from FemtoEvent::bField()
  179. TVector3 gMom() const { return TVector3(mGlobMomX, mGlobMomY, mGlobMomZ); }
  180. /// Return magnitude of global track 3-momentum
  181. Float_t gP() const { return gPtot(); }
  182. /// Return magnitude of global track 3-momentum
  183. Float_t gPtot() const { return gMom().Mag(); }
  184. /// Return squared magnitude of global track 3-momentum
  185. Float_t gPtot2() const { return gMom().Mag2(); }
  186. /// Return transverse momentum of global track
  187. Float_t gPt() const { return gMom().Perp(); }
  188. /// Return squared transverse momentum of global track
  189. Float_t gPt2() const { return gMom().Perp2(); }
  190. /// Return pseudorapidity of global track
  191. Float_t gEta() const { return gMom().PseudoRapidity(); }
  192. /// Return beta measured in TOF
  193. Float_t beta() const { return ( isTofTrack() ) ? (Float_t)mBeta / 20000. : -999.; }
  194. /// Return 1./beta
  195. Float_t invBeta() const { return ( isTofTrack() ) ? 1./beta() : -999.; }
  196. /// Return 1./beta^2
  197. Float_t invBeta2() const { return ( isTofTrack() ) ? invBeta()*invBeta() : -999.; }
  198. /// Return squared mass estimated by TOF beta
  199. Float_t massSqr() const;
  200. /// Check if track had signal in TOF
  201. Bool_t isTofTrack() const { return (mBeta>0) ? true : false; }
  202. /// Global momentum at point of DCA to pVtx, B should be in kilogauss
  203. TVector3 gMom(const TVector3 pVtx, const Float_t B) const;
  204. /// Return helix of the global track
  205. StFemtoPhysicalHelix helix(const Float_t& B) const;
  206. private:
  207. /// Track unique ID
  208. UShort_t mId;
  209. /// Number of hits * charge
  210. Char_t mNHits;
  211. /// Number of possible hits
  212. UChar_t mNHitsPoss;
  213. /// Number of hits used for dE/dx estimation
  214. UChar_t mNHitsDedx;
  215. /// dEdx * 10^9 (in eV/cm)
  216. UShort_t mDedx;
  217. /// Topology map (word0 and word1)
  218. UInt_t mMap[2];
  219. /// chi2 of the track reconstruction (covnersion factor = *10)
  220. UChar_t mChi2;
  221. /// nSigma(electron) via dE/dx (conversion factor = *1000)
  222. Short_t mNSigmaElectron;
  223. /// nSigma(pion) via dE/dx (conversion factor = *1000)
  224. Short_t mNSigmaPion;
  225. /// nSigma(kaon) via dE/dx (conversion factor = *1000)
  226. Short_t mNSigmaKaon;
  227. /// nSigma(proton) via dE/dx (conversion factor = *1000)
  228. Short_t mNSigmaProton;
  229. /// Beta in TOF (0 if no signal) (conversion factor = *20000)
  230. UShort_t mBeta;
  231. /// X component of track momentum fitted to the primary vertex
  232. Float_t mPrimMomX;
  233. /// Y component of track momentum fitted to the primary vertex
  234. Float_t mPrimMomY;
  235. /// Z component of track momentum fitted to the primary vertex
  236. Float_t mPrimMomZ;
  237. /// X component of global track momentum at DCA to primary vertex
  238. Float_t mGlobMomX;
  239. /// Y component of global track momentum at DCA to primary vertex
  240. Float_t mGlobMomY;
  241. /// Z component of global track momentum at DCA to primary vertex
  242. Float_t mGlobMomZ;
  243. /// X position of global track at DCA to primary vertex
  244. Float_t mOriginX;
  245. /// Y position of global track at DCA to primary vertex
  246. Float_t mOriginY;
  247. /// Z position of global track at DCA to primary vertex
  248. Float_t mOriginZ;
  249. ClassDef(StFemtoTrack, 2)
  250. };
  251. //_________________
  252. inline TVector3 StFemtoTrack::gMom(const TVector3 pVtx, const Float_t B) const {
  253. StFemtoPhysicalHelix gHelix = helix(B);
  254. return gHelix.momentumAt( gHelix.pathLength(pVtx), B * kilogauss );
  255. }
  256. //_________________
  257. inline StFemtoPhysicalHelix StFemtoTrack::helix(const Float_t& B) const {
  258. return StFemtoPhysicalHelix( TVector3(mGlobMomX, mGlobMomY, mGlobMomZ),
  259. TVector3(mOriginX, mOriginY, mOriginZ),
  260. B * kilogauss,
  261. static_cast<float>( charge() ) );
  262. }
  263. #endif // StFemtoTrack_h