StFemtoTrack.h 5.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134
  1. //
  2. // The version of the class should be changed every time
  3. // when any changes in the codes are done
  4. // Grigory Nigmatkulov: 2016/12/15
  5. //
  6. #ifndef STFEMTOTRACK_HH
  7. #define STFEMTOTRACK_HH
  8. #include "Rtypes.h"
  9. #include "TObject.h"
  10. #include "TVector3.h"
  11. #include "TMath.h"
  12. const Float_t __SCALE__ = 1000.;
  13. const Float_t __INV_SCALE__ = 0.001;
  14. const UShort_t USHORTMAX = 65535;
  15. const Short_t SHORTMAX = 32767;
  16. //_________________
  17. class StFemtoTrack : public TObject {
  18. public:
  19. StFemtoTrack();
  20. StFemtoTrack(Float_t px, Float_t py, Float_t pz);
  21. ~StFemtoTrack();
  22. //Setters
  23. void SetId(Short_t id) { mId = (UShort_t)id; }
  24. void SetNHits(Char_t nHits) { mNHits = nHits; } //Should give Char_t
  25. void SetNHitsPoss(UShort_t nHits) { mNHitsPoss = (UChar_t)nHits; }
  26. //void SetNHitsDedx(UShort_t nHits) { mNHitsDedx = (UChar_t)nHits; }
  27. void SetNSigmaElectron(Float_t ns);
  28. void SetNSigmaPion(Float_t ns);
  29. void SetNSigmaKaon(Float_t ns);
  30. void SetNSigmaProton(Float_t ns);
  31. void SetDedx(Float_t dEdx);
  32. void SetDedx(Double_t dEdx);
  33. void SetMap0(Int_t map0) { mMap0 = (UInt_t)map0; }
  34. void SetMap0(UInt_t map0) { mMap0 = map0; }
  35. void SetMap1(Int_t map1) { mMap1 = (UInt_t)map1; }
  36. void SetMap1(UInt_t map1) { mMap1 = map1; }
  37. void SetPx(Float_t px) { mMomX = px; }
  38. void SetPy(Float_t py) { mMomY = py; }
  39. void SetPz(Float_t pz) { mMomZ = pz; }
  40. void SetP(Float_t px, Float_t py, Float_t pz);
  41. void SetDCAxGlobal(Float_t x);
  42. void SetDCAyGlobal(Float_t y);
  43. void SetDCAzGlobal(Float_t z);
  44. void SetDCAGlobal(Float_t x, Float_t y, Float_t z);
  45. void SetGlobalPx(Float_t px) { mGlobMomX = px; }
  46. void SetGlobalPy(Float_t py) { mGlobMomY = py; }
  47. void SetGlobalPz(Float_t pz) { mGlobMomZ = pz; }
  48. void SetGlobalP(Float_t px, Float_t py, Float_t pz);
  49. void SetBeta(Float_t beta) { mBeta = beta; }
  50. //Getters
  51. Short_t GetId() const { return (Short_t)mId; }
  52. Short_t GetType() const;
  53. Short_t GetFlag() const { return 300; }
  54. UShort_t GetNHits() const { return (UShort_t)TMath::Abs(mNHits); }
  55. UShort_t GetNHitsFit() const;
  56. UShort_t GetNHitsPoss() const { return (UShort_t)mNHitsPoss; }
  57. //UShort_t GetNHitsDedx() const { return (UShort_t)mNHitsDedx; }
  58. UShort_t GetNHitsDedx() const { return (UShort_t)TMath::Abs(mNHits); } //Temporary
  59. Float_t GetNSigmaElectron() const { return (Float_t)mNSigmaElectron*__INV_SCALE__; }
  60. Float_t GetNSigmaPion() const { return (Float_t)mNSigmaPion*__INV_SCALE__; }
  61. Float_t GetNSigmaKaon() const { return (Float_t)mNSigmaKaon*__INV_SCALE__; }
  62. Float_t GetNSigmaProton() const { return (Float_t)mNSigmaProton*__INV_SCALE__; }
  63. Float_t GetPidProbElectron() const { return 0.5; }
  64. Float_t GetPidProbPion() const { return 0.5; }
  65. Float_t GetPidProbKaon() const { return 0.5; }
  66. Float_t GetPidProbProton() const { return 0.5; }
  67. Double_t GetDedx() const { return (Double_t)mDedx * 1e-9; }
  68. UInt_t GetMap0() const { return (UInt_t)mMap0; }
  69. UInt_t GetMap1() const { return (UInt_t)mMap1; }
  70. Float_t GetChi2() const { return 1.; }
  71. Short_t GetCharge() const { return (mNHits>0) ? +1 : -1; }
  72. Float_t GetPx() const { return mMomX; }
  73. Float_t GetPy() const { return mMomY; }
  74. Float_t GetPz() const { return mMomZ; }
  75. TVector3 GetP() const;
  76. Float_t GetPt() const;
  77. Float_t GetPtot() const;
  78. Float_t GetEta() const;
  79. Float_t GetDCAxy() const { return 0.; }
  80. Float_t GetDCAx() const { return 0.; }
  81. Float_t GetDCAy() const { return 0.; }
  82. Float_t GetDCAz() const { return 0.; }
  83. Float_t GetDCAxyGlobal() const;
  84. Float_t GetDCAxGlobal() const { return (Float_t)mDCAxGlobal * 0.0001; }
  85. Float_t GetDCAyGlobal() const { return (Float_t)mDCAyGlobal * 0.0001; }
  86. Float_t GetDCAzGlobal() const { return (Float_t)mDCAzGlobal * 0.0001; }
  87. Float_t GetDCAGlobal() const;
  88. Float_t GetGlobalPx() const { return mGlobMomX; }
  89. Float_t GetGlobalPy() const { return mGlobMomY; }
  90. Float_t GetGlobalPz() const { return mGlobMomZ; }
  91. TVector3 GetGlobalP() const;
  92. Float_t GetPtGlobal() const;
  93. Float_t GetPtotGlobal() const;
  94. Float_t GetBeta() const { return mBeta; }
  95. Float_t GetMassSqr() const;
  96. Bool_t GetIsTofTrack();
  97. private:
  98. UShort_t mId;
  99. Char_t mNHits; //q * NHits ; nHitsFit=nHits+1 for primaries
  100. UChar_t mNHitsPoss;
  101. //UChar_t mNHitsDedx;
  102. Short_t mNSigmaElectron; // ns*1000 !
  103. Short_t mNSigmaPion; // ns*1000 !
  104. Short_t mNSigmaKaon; // ns*1000 !
  105. Short_t mNSigmaProton; // ns*1000 !
  106. UShort_t mDedx; // dEdx*1e9 !
  107. UInt_t mMap0; // TopologyMap data0
  108. UInt_t mMap1; // TopologyMap data1
  109. Float_t mMomX; // Primary Px at pVtx
  110. Float_t mMomY; // Primary Py at pVtx
  111. Float_t mMomZ; // Primary Pz at pVtx
  112. Short_t mDCAxGlobal; // gDCA * 10000 (IMPORTANT: for primaries)
  113. Short_t mDCAyGlobal; // gDCA * 10000 (IMPORTANT: for primaries)
  114. Short_t mDCAzGlobal; // gDCA * 10000 (IMPORTANT: for primaries)
  115. Float_t mGlobMomX; // Px Global at DCA to pVtx
  116. Float_t mGlobMomY; // Py Global at DCA to pVtx
  117. Float_t mGlobMomZ; // Pz Global at DCA to pVtx
  118. //To get the origin at gDCA (3D): mOrigin = PrimVtxPosition + gDCA
  119. Float_t mBeta; // ToF beta
  120. ClassDef(StFemtoTrack, 2)
  121. };
  122. #endif