Belle II Software  release-06-02-00
TrackFitResult.h
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 #pragma once
9 
10 #include <framework/gearbox/Const.h>
11 #include <framework/datastore/RelationsObject.h>
12 #include <framework/dataobjects/Helix.h>
13 #include <framework/dataobjects/UncertainHelix.h>
14 #include <framework/geometry/BFieldManager.h>
15 
16 #include <TVector3.h>
17 #include <TMatrixDSym.h>
18 #include <TLorentzVector.h>
19 
20 #include <stdint.h>
21 #include <vector>
22 #include <cmath>
23 
24 namespace Belle2 {
29  class HitPatternCDC;
30  class HitPatternVXD;
31 
50  public:
53 
72  TrackFitResult(const TVector3& position, const TVector3& momentum,
73  const TMatrixDSym& covariance, const short int charge,
74  const Const::ParticleType& particleType, const float pValue,
75  const float bField,
76  const uint64_t hitPatternCDCInitializer,
77  const uint32_t hitPatternVXDInitializer,
78  const uint16_t NDF);
79 
91  TrackFitResult(const std::vector<float>& tau, const std::vector<float>& cov5,
92  const Const::ParticleType& particleType, const float pValue,
93  const uint64_t hitPatternCDCInitializer,
94  const uint32_t hitPatternVXDInitializer,
95  const uint16_t NDF
96  );
97 
99  TVector3 getPosition() const { return getHelix().getPerigee(); }
100 
106  TVector3 getMomentum() const
107  {
108  const double bField = BFieldManager::getField(getPosition()).Z() / Unit::T;
109  return getHelix().getMomentum(bField);
110  }
111 
115  TLorentzVector get4Momentum() const
116  {
117  return TLorentzVector(getMomentum(), getEnergy());
118  }
119 
123  double getEnergy() const
124  {
125  return std::sqrt(getMomentum().Mag2() + getParticleType().getMass() * getParticleType().getMass());
126  }
127 
130  double getTransverseMomentum() const
131  {
132  const double bField = BFieldManager::getField(getPosition()).Z() / Unit::T;
133  return getHelix().getTransverseMomentum(bField);
134  }
135 
137  TMatrixDSym getCovariance6() const
138  {
139  const double bField = BFieldManager::getField(getPosition()).Z() / Unit::T;
140  return getUncertainHelix().getCartesianCovariance(bField);
141  }
142 
145 
150  short getChargeSign() const { return (getOmega() > 0) - (getOmega() < 0); }
151 
153  double getPValue() const { return m_pValue; }
154 
156  int getNDF() const;
157 
159  double getChi2() const;
160 
161  //------------------------------------------------------------------------
162  // --- Getters for perigee helix parameters
163  //------------------------------------------------------------------------
168  double getD0() const { return m_tau[iD0]; }
169 
174  double getPhi0() const { return m_tau[iPhi0]; }
175 
177  double getPhi() const { return getPhi0(); }
178 
184  double getOmega() const { return m_tau[iOmega]; }
185 
190  double getZ0() const { return m_tau[iZ0]; }
191 
196  double getTanLambda() const { return m_tau[iTanLambda]; }
197 
199  double getCotTheta() const { return getTanLambda(); }
200 
205  std::vector<float> getTau() const { return std::vector<float>(m_tau, m_tau + c_NPars); }
206 
211  std::vector<float> getCov() const { return std::vector<float>(m_cov5, m_cov5 + c_NCovEntries); }
212 
217  TMatrixDSym getCovariance5() const;
218 
220  Helix getHelix() const
221  { return Helix(getD0(), getPhi0(), getOmega(), getZ0(), getTanLambda()); }
222 
224  Helix getHelix(float momentumScale) const
225  { return Helix(getD0(), getPhi0(), getOmega() / momentumScale, getZ0(), getTanLambda()); }
226 
229  {
230  return UncertainHelix(getD0(), getPhi0(), getOmega(), getZ0(),
232  }
233 
236 
239 
241  virtual std::string getInfoHTML() const override;
242 
244  private:
245 
246  //---------------------------------------------------------------------------------------------------------------------------
248  const unsigned int m_pdg;
249 
251  const Double32_t m_pValue;
252 
258  static const unsigned int c_NPars = 5;
259  static const unsigned int c_NCovEntries = 5 * 6 / 2;
265  static const unsigned int iD0 = 0;
266  static const unsigned int iPhi0 = 1;
267  static const unsigned int iOmega = 2;
268  static const unsigned int iZ0 = 3;
269  static const unsigned int iTanLambda = 4;
273  Double32_t m_tau[c_NPars];
274 
279  Double32_t m_cov5[c_NCovEntries];
280 
291 
293  static const uint16_t c_NDFFlag = 0xFFFF;
294 
296  uint16_t m_NDF;
297 
299  /* Version history:
300  ver 8: add NDF
301  ver 7: fixed sign errors in the translation of position and momentum covariances.
302  ver 6: use fixed size arrays instead of vectors (add schema evolution rule), use Double32_t.
303  ver 5: CDC Hit Pattern now a single variable (add schema evolution rule).
304  ver 4: added hit pattern
305  ver 3: back to vectors for the elements
306  ver 4: fixed-size arrays for the elements <------- incompatible with later version 4
307  ver 3: add data fields for hit patterns
308  ver 2: re-arranging of members
309  ver 1:
310  ver 2: <------- incompatible with later version 2
311  ver 1:
312  */
313  };
315 }
Helix parameter class.
Definition: Helix.h:48
The ParticleType class for identifying different particle types.
Definition: Const.h:289
Hit pattern of CDC hits within a track.
Definition: HitPatternCDC.h:35
Hit pattern of the VXD within a track.
Definition: HitPatternVXD.h:37
Defines interface for accessing relations of objects in StoreArray.
Values of the result of a track fit with a given particle hypothesis.
std::vector< float > getTau() const
Getter for all perigee parameters.
ClassDefOverride(TrackFitResult, 8)
Values of the result of a track fit with a given particle hypothesis.
static const unsigned int c_NPars
Number of helix parameters.
Helix getHelix() const
Conversion to framework Helix (without covariance).
TMatrixDSym getCovariance5() const
Getter for covariance matrix of perigee parameters in matrix form.
static const unsigned int iZ0
Index for z0.
static const unsigned int c_NCovEntries
Number of covariance entries.
Helix getHelix(float momentumScale) const
Conversion to framework Helix with momentum scaling (without covariance).
double getPhi() const
Getter for phi0 with CDF naming convention.
Double32_t m_tau[c_NPars]
perigee helix parameters; tau = d0, phi0, omega, z0, tanLambda.
double getChi2() const
Get chi2 given NDF and p-value.
short getChargeSign() const
Return track charge (1 or -1).
double getPValue() const
Getter for Chi2 Probability of the track fit.
TLorentzVector get4Momentum() const
Getter for the 4Momentum at the closest approach of the track in the r/phi projection.
double getEnergy() const
Getter for the Energy at the closest approach of the track in the r/phi projection.
TMatrixDSym getCovariance6() const
Position and Momentum Covariance Matrix.
const uint32_t m_hitPatternVXDInitializer
Member for initializing the information about hits in the VXD.
double getCotTheta() const
Getter for tanLambda with CDF naming convention.
virtual std::string getInfoHTML() const override
Return a short summary of this object's contents in HTML format.
double getOmega() const
Getter for omega.
TVector3 getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
int getNDF() const
Getter for number of degrees of freedom of the track fit.
const Double32_t m_pValue
Chi2 Probability of the fit.
double getD0() const
Getter for d0.
double getTransverseMomentum() const
Getter for the absolute value of the transverse momentum at the perigee.
double getTanLambda() const
Getter for tanLambda.
static const unsigned int iTanLambda
Index tan lambda.
TrackFitResult()
Constructor initializing everything to zero.
std::vector< float > getCov() const
Getter for all covariance matrix elements of perigee parameters.
TVector3 getPosition() const
Getter for vector of position at closest approach of track in r/phi projection.
static const unsigned int iD0
Index for d0.
uint64_t m_hitPatternCDCInitializer
Member for initializing the information about hits in the CDC.
static const unsigned int iOmega
Index for omega.
Const::ParticleType getParticleType() const
Getter for ParticleType of the mass hypothesis of the track fit.
static const uint16_t c_NDFFlag
backward compatibility initialisation for NDF
double getZ0() const
Getter for z0.
uint16_t m_NDF
Memeber for number of degrees of freedom.
Double32_t m_cov5[c_NCovEntries]
The 15 = 5*(5+1)/2 covariance matrix elements.
static const unsigned int iPhi0
Index for phi0.
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
UncertainHelix getUncertainHelix() const
Conversion to framework Uncertain Helix (i.e., with covariance).
const unsigned int m_pdg
PDG Code for hypothesis with which the corresponding fit was performed.
double getPhi0() const
Getter for phi0.
HitPatternVXD getHitPatternVXD() const
Getter for the hit pattern in the VXD;.
This class represents an ideal helix in perigee parameterization including the covariance matrix of t...
TMatrixDSym getCartesianCovariance(const double bZ_tesla=1.5) const
Getter for the position and momentum covariance matrix.
static const double T
[tesla]
Definition: Unit.h:120
static void getField(const double *pos, double *field)
return the magnetic field at a given position.
Abstract base class for different kinds of events.