Belle II Software  release-05-02-19
TrackFitResult.cc
1 /**************************************************************************
2 x * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2013, 2014, 2015 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Martin Heck, Tobias Schlüter *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 #include <mdst/dataobjects/TrackFitResult.h>
11 #include <framework/dataobjects/UncertainHelix.h>
12 #include <mdst/dataobjects/HitPatternCDC.h>
13 #include <mdst/dataobjects/HitPatternVXD.h>
14 
15 #include <framework/utilities/HTML.h>
16 
17 #include <boost/math/special_functions/gamma.hpp>
18 
19 #include <sstream>
20 #include <assert.h>
21 
22 using namespace Belle2;
23 
25  m_pdg(0), m_pValue(0),
26  m_hitPatternCDCInitializer(0),
27  m_hitPatternVXDInitializer(0),
28  m_NDF(c_NDFFlag)
29 {
30  memset(m_tau, 0, sizeof(m_tau));
31  memset(m_cov5, 0, sizeof(m_cov5));
32 }
33 
34 TrackFitResult::TrackFitResult(const TVector3& position, const TVector3& momentum,
35  const TMatrixDSym& covariance, const short int charge,
36  const Const::ParticleType& particleType, const float pValue,
37  const float bField,
38  const uint64_t hitPatternCDCInitializer,
39  const uint32_t hitPatternVXDInitializer,
40  const uint16_t NDF) :
41  m_pdg(std::abs(particleType.getPDGCode())), m_pValue(pValue),
42  m_hitPatternCDCInitializer(hitPatternCDCInitializer),
43  m_hitPatternVXDInitializer(hitPatternVXDInitializer),
44  m_NDF(NDF)
45 {
46  UncertainHelix h(position, momentum, charge, bField, covariance, pValue);
47 
48  m_tau[iD0] = h.getD0();
49  m_tau[iPhi0] = h.getPhi0();
50  m_tau[iOmega] = h.getOmega();
51  m_tau[iZ0] = h.getZ0();
52  m_tau[iTanLambda] = h.getTanLambda();
53 
54  // Upper half of the covariance matrix goes into m_cov5.
55  const TMatrixDSym& cov = h.getCovariance();
56  unsigned int counter = 0;
57  for (unsigned int i = 0; i < c_NPars; ++i) {
58  for (unsigned int j = i; j < c_NPars; ++j) {
59  m_cov5[counter++] = cov(i, j);
60  }
61  }
62  assert(counter == c_NCovEntries);
63 }
64 
65 TrackFitResult::TrackFitResult(const std::vector<float>& tau, const std::vector<float>& cov5,
66  const Const::ParticleType& particleType, const float pValue,
67  const uint64_t hitPatternCDCInitializer,
68  const uint32_t hitPatternVXDInitializer,
69  const uint16_t NDF) :
70  m_pdg(std::abs(particleType.getPDGCode())), m_pValue(pValue),
71  m_hitPatternCDCInitializer(hitPatternCDCInitializer),
72  m_hitPatternVXDInitializer(hitPatternVXDInitializer),
73  m_NDF(NDF)
74 {
75  if (tau.size() != c_NPars
76  || cov5.size() != c_NCovEntries)
77  B2FATAL("Invalid initializer for TrackFitResult.");
78 
79  for (unsigned int i = 0; i < c_NPars; ++i)
80  m_tau[i] = tau[i];
81  for (unsigned int i = 0; i < c_NCovEntries; ++i)
82  m_cov5[i] = cov5[i];
83 }
84 
86 {
87  if (m_NDF == c_NDFFlag) {
88  return -1;
89  }
90  return m_NDF;
91 }
92 
94 {
95  double pValue = getPValue();
96  int nDF = getNDF();
97  if (pValue == 0) {
98  return std::numeric_limits<double>::infinity();
99  }
100  if (nDF < 0) {
101  return std::numeric_limits<double>::quiet_NaN();
102  }
103  return 2 * boost::math::gamma_q_inv(nDF / 2., pValue);
104 }
105 
107 {
108  TMatrixDSym cov5(c_NPars);
109  unsigned int counter = 0;
110  for (unsigned int i = 0; i < c_NPars; ++i) {
111  for (unsigned int j = i; j < c_NPars; ++j) {
112  cov5(i, j) = cov5(j, i) = m_cov5[counter];
113  ++counter;
114  }
115  }
116  assert(counter == c_NCovEntries);
117  return cov5;
118 }
119 
121 {
123 }
124 
126 {
128 }
129 std::string TrackFitResult::getInfoHTML() const
130 {
131  std::stringstream out;
132  out << "<b>Fit Hypothesis (PDG)</b>: " << m_pdg << "<br>";
133  out << "<b>nPXDHits</b>: " << getHitPatternVXD().getNPXDHits() << "<br>";
134  out << "<b>nSVDHits</b>: " << getHitPatternVXD().getNSVDHits() << "<br>";
135  out << "<b>nCDCHits</b>: " << getHitPatternCDC().getNHits() << "<br>";
136  out << "<b>NDF</b>: " << getNDF() << "<br>";
137  out << " <br>";
138  out << "<b>d0</b>: " << m_tau[iD0] << " cm <br>";
139  out << "<b>phi0</b>: " << m_tau[iPhi0] << " rad <br>";
140  out << "<b>omega</b>: " << m_tau[iOmega] << " 1/GeV<br>";
141  out << "<b>z0</b>: " << m_tau[iZ0] << " cm <br>";
142  out << "<b>tanLambda</b>: " << m_tau[iTanLambda] << "<br>";
143  out << " <br>";
144  out << "<b>Covariance Matrix (d0, phi0, omega, z0, tanLambda)</b>: " << HTML::getString(getCovariance5()) << "<br>";
145  out << "<b>pValue</b>: " << m_pValue << "<br>";
146  return out.str();
147 }
Belle2::TrackFitResult::TrackFitResult
TrackFitResult()
Constructor initializing everything to zero.
Definition: TrackFitResult.cc:24
Belle2::TrackFitResult::getNDF
int getNDF() const
Getter for number of degrees of freedom of the track fit.
Definition: TrackFitResult.cc:85
Belle2::TrackFitResult::iTanLambda
static const unsigned int iTanLambda
Index tan lambda.
Definition: TrackFitResult.h:275
Belle2::UncertainHelix
This class represents an ideal helix in perigee parameterization including the covariance matrix of t...
Definition: UncertainHelix.h:40
Belle2::TrackFitResult::c_NDFFlag
static const uint16_t c_NDFFlag
backward compatibility initialisation for NDF
Definition: TrackFitResult.h:299
Belle2::HitPatternVXD::getNPXDHits
unsigned short getNPXDHits() const
Get total number of hits in the PXD.
Definition: HitPatternVXD.cc:147
Belle2::TrackFitResult::m_cov5
Double32_t m_cov5[c_NCovEntries]
The 15 = 5*(5+1)/2 covariance matrix elements.
Definition: TrackFitResult.h:285
Belle2::TrackFitResult::iD0
static const unsigned int iD0
Index for d0.
Definition: TrackFitResult.h:271
Belle2::TrackFitResult::getPValue
double getPValue() const
Getter for Chi2 Probability of the track fit.
Definition: TrackFitResult.h:163
Belle2::TrackFitResult::c_NCovEntries
static const unsigned int c_NCovEntries
Number of covariance entries.
Definition: TrackFitResult.h:265
Belle2::TrackFitResult::m_hitPatternVXDInitializer
const uint32_t m_hitPatternVXDInitializer
Member for initializing the information about hits in the VXD.
Definition: TrackFitResult.h:296
Belle2::TrackFitResult::m_pValue
const Double32_t m_pValue
Chi2 Probability of the fit.
Definition: TrackFitResult.h:257
Belle2::TrackFitResult::getHitPatternCDC
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
Definition: TrackFitResult.cc:120
Belle2::TrackFitResult::m_NDF
uint16_t m_NDF
Memeber for number of degrees of freedom.
Definition: TrackFitResult.h:302
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::TrackFitResult::iPhi0
static const unsigned int iPhi0
Index for phi0.
Definition: TrackFitResult.h:272
Belle2::HTML::getString
std::string getString(const TMatrixFBase &matrix, int precision=2, bool color=true)
get HTML table representing a matrix.
Definition: HTML.cc:18
Belle2::TrackFitResult::m_hitPatternCDCInitializer
uint64_t m_hitPatternCDCInitializer
Member for initializing the information about hits in the CDC.
Definition: TrackFitResult.h:291
Belle2::TrackFitResult::getChi2
double getChi2() const
Get chi2 given NDF and p-value.
Definition: TrackFitResult.cc:93
Belle2::TrackFitResult::getHitPatternVXD
HitPatternVXD getHitPatternVXD() const
Getter for the hit pattern in the VXD;.
Definition: TrackFitResult.cc:125
Belle2::TrackFitResult::m_pdg
const unsigned int m_pdg
PDG Code for hypothesis with which the corresponding fit was performed.
Definition: TrackFitResult.h:254
Belle2::TrackFitResult::getCovariance5
TMatrixDSym getCovariance5() const
Getter for covariance matrix of perigee parameters in matrix form.
Definition: TrackFitResult.cc:106
Belle2::Const::ParticleType
The ParticleType class for identifying different particle types.
Definition: Const.h:284
Belle2::TrackFitResult::iOmega
static const unsigned int iOmega
Index for omega.
Definition: TrackFitResult.h:273
Belle2::TrackFitResult::getInfoHTML
virtual std::string getInfoHTML() const override
Return a short summary of this object's contents in HTML format.
Definition: TrackFitResult.cc:129
Belle2::HitPatternCDC::getNHits
unsigned short getNHits() const
Get the total Number of CDC hits in the fit.
Definition: HitPatternCDC.cc:49
Belle2::HitPatternVXD::getNSVDHits
unsigned short getNSVDHits() const
Get total number of hits in the SVD.
Definition: HitPatternVXD.cc:83
Belle2::TrackFitResult::c_NPars
static const unsigned int c_NPars
Number of helix parameters.
Definition: TrackFitResult.h:264
Belle2::HitPatternVXD
Hit pattern of the VXD within a track.
Definition: HitPatternVXD.h:44
Belle2::HitPatternCDC
Hit pattern of CDC hits within a track.
Definition: HitPatternCDC.h:45
Belle2::TrackFitResult::m_tau
Double32_t m_tau[c_NPars]
perigee helix parameters; tau = d0, phi0, omega, z0, tanLambda.
Definition: TrackFitResult.h:279
Belle2::TrackFitResult::iZ0
static const unsigned int iZ0
Index for z0.
Definition: TrackFitResult.h:274