CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private 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 (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::Handle< 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::Handle< TrackingParticleCollection > &) const
 
 TrackToTrackingParticleAssociatorBaseImpl ()
 Constructor. More...
 
virtual ~TrackToTrackingParticleAssociatorBaseImpl ()
 Destructor. More...
 

Private Member Functions

double getChi2 (const reco::TrackBase::ParameterVector &rParameters, const reco::TrackBase::CovarianceMatrix &recoTrackCovMatrix, const Basic3DVector< double > &momAtVtx, const Basic3DVector< double > &vert, int charge, const reco::BeamSpot &) const
 basic method where chi2 is computed More...
 

Private Attributes

double chi2cut
 
bool onlyDiagonal
 
const reco::BeamSpottheBeamSpot
 
const MagneticFieldtheMF
 

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 34 of file TrackAssociatorByChi2Impl.h.

Member Typedef Documentation

Definition at line 36 of file TrackAssociatorByChi2Impl.h.

Definition at line 37 of file TrackAssociatorByChi2Impl.h.

Definition at line 38 of file TrackAssociatorByChi2Impl.h.

Constructor & Destructor Documentation

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

Member Function Documentation

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 21 of file TrackAssociatorByChi2Impl.cc.

References edm::RefToBaseVector< T >::begin(), ALCARECOTkAlJpsiMuMu_cff::charge, hltPixelTracks_cff::chi2, trackAssociatorByChi2_cfi::chi2cut, edm::RefToBaseVector< T >::end(), mps_fire::i, edm::AssociationMap< Tag >::insert(), dqmiolumiharvest::j, LogDebug, trackAssociatorByChi2_cfi::onlyDiagonal, edm::AssociationMap< Tag >::post_insert(), and edm::RefVector< C, T, F >::size().

22  {
23  const reco::BeamSpot& bs = *theBeamSpot;
24 
25  RecoToSimCollection outputCollection;
26 
27  //dereference the edm::Refs only once
28  std::vector<TrackingParticle const*> tPC;
29  tPC.reserve(tPCH.size());
30  for (auto const& ref : tPCH) {
31  tPC.push_back(&(*ref));
32  }
33 
34  int tindex = 0;
35  for (RefToBaseVector<reco::Track>::const_iterator rt = tC.begin(); rt != tC.end(); rt++, tindex++) {
36  LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION==========="
37  << "\n"
38  << "rec::Track #" << tindex << " with pt=" << (*rt)->pt() << "\n"
39  << "==========================================="
40  << "\n";
41 
42  TrackBase::ParameterVector rParameters = (*rt)->parameters();
43 
44  TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
45  if (onlyDiagonal) {
46  for (unsigned int i = 0; i < 5; i++) {
47  for (unsigned int j = 0; j < 5; j++) {
48  if (i != j)
49  recoTrackCovMatrix(i, j) = 0;
50  }
51  }
52  }
53 
54  recoTrackCovMatrix.Invert();
55 
56  int tpindex = 0;
57  for (auto tp = tPC.begin(); tp != tPC.end(); tp++, ++tpindex) {
58  //skip tps with a very small pt
59  //if (sqrt((*tp)->momentum().perp2())<0.5) continue;
60  int charge = (*tp)->charge();
61  if (charge == 0)
62  continue;
63  Basic3DVector<double> momAtVtx((*tp)->momentum().x(), (*tp)->momentum().y(), (*tp)->momentum().z());
64  Basic3DVector<double> vert = (Basic3DVector<double>)(*tp)->vertex();
65 
66  double chi2 = getChi2(rParameters, recoTrackCovMatrix, momAtVtx, vert, charge, bs);
67 
68  if (chi2 < chi2cut) {
69  outputCollection.insert(
70  tC[tindex],
71  std::make_pair(tPCH[tpindex],
72  -chi2)); //-chi2 because the Association Map is ordered using std::greater
73  }
74  }
75  }
76  outputCollection.post_insert();
77  return outputCollection;
78 }
#define LogDebug(id)
const reco::BeamSpot * theBeamSpot
double getChi2(const reco::TrackBase::ParameterVector &rParameters, const reco::TrackBase::CovarianceMatrix &recoTrackCovMatrix, const Basic3DVector< double > &momAtVtx, const Basic3DVector< double > &vert, int charge, const reco::BeamSpot &) const
basic method where chi2 is computed
const_iterator end() const
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:71
void post_insert()
post insert action
void insert(const key_type &k, const data_type &v)
insert an association
const_iterator begin() const
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:74
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 69 of file TrackAssociatorByChi2Impl.h.

References reco::TrackToTrackingParticleAssociatorBaseImpl::associateRecoToSim().

