CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Attributes
TrackAssociatorByChi2Impl Class Reference

#include <TrackAssociatorByChi2Impl.h>

Inheritance diagram for TrackAssociatorByChi2Impl:
reco::TrackToTrackingParticleAssociatorBaseImpl

Public Types

typedef std::map< double, SimTrackChi2SimMap
 
typedef std::pair< reco::Track, Chi2SimMapRecoToSimPair
 
typedef std::vector< RecoToSimPairRecoToSimPairAssociation
 

Public Member Functions

reco::RecoToSimCollection associateRecoToSim (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &) const override
 Association Reco To Sim with Collections. More...
 
reco::RecoToSimCollection associateRecoToSim (const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const override
 compare reco to sim the handle of reco::Track and TrackingParticle collections More...
 
reco::SimToRecoCollection associateSimToReco (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &) const override
 Association Sim To Reco with Collections. More...
 
reco::SimToRecoCollection associateSimToReco (const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const override
 compare reco to sim the handle of reco::Track and TrackingParticle collections More...
 
 TrackAssociatorByChi2Impl (edm::EDProductGetter const &productGetter, const MagneticField &mF, const reco::BeamSpot &bs, double chi2Cut, bool onlyDiag)
 Constructor. More...
 
- Public Member Functions inherited from reco::TrackToTrackingParticleAssociatorBaseImpl
virtual reco::RecoToSimCollection associateRecoToSim (const edm::Handle< edm::View< reco::Track >> &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
 
virtual reco::RecoToSimCollectionSeed associateRecoToSim (const edm::Handle< edm::View< TrajectorySeed >> &, const edm::Handle< TrackingParticleCollection > &) const
 
virtual reco::RecoToSimCollectionTCandidate associateRecoToSim (const edm::Handle< TrackCandidateCollection > &, const edm::RefVector< TrackingParticleCollection > &) const
 
virtual reco::SimToRecoCollection associateSimToReco (const edm::Handle< edm::View< reco::Track >> &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
 
virtual reco::SimToRecoCollectionSeed associateSimToReco (const edm::Handle< edm::View< TrajectorySeed >> &, const edm::Handle< TrackingParticleCollection > &) const
 
virtual reco::SimToRecoCollectionTCandidate associateSimToReco (const edm::Handle< TrackCandidateCollection > &, const edm::RefVector< TrackingParticleCollection > &) const
 
 TrackToTrackingParticleAssociatorBaseImpl ()
 Constructor. More...
 
virtual ~TrackToTrackingParticleAssociatorBaseImpl ()
 Destructor. More...
 

Private Attributes

const reco::BeamSpotbeamSpot_
 
double chi2cut_
 
const MagneticFieldmF_
 
bool onlyDiagonal_
 
edm::EDProductGetter const * productGetter_
 

Detailed Description

Class that performs the association of reco::Tracks and TrackingParticles evaluating the chi2 of reco tracks parameters and sim tracks parameters. The cut can be tuned from the config file: see data/TrackAssociatorByChi2.cfi. Note that the Association Map is filled with -ch2 and not chi2 because it is ordered using std::greater: the track with the lowest association chi2 will be the first in the output map.It is possible to use only diagonal terms (associator by pulls) seeting onlyDiagonal = true in the PSet

Author
cerati, magni

Definition at line 38 of file TrackAssociatorByChi2Impl.h.

Member Typedef Documentation

◆ Chi2SimMap

Definition at line 40 of file TrackAssociatorByChi2Impl.h.

◆ RecoToSimPair

Definition at line 41 of file TrackAssociatorByChi2Impl.h.

◆ RecoToSimPairAssociation

Definition at line 42 of file TrackAssociatorByChi2Impl.h.

Constructor & Destructor Documentation

◆ TrackAssociatorByChi2Impl()

TrackAssociatorByChi2Impl::TrackAssociatorByChi2Impl ( edm::EDProductGetter const &  productGetter,
const MagneticField mF,
const reco::BeamSpot bs,
double  chi2Cut,
bool  onlyDiag 
)
inline

Constructor.

Definition at line 59 of file TrackAssociatorByChi2Impl.h.

Member Function Documentation

◆ associateRecoToSim() [1/2]

RecoToSimCollection TrackAssociatorByChi2Impl::associateRecoToSim ( const edm::RefToBaseVector< reco::Track > &  tC,
const edm::RefVector< TrackingParticleCollection > &  tPCH 
) const
overridevirtual

Association Reco To Sim with Collections.

Implements reco::TrackToTrackingParticleAssociatorBaseImpl.

Definition at line 14 of file TrackAssociatorByChi2Impl.cc.

References edm::RefToBaseVector< T >::begin(), cms::cuda::bs, ALCARECOTkAlJpsiMuMu_cff::charge, hltPixelTracks_cff::chi2, edm::RefToBaseVector< T >::end(), edm::first(), mps_fire::i, edm::AssociationMap< Tag >::insert(), dqmiolumiharvest::j, LogDebug, edm::AssociationMap< Tag >::post_insert(), hcal_runs::rt, edm::second(), edm::RefVector< C, T, F >::size(), cmsswSequenceInfo::tp, track_associator::trackAssociationChi2(), and reco::trackingParametersAtClosestApproachToBeamSpot().

15  {
16  const BeamSpot& bs = *beamSpot_;
17 
18  RecoToSimCollection outputCollection(productGetter_);
19 
20  //dereference the Refs only once and precompute params
21  std::vector<TrackingParticle const*> tPC;
22  std::vector<std::pair<bool, TrackBase::ParameterVector>> tpParams;
23  tPC.reserve(tPCH.size());
24  tpParams.reserve(tPCH.size());
25  for (auto const& ref : tPCH) {
26  auto const& tp = *ref;
27  tPC.push_back(&tp);
28 
29  int charge = tp.charge();
30  if (charge == 0)
31  tpParams.emplace_back(false, TrackBase::ParameterVector());
32  else {
33  using BVec = Basic3DVector<double>;
34  tpParams.emplace_back(
35  trackingParametersAtClosestApproachToBeamSpot(BVec(tp.vertex()), BVec(tp.momentum()), charge, *mF_, bs));
36  }
37  }
38 
39  int tindex = 0;
40  for (RefToBaseVector<Track>::const_iterator rt = tC.begin(); rt != tC.end(); rt++, tindex++) {
41  LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION==========="
42  << "\n"
43  << "rec::Track #" << tindex << " with pt=" << (*rt)->pt() << "\n"
44  << "==========================================="
45  << "\n";
46 
47  TrackBase::ParameterVector rParameters = (*rt)->parameters();
48 
49  TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
50  if (onlyDiagonal_) {
51  for (unsigned int i = 0; i < 5; i++) {
52  for (unsigned int j = 0; j < 5; j++) {
53  if (i != j)
54  recoTrackCovMatrix(i, j) = 0;
55  }
56  }
57  }
58 
59  recoTrackCovMatrix.Invert();
60 
61  int tpindex = 0;
62  for (auto tp = tPC.begin(); tp != tPC.end(); tp++, ++tpindex) {
63  //skip tps with a very small pt
64  //if (sqrt((*tp)->momentum().perp2())<0.5) continue;
65  if (!tpParams[tpindex].first)
66  continue;
67 
68  double chi2 = trackAssociationChi2(rParameters, recoTrackCovMatrix, tpParams[tpindex].second);
69 
70  if (chi2 < chi2cut_) {
71  //-chi2 because the Association Map is ordered using std::greater
72  outputCollection.insert(tC[tindex], std::make_pair(tPCH[tpindex], -chi2));
73  }
74  }
75  }
76  outputCollection.post_insert();
77  return outputCollection;
78 }
edm::EDProductGetter const * productGetter_
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:72
U second(std::pair< T, U > const &p)
const_iterator begin() const
const_iterator end() const
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
std::pair< bool, reco::TrackBase::ParameterVector > trackingParametersAtClosestApproachToBeamSpot(const Basic3DVector< double > &vertex, const Basic3DVector< double > &momAtVtx, float charge, const MagneticField &magField, const BeamSpot &bs)
double trackAssociationChi2(const reco::TrackBase::ParameterVector &rParameters, const reco::TrackBase::CovarianceMatrix &recoTrackCovMatrix, const reco::TrackBase::ParameterVector &sParameters)
basic method where chi2 is computed
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:75
#define LogDebug(id)
const reco::BeamSpot * beamSpot_

◆ associateRecoToSim() [2/2]

reco::RecoToSimCollection TrackAssociatorByChi2Impl::associateRecoToSim ( const edm::Handle< edm::View< reco::Track > > &  tCH,
const edm::Handle< TrackingParticleCollection > &  tPCH 
) const
inlineoverride

compare reco to sim the handle of reco::Track and TrackingParticle collections

Definition at line 77 of file TrackAssociatorByChi2Impl.h.

78  {
79  return TrackToTrackingParticleAssociatorBaseImpl::associateRecoToSim(tCH, tPCH);
80  }

◆ associateSimToReco() [1/2]

SimToRecoCollection TrackAssociatorByChi2Impl::associateSimToReco ( const edm::RefToBaseVector< reco::Track > &  tC,
const edm::RefVector< TrackingParticleCollection > &  tPCH 
) const
overridevirtual

Association Sim To Reco with Collections.

Implements reco::TrackToTrackingParticleAssociatorBaseImpl.

Definition at line 80 of file TrackAssociatorByChi2Impl.cc.

References edm::RefVector< C, T, F >::begin(), cms::cuda::bs, ALCARECOTkAlJpsiMuMu_cff::charge, hltPixelTracks_cff::chi2, edm::RefVector< C, T, F >::end(), mps_fire::i, edm::AssociationMap< Tag >::insert(), dqmiolumiharvest::j, LogDebug, edm::AssociationMap< Tag >::post_insert(), edm::RefToBaseVector< T >::size(), mathSSE::sqrt(), cmsswSequenceInfo::tp, track_associator::trackAssociationChi2(), and reco::trackingParametersAtClosestApproachToBeamSpot().

81  {
82  const BeamSpot& bs = *beamSpot_;
83 
84  SimToRecoCollection outputCollection(productGetter_);
85 
86  //compute track parameters only once
87  std::vector<TrackBase::ParameterVector> tPars;
88  tPars.reserve(tC.size());
89  std::vector<TrackBase::CovarianceMatrix> tCovs;
90  tCovs.reserve(tC.size());
91  for (auto const& ref : tC) {
92  auto const& aTk = *ref;
93  tPars.emplace_back(aTk.parameters());
94 
95  TrackBase::CovarianceMatrix recoTrackCovMatrix = aTk.covariance();
96  if (onlyDiagonal_) {
97  for (unsigned int i = 0; i < 5; i++) {
98  for (unsigned int j = 0; j < 5; j++) {
99  if (i != j)
100  recoTrackCovMatrix(i, j) = 0;
101  }
102  }
103  }
104  recoTrackCovMatrix.Invert();
105  tCovs.emplace_back(recoTrackCovMatrix);
106  }
107 
108  int tpindex = 0;
109  for (auto tp = tPCH.begin(); tp != tPCH.end(); tp++, ++tpindex) {
110  //skip tps with a very small pt
111  //if (sqrt(tp->momentum().perp2())<0.5) continue;
112  auto const& aTP = **tp;
113  int charge = aTP.charge();
114  if (charge == 0)
115  continue;
116 
117  LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION==========="
118  << "\n"
119  << "TrackingParticle #" << tpindex << " with pt=" << sqrt(aTP.momentum().perp2())
120  << "\n"
121  << "==========================================="
122  << "\n";
123 
124  using BVec = Basic3DVector<double>;
125  auto const tpBoolParams =
126  trackingParametersAtClosestApproachToBeamSpot(BVec(aTP.vertex()), BVec(aTP.momentum()), charge, *mF_, bs);
127  if (!tpBoolParams.first)
128  continue;
129 
130  for (unsigned int tindex = 0; tindex < tC.size(); tindex++) {
131  TrackBase::ParameterVector const& rParameters = tPars[tindex];
132  TrackBase::CovarianceMatrix const& recoTrackCovMatrix = tCovs[tindex];
133 
134  double chi2 = trackAssociationChi2(rParameters, recoTrackCovMatrix, tpBoolParams.second);
135 
136  if (chi2 < chi2cut_) {
137  //-chi2 because the Association Map is ordered using std::greater
138  outputCollection.insert(*tp, std::make_pair(tC[tindex], -chi2));
139  }
140  }
141  }
142  outputCollection.post_insert();
143  return outputCollection;
144 }
edm::EDProductGetter const * productGetter_
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:72
size_type size() const
T sqrt(T t)
Definition: SSEVec.h:19
std::pair< bool, reco::TrackBase::ParameterVector > trackingParametersAtClosestApproachToBeamSpot(const Basic3DVector< double > &vertex, const Basic3DVector< double > &momAtVtx, float charge, const MagneticField &magField, const BeamSpot &bs)
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223
double trackAssociationChi2(const reco::TrackBase::ParameterVector &rParameters, const reco::TrackBase::CovarianceMatrix &recoTrackCovMatrix, const reco::TrackBase::ParameterVector &sParameters)
basic method where chi2 is computed
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:75
#define LogDebug(id)
const reco::BeamSpot * beamSpot_

◆ associateSimToReco() [2/2]

reco::SimToRecoCollection TrackAssociatorByChi2Impl::associateSimToReco ( const edm::Handle< edm::View< reco::Track > > &  tCH,
const edm::Handle< TrackingParticleCollection > &  tPCH 
) const
inlineoverride

compare reco to sim the handle of reco::Track and TrackingParticle collections

Definition at line 84 of file TrackAssociatorByChi2Impl.h.

85  {
86  return TrackToTrackingParticleAssociatorBaseImpl::associateSimToReco(tCH, tPCH);
87  }

Member Data Documentation

◆ beamSpot_

const reco::BeamSpot* TrackAssociatorByChi2Impl::beamSpot_
private

Definition at line 92 of file TrackAssociatorByChi2Impl.h.

◆ chi2cut_

double TrackAssociatorByChi2Impl::chi2cut_
private

Definition at line 93 of file TrackAssociatorByChi2Impl.h.

◆ mF_

const MagneticField* TrackAssociatorByChi2Impl::mF_
private

Definition at line 91 of file TrackAssociatorByChi2Impl.h.

◆ onlyDiagonal_

bool TrackAssociatorByChi2Impl::onlyDiagonal_
private

Definition at line 94 of file TrackAssociatorByChi2Impl.h.

◆ productGetter_

edm::EDProductGetter const* TrackAssociatorByChi2Impl::productGetter_
private

Definition at line 90 of file TrackAssociatorByChi2Impl.h.