CMS 3D CMS Logo

trackAssociationChi2.cc
Go to the documentation of this file.
5 
6 namespace track_associator {
7  constexpr double invalidChi2 = 10000000000.;
8 
10  const reco::TrackBase::CovarianceMatrix &recoTrackCovMatrix,
11  const reco::TrackBase::ParameterVector &sParameters) {
12  double chi2 = invalidChi2;
13 
14  reco::TrackBase::ParameterVector diffParameters = rParameters - sParameters;
15  diffParameters[2] = reco::deltaPhi(diffParameters[2], 0.f);
16  chi2 = ROOT::Math::Dot(diffParameters * recoTrackCovMatrix, diffParameters);
17  chi2 /= 5;
18 
19  LogDebug("TrackAssociator") << "====NEW RECO TRACK WITH PT=" << sin(rParameters[1]) / rParameters[0] << "====\n"
20  << "qoverp sim: " << sParameters[0] << "\n"
21  << "lambda sim: " << sParameters[1] << "\n"
22  << "phi sim: " << sParameters[2] << "\n"
23  << "dxy sim: " << sParameters[3] << "\n"
24  << "dsz sim: " << sParameters[4] << "\n"
25  << ": " /*<< */ << "\n"
26  << "qoverp rec: " << rParameters[0] << "\n"
27  << "lambda rec: " << rParameters[1] << "\n"
28  << "phi rec: " << rParameters[2] << "\n"
29  << "dxy rec: " << rParameters[3] << "\n"
30  << "dsz rec: " << rParameters[4] << "\n"
31  << ": " /*<< */ << "\n"
32  << "chi2: " << chi2 << "\n";
33  return chi2;
34  }
35 
37  const reco::TrackBase::CovarianceMatrix &recoTrackCovMatrix,
38  const Basic3DVector<double> &momAtVtx,
39  const Basic3DVector<double> &vert,
40  int charge,
41  const MagneticField &magfield,
42  const reco::BeamSpot &bs) {
43  double chi2 = invalidChi2;
44 
45  std::pair<bool, reco::TrackBase::ParameterVector> params =
47  if (params.first) {
48  return trackAssociationChi2(rParameters, recoTrackCovMatrix, params.second);
49  }
50  return chi2;
51  }
52 
54  const reco::TrackBase::CovarianceMatrix &recoTrackCovMatrix,
55  const TrackingParticle &trackingParticle,
56  const MagneticField &magfield,
57  const reco::BeamSpot &bs) {
58  const int charge = trackingParticle.charge();
59  if (charge == 0)
60  return invalidChi2;
61 
62  const auto tpMom = trackingParticle.momentum();
63  Basic3DVector<double> momAtVtx(tpMom.x(), tpMom.y(), tpMom.z());
64  Basic3DVector<double> vert(trackingParticle.vertex());
65 
66  return trackAssociationChi2(rParameters, recoTrackCovMatrix, momAtVtx, vert, charge, magfield, bs);
67  }
68 
70  const TrackingParticle &trackingParticle,
71  const MagneticField &magfield,
72  const reco::BeamSpot &bs) {
73  auto cov = track.covariance();
74  cov.Invert();
75 
76  return trackAssociationChi2(track.parameters(), cov, trackingParticle, magfield, bs);
77  }
78 } // namespace track_associator
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11724
MessageLogger.h
f
double f[11][100]
Definition: MuScleFitUtils.cc:78
CalibrationSummaryClient_cfi.params
params
Definition: CalibrationSummaryClient_cfi.py:14
reco::deltaPhi
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
deltaPhi.h
hltPixelTracks_cff.chi2
chi2
Definition: hltPixelTracks_cff.py:25
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
TrackingParticle::charge
float charge() const
Electric charge. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:98
cms::cuda::bs
bs
Definition: HistoContainer.h:76
TrackingParticle
Monte Carlo truth information used for tracking validation.
Definition: TrackingParticle.h:29
reco::BeamSpot
Definition: BeamSpot.h:21
reco::TrackBase::ParameterVector
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:71
TrackingParticle::momentum
Vector momentum() const
spatial momentum vector
Definition: TrackingParticle.h:106
TrackingParticle::vertex
Point vertex() const
Parent vertex position.
Definition: TrackingParticle.h:169
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
reco::TrackBase
Definition: TrackBase.h:62
track_associator::trackAssociationChi2
double trackAssociationChi2(const reco::TrackBase::ParameterVector &rParameters, const reco::TrackBase::CovarianceMatrix &recoTrackCovMatrix, const reco::TrackBase::ParameterVector &sParameters)
basic method where chi2 is computed
Definition: trackAssociationChi2.cc:9
reco::trackingParametersAtClosestApproachToBeamSpot
std::pair< bool, reco::TrackBase::ParameterVector > trackingParametersAtClosestApproachToBeamSpot(const Basic3DVector< double > &vertex, const Basic3DVector< double > &momAtVtx, float charge, const MagneticField &magField, const BeamSpot &bs)
Definition: trackingParametersAtClosestApproachToBeamSpot.cc:21
track_associator
Definition: trackAssociationChi2.h:10
track_associator::invalidChi2
constexpr double invalidChi2
Definition: trackAssociationChi2.cc:7
trackAssociationChi2.h
reco::TrackBase::CovarianceMatrix
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:74
MagneticField
Definition: MagneticField.h:19
volumeBasedMagneticField_160812_cfi.magfield
magfield
Definition: volumeBasedMagneticField_160812_cfi.py:11
Basic3DVector
Definition: extBasic3DVector.h:30
trackingParametersAtClosestApproachToBeamSpot.h