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 Basic3DVector<double> &momAtVtx,
12  const Basic3DVector<double> &vert,
13  int charge,
14  const MagneticField &magfield,
15  const reco::BeamSpot &bs) {
16  double chi2 = invalidChi2;
17 
18  std::pair<bool, reco::TrackBase::ParameterVector> params =
20  if (params.first) {
21  reco::TrackBase::ParameterVector sParameters = params.second;
22 
23  reco::TrackBase::ParameterVector diffParameters = rParameters - sParameters;
24  diffParameters[2] = reco::deltaPhi(diffParameters[2], 0.f);
25  chi2 = ROOT::Math::Dot(diffParameters * recoTrackCovMatrix, diffParameters);
26  chi2 /= 5;
27 
28  LogDebug("TrackAssociator") << "====NEW RECO TRACK WITH PT="
29  << sin(rParameters[1]) * float(charge) / rParameters[0] << "====\n"
30  << "qoverp sim: " << sParameters[0] << "\n"
31  << "lambda sim: " << sParameters[1] << "\n"
32  << "phi sim: " << sParameters[2] << "\n"
33  << "dxy sim: " << sParameters[3] << "\n"
34  << "dsz sim: " << sParameters[4] << "\n"
35  << ": " /*<< */ << "\n"
36  << "qoverp rec: " << rParameters[0] << "\n"
37  << "lambda rec: " << rParameters[1] << "\n"
38  << "phi rec: " << rParameters[2] << "\n"
39  << "dxy rec: " << rParameters[3] << "\n"
40  << "dsz rec: " << rParameters[4] << "\n"
41  << ": " /*<< */ << "\n"
42  << "chi2: " << chi2 << "\n";
43  }
44  return chi2;
45  }
46 
48  const reco::TrackBase::CovarianceMatrix &recoTrackCovMatrix,
49  const TrackingParticle &trackingParticle,
50  const MagneticField &magfield,
51  const reco::BeamSpot &bs) {
52  const int charge = trackingParticle.charge();
53  if (charge == 0)
54  return invalidChi2;
55 
56  const auto tpMom = trackingParticle.momentum();
57  Basic3DVector<double> momAtVtx(tpMom.x(), tpMom.y(), tpMom.z());
58  Basic3DVector<double> vert(trackingParticle.vertex());
59 
60  return trackAssociationChi2(rParameters, recoTrackCovMatrix, momAtVtx, vert, charge, magfield, bs);
61  }
62 
64  const TrackingParticle &trackingParticle,
65  const MagneticField &magfield,
66  const reco::BeamSpot &bs) {
67  auto cov = track.covariance();
68  cov.Invert();
69 
70  return trackAssociationChi2(track.parameters(), cov, trackingParticle, magfield, bs);
71  }
72 } // namespace track_associator
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11779
MessageLogger.h
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
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
track_associator::trackAssociationChi2
double trackAssociationChi2(const reco::TrackBase::ParameterVector &rParameters, const reco::TrackBase::CovarianceMatrix &recoTrackCovMatrix, const Basic3DVector< double > &momAtVtx, const Basic3DVector< double > &vert, int charge, const MagneticField &magfield, const reco::BeamSpot &bs)
basic method where chi2 is computed
Definition: trackAssociationChi2.cc:9
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:127
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:160
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:223
reco::TrackBase
Definition: TrackBase.h:62
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