Belle II Software  light-2205-abys
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 <Math/Vector4D.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  ROOT::Math::PxPyPzEVector get4Momentum() const
116  {
117  const B2Vector3D momentum = getMomentum();
118  return ROOT::Math::PxPyPzEVector(momentum.x(), momentum.y(), momentum.z(), getEnergy());
119  }
120 
124  double getEnergy() const
125  {
126  return std::sqrt(getMomentum().Mag2() + getParticleType().getMass() * getParticleType().getMass());
127  }
128 
131  double getTransverseMomentum() const
132  {
133  const double bField = BFieldManager::getField(getPosition()).Z() / Unit::T;
134  return getHelix().getTransverseMomentum(bField);
135  }
136 
138  TMatrixDSym getCovariance6() const
139  {
140  const double bField = BFieldManager::getField(getPosition()).Z() / Unit::T;
141  return getUncertainHelix().getCartesianCovariance(bField);
142  }
143 
146 
151  short getChargeSign() const { return (getOmega() > 0) - (getOmega() < 0); }
152 
154  double getPValue() const { return m_pValue; }
155 
157  int getNDF() const;
158 
160  double getChi2() const;
161 
162  //------------------------------------------------------------------------
163  // --- Getters for perigee helix parameters
164  //------------------------------------------------------------------------
169  double getD0() const { return m_tau[iD0]; }
170 
175  double getPhi0() const { return m_tau[iPhi0]; }
176 
178  double getPhi() const { return getPhi0(); }
179 
185  double getOmega() const { return m_tau[iOmega]; }
186 
191  double getZ0() const { return m_tau[iZ0]; }
192 
197  double getTanLambda() const { return m_tau[iTanLambda]; }
198 
200  double getCotTheta() const { return getTanLambda(); }
201 
206  std::vector<float> getTau() const { return std::vector<float>(m_tau, m_tau + c_NPars); }
207 
212  std::vector<float> getCov() const { return std::vector<float>(m_cov5, m_cov5 + c_NCovEntries); }
213 
218  TMatrixDSym getCovariance5() const;
219 
221  Helix getHelix() const
222  { return Helix(getD0(), getPhi0(), getOmega(), getZ0(), getTanLambda()); }
223 
225  Helix getHelix(float momentumScale) const
226  { return Helix(getD0(), getPhi0(), getOmega() / momentumScale, getZ0(), getTanLambda()); }
227 
230  {
231  return UncertainHelix(getD0(), getPhi0(), getOmega(), getZ0(),
233  }
234 
237 
240 
242  virtual std::string getInfoHTML() const override;
243 
245  private:
246 
247  //---------------------------------------------------------------------------------------------------------------------------
249  const unsigned int m_pdg;
250 
252  const Double32_t m_pValue;
253 
259  static const unsigned int c_NPars = 5;
260  static const unsigned int c_NCovEntries = 5 * 6 / 2;
266  static const unsigned int iD0 = 0;
267  static const unsigned int iPhi0 = 1;
268  static const unsigned int iOmega = 2;
269  static const unsigned int iZ0 = 3;
270  static const unsigned int iTanLambda = 4;
274  Double32_t m_tau[c_NPars];
275 
280  Double32_t m_cov5[c_NCovEntries];
281 
292 
294  static const uint16_t c_NDFFlag = 0xFFFF;
295 
297  uint16_t m_NDF;
298 
300  /* Version history:
301  ver 8: add NDF
302  ver 7: fixed sign errors in the translation of position and momentum covariances.
303  ver 6: use fixed size arrays instead of vectors (add schema evolution rule), use Double32_t.
304  ver 5: CDC Hit Pattern now a single variable (add schema evolution rule).
305  ver 4: added hit pattern
306  ver 3: back to vectors for the elements
307  ver 4: fixed-size arrays for the elements <------- incompatible with later version 4
308  ver 3: add data fields for hit patterns
309  ver 2: re-arranging of members
310  ver 1:
311  ver 2: <------- incompatible with later version 2
312  ver 1:
313  */
314  };
316 }
DataType y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:421
DataType z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:423
DataType x() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:419
The ParticleType class for identifying different particle types.
Definition: Const.h:289
This class represents an ideal helix in perigee parameterization.
Definition: Helix.h:59
TVector3 getPerigee() const
Getter for the perigee position.
Definition: Helix.cc:68
double getTransverseMomentum(const double bZ) const
Getter for the absolute value of the transverse momentum at the perigee.
Definition: Helix.cc:98
TVector3 getMomentum(const double bZ) const
Getter for vector of momentum at the perigee position.
Definition: Helix.cc:88
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.
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.
ROOT::Math::PxPyPzEVector get4Momentum() const
Getter for the 4Momentum at the closest approach of the track in the r/phi projection.
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.
Definition: ClusterUtils.h:23