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 
17  double chi2 = invalidChi2;
18 
19  std::pair<bool,reco::TrackBase::ParameterVector> params = reco::trackingParametersAtClosestApproachToBeamSpot(vert, momAtVtx, charge, magfield, bs);
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=" << sin(rParameters[1])*float(charge)/rParameters[0] << "====\n"
29  << "qoverp sim: " << sParameters[0] << "\n"
30  << "lambda sim: " << sParameters[1] << "\n"
31  << "phi sim: " << sParameters[2] << "\n"
32  << "dxy sim: " << sParameters[3] << "\n"
33  << "dsz sim: " << sParameters[4] << "\n"
34  << ": " /*<< */ << "\n"
35  << "qoverp rec: " << rParameters[0] << "\n"
36  << "lambda rec: " << rParameters[1] << "\n"
37  << "phi rec: " << rParameters[2] << "\n"
38  << "dxy rec: " << rParameters[3] << "\n"
39  << "dsz rec: " << rParameters[4] << "\n"
40  << ": " /*<< */ << "\n"
41  << "chi2: " << chi2 << "\n";
42 
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 
63 
65  const TrackingParticle& trackingParticle,
66  const MagneticField& magfield,
67  const reco::BeamSpot& bs) {
68  auto cov = track.covariance();
69  cov.Invert();
70 
71  return trackAssociationChi2(track.parameters(), cov, trackingParticle, magfield, bs);
72  }
73 }
#define LogDebug(id)
Vector momentum() const
spatial momentum vector
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:74
#define constexpr
float charge() const
Electric charge. Note this is taken from the first SimTrack only.
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:731
double f[11][100]
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
std::pair< bool, reco::TrackBase::ParameterVector > trackingParametersAtClosestApproachToBeamSpot(const Basic3DVector< double > &vertex, const Basic3DVector< double > &momAtVtx, float charge, const MagneticField &magField, const BeamSpot &bs)
Point vertex() const
Parent vertex position.
ParameterVector parameters() const
Track parameters with one-to-one correspondence to the covariance matrix.
Definition: TrackBase.h:725
Monte Carlo truth information used for tracking validation.
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
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:77