70  {
71  return TrackToTrackingParticleAssociatorBaseImpl::associateRecoToSim(tCH, tPCH);
72  }
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::RefToBaseVector< T >::begin(), edm::RefVector< C, T, F >::begin(), ALCARECOTkAlJpsiMuMu_cff::charge, hltPixelTracks_cff::chi2, trackAssociatorByChi2_cfi::chi2cut, edm::RefToBaseVector< T >::end(), edm::RefVector< C, T, F >::end(), mps_fire::i, edm::AssociationMap< Tag >::insert(), dqmiolumiharvest::j, LogDebug, trackAssociatorByChi2_cfi::onlyDiagonal, edm::AssociationMap< Tag >::post_insert(), and mathSSE::sqrt().

81  {
82  const reco::BeamSpot& bs = *theBeamSpot;
83 
84  SimToRecoCollection outputCollection;
85 
86  int tpindex = 0;
87  for (auto tp = tPCH.begin(); tp != tPCH.end(); tp++, ++tpindex) {
88  //skip tps with a very small pt
89  //if (sqrt(tp->momentum().perp2())<0.5) continue;
90  int charge = (*tp)->charge();
91  if (charge == 0)
92  continue;
93 
94  LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION==========="
95  << "\n"
96  << "TrackingParticle #" << tpindex << " with pt=" << sqrt((*tp)->momentum().perp2())
97  << "\n"
98  << "==========================================="
99  << "\n";
100 
101  Basic3DVector<double> momAtVtx((*tp)->momentum().x(), (*tp)->momentum().y(), (*tp)->momentum().z());
102  Basic3DVector<double> vert((*tp)->vertex().x(), (*tp)->vertex().y(), (*tp)->vertex().z());
103 
104  int tindex = 0;
105  for (RefToBaseVector<reco::Track>::const_iterator rt = tC.begin(); rt != tC.end(); rt++, tindex++) {
106  TrackBase::ParameterVector rParameters = (*rt)->parameters();
107  TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
108  if (onlyDiagonal) {
109  for (unsigned int i = 0; i < 5; i++) {
110  for (unsigned int j = 0; j < 5; j++) {
111  if (i != j)
112  recoTrackCovMatrix(i, j) = 0;
113  }
114  }
115  }
116  recoTrackCovMatrix.Invert();
117 
118  double chi2 = getChi2(rParameters, recoTrackCovMatrix, momAtVtx, vert, charge, bs);
119 
120  if (chi2 < chi2cut) {
121  outputCollection.insert(
122  *tp,
123  std::make_pair(tC[tindex],
124  -chi2)); //-chi2 because the Association Map is ordered using std::greater
125  }
126  }
127  }
128  outputCollection.post_insert();
129  return outputCollection;
130 }
#define LogDebug(id)
const reco::BeamSpot * theBeamSpot
double getChi2(const reco::TrackBase::ParameterVector &rParameters, const reco::TrackBase::CovarianceMatrix &recoTrackCovMatrix, const Basic3DVector< double > &momAtVtx, const Basic3DVector< double > &vert, int charge, const reco::BeamSpot &) const
basic method where chi2 is computed
const_iterator end() const
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:71
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
void post_insert()
post insert action
T sqrt(T t)
Definition: SSEVec.h:19
void insert(const key_type &k, const data_type &v)
insert an association
const_iterator begin() const
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:74
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 76 of file TrackAssociatorByChi2Impl.h.

References reco::TrackToTrackingParticleAssociatorBaseImpl::associateSimToReco(), and ALCARECOTkAlJpsiMuMu_cff::charge.

77  {
78  return TrackToTrackingParticleAssociatorBaseImpl::associateSimToReco(tCH, tPCH);
79  }
double TrackAssociatorByChi2Impl::getChi2 ( const reco::TrackBase::ParameterVector rParameters,
const reco::TrackBase::CovarianceMatrix recoTrackCovMatrix,
const Basic3DVector< double > &  momAtVtx,
const Basic3DVector< double > &  vert,
int  charge,
const reco::BeamSpot bs 
) const
private

basic method where chi2 is computed

Definition at line 12 of file TrackAssociatorByChi2Impl.cc.

References track_associator::trackAssociationChi2().

17  {
18  return track_associator::trackAssociationChi2(rParameters, recoTrackCovMatrix, momAtVtx, vert, charge, *theMF, bs);
19 }
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

Member Data Documentation

double TrackAssociatorByChi2Impl::chi2cut
private

Definition at line 92 of file TrackAssociatorByChi2Impl.h.

bool TrackAssociatorByChi2Impl::onlyDiagonal
private

Definition at line 93 of file TrackAssociatorByChi2Impl.h.

const reco::BeamSpot* TrackAssociatorByChi2Impl::theBeamSpot
private

Definition at line 91 of file TrackAssociatorByChi2Impl.h.

const MagneticField* TrackAssociatorByChi2Impl::theMF
private

Definition at line 90 of file TrackAssociatorByChi2Impl.h.