CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes
TrackGenAssociatorByChi2Impl Class Reference

#include <TrackGenAssociatorByChi2Impl.h>

Inheritance diagram for TrackGenAssociatorByChi2Impl:
reco::TrackToGenParticleAssociatorBaseImpl

Public Types

typedef std::map< double,
SimTrack
Chi2SimMap
 
typedef std::pair< reco::Track,
Chi2SimMap
RecoToSimPair
 
typedef std::vector
< RecoToSimPair
RecoToSimPairAssociation
 

Public Member Functions

reco::GenToRecoCollection associateGenToReco (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< reco::GenParticleCollection > &) const override
 Association Sim To Reco with Collections (Gen Particle version) More...
 
virtual reco::GenToRecoCollection associateGenToReco (const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< reco::GenParticleCollection > &tPCH) const override
 compare reco to sim the handle of reco::Track and GenParticle collections More...
 
reco::RecoToGenCollection associateRecoToGen (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< reco::GenParticleCollection > &) const override
 Association Sim To Reco with Collections (Gen Particle version) More...
 
virtual reco::RecoToGenCollection associateRecoToGen (const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< reco::GenParticleCollection > &tPCH) const override
 compare reco to sim the handle of reco::Track and GenParticle collections More...
 
 TrackGenAssociatorByChi2Impl (const MagneticField &mF, const reco::BeamSpot &bs, double chi2Cut, bool onlyDiag)
 Constructor with PSet. More...
 
- Public Member Functions inherited from reco::TrackToGenParticleAssociatorBaseImpl
 TrackToGenParticleAssociatorBaseImpl ()
 Constructor. More...
 
virtual ~TrackToGenParticleAssociatorBaseImpl ()
 

Private Member Functions

double associateRecoToSim (reco::TrackCollection::const_iterator, TrackingParticleCollection::const_iterator, const reco::BeamSpot &) const
 compare reco::TrackCollection and TrackingParticleCollection iterators: returns the chi2 More...
 
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/TrackGenAssociatorByChi2.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 36 of file TrackGenAssociatorByChi2Impl.h.

Member Typedef Documentation

Definition at line 39 of file TrackGenAssociatorByChi2Impl.h.

Definition at line 40 of file TrackGenAssociatorByChi2Impl.h.

Definition at line 41 of file TrackGenAssociatorByChi2Impl.h.

Constructor & Destructor Documentation

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

Constructor with PSet.

Constructor with magnetic field, double, bool and InputTag

Definition at line 46 of file TrackGenAssociatorByChi2Impl.h.

48  :
49  theMF(&mF),
50  theBeamSpot(&bs),
51  chi2cut(chi2Cut),
52  onlyDiagonal(onlyDiag)
53  {
54  }

Member Function Documentation

GenToRecoCollection TrackGenAssociatorByChi2Impl::associateGenToReco ( const edm::RefToBaseVector< reco::Track > &  tC,
const edm::RefVector< reco::GenParticleCollection > &  tPCH 
) const
overridevirtual

Association Sim To Reco with Collections (Gen Particle version)

Implements reco::TrackToGenParticleAssociatorBaseImpl.

Definition at line 79 of file TrackGenAssociatorByChi2Impl.cc.

References edm::RefToBaseVector< T >::begin(), edm::RefVector< C, T, F >::begin(), RecoTauCleanerPlugins::charge, HLT_FULL_cff::chi2, edm::RefToBaseVector< T >::end(), edm::RefVector< C, T, F >::end(), track_associator::getChi2(), i, edm::AssociationMap< Tag >::insert(), j, LogDebug, edm::AssociationMap< Tag >::post_insert(), and mathSSE::sqrt().

Referenced by associateGenToReco().

80  {
81  reco::BeamSpot const& bs = *theBeamSpot;
82 
83  GenToRecoCollection outputCollection;
84 
85  int tpindex =0;
86  for (auto tp=tPCH.begin(); tp!=tPCH.end(); tp++, ++tpindex){
87 
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) continue;
92 
93  LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION===========" << "\n"
94  << "TrackingParticle #"<<tpindex<<" with pt=" << sqrt((*tp)->momentum().perp2()) << "\n"
95  << "===========================================" << "\n";
96 
97  Basic3DVector<double> momAtVtx((*tp)->momentum().x(),(*tp)->momentum().y(),(*tp)->momentum().z());
98  Basic3DVector<double> vert((*tp)->vertex().x(),(*tp)->vertex().y(),(*tp)->vertex().z());
99 
100  int tindex=0;
101  for (RefToBaseVector<reco::Track>::const_iterator rt=tC.begin(); rt!=tC.end(); rt++, tindex++){
102 
103  TrackBase::ParameterVector rParameters = (*rt)->parameters();
104  TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
105  if (onlyDiagonal) {
106  for (unsigned int i=0;i<5;i++){
107  for (unsigned int j=0;j<5;j++){
108  if (i!=j) recoTrackCovMatrix(i,j)=0;
109  }
110  }
111  }
112  recoTrackCovMatrix.Invert();
113 
114  double chi2 = getChi2(rParameters,recoTrackCovMatrix,momAtVtx,vert,charge,bs);
115 
116  if (chi2<chi2cut) {
117  outputCollection.insert(*tp,
118  std::make_pair(tC[tindex],
119  -chi2));//-chi2 because the Association Map is ordered using std::greater
120  }
121  }
122  }
123  outputCollection.post_insert();
124  return outputCollection;
125 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
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:74
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:255
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:250
void post_insert()
post insert action
T sqrt(T t)
Definition: SSEVec.h:18
int j
Definition: DBlmapReader.cc:9
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:77
virtual reco::GenToRecoCollection TrackGenAssociatorByChi2Impl::associateGenToReco ( const edm::Handle< edm::View< reco::Track > > &  tCH,
const edm::Handle< reco::GenParticleCollection > &  tPCH 
) const
inlineoverridevirtual

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

Implements reco::TrackToGenParticleAssociatorBaseImpl.

Definition at line 79 of file TrackGenAssociatorByChi2Impl.h.

References associateGenToReco(), edm::HandleBase::id(), j, edm::RefVector< C, T, F >::push_back(), and edm::RefToBaseVector< T >::push_back().

80  {
82  for (unsigned int j=0; j<tCH->size();j++)
83  tc.push_back(edm::RefToBase<reco::Track>(tCH,j));
84 
86  for (unsigned int j=0; j<tPCH->size();j++)
88 
89  return associateGenToReco(tc,tpc);
90  }
ProductID id() const
Definition: HandleBase.cc:15
int j
Definition: DBlmapReader.cc:9
reco::GenToRecoCollection associateGenToReco(const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< reco::GenParticleCollection > &) const override
Association Sim To Reco with Collections (Gen Particle version)
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:62
RecoToGenCollection TrackGenAssociatorByChi2Impl::associateRecoToGen ( const edm::RefToBaseVector< reco::Track > &  tC,
const edm::RefVector< reco::GenParticleCollection > &  tPCH 
) const
overridevirtual

Association Sim To Reco with Collections (Gen Particle version)

Implements reco::TrackToGenParticleAssociatorBaseImpl.

Definition at line 21 of file TrackGenAssociatorByChi2Impl.cc.

References edm::RefToBaseVector< T >::begin(), RecoTauCleanerPlugins::charge, HLT_FULL_cff::chi2, edm::RefToBaseVector< T >::end(), track_associator::getChi2(), i, edm::AssociationMap< Tag >::insert(), j, LogDebug, edm::AssociationMap< Tag >::post_insert(), and edm::RefVector< C, T, F >::size().

Referenced by associateRecoToGen().

22  {
23  reco::BeamSpot const& bs = *theBeamSpot;
24 
25  RecoToGenCollection outputCollection;
26 
27  //dereference the edm::Refs only once
28  std::vector<GenParticle 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 
37  LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION===========" << "\n"
38  << "rec::Track #"<<tindex<<" with pt=" << (*rt)->pt() << "\n"
39  << "===========================================" << "\n";
40 
41  TrackBase::ParameterVector rParameters = (*rt)->parameters();
42 
43  TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
44  if (onlyDiagonal){
45  for (unsigned int i=0;i<5;i++){
46  for (unsigned int j=0;j<5;j++){
47  if (i!=j) recoTrackCovMatrix(i,j)=0;
48  }
49  }
50  }
51 
52  recoTrackCovMatrix.Invert();
53 
54  int tpindex =0;
55  for (auto tp=tPC.begin(); tp!=tPC.end(); tp++, ++tpindex){
56 
57  //skip tps with a very small pt
58  //if (sqrt((*tp)->momentum().perp2())<0.5) continue;
59  int charge = (*tp)->charge();
60  if (charge==0) continue;
61  Basic3DVector<double> momAtVtx((*tp)->momentum().x(),(*tp)->momentum().y(),(*tp)->momentum().z());
62  Basic3DVector<double> vert=(Basic3DVector<double>) (*tp)->vertex();
63 
64  double chi2 = getChi2(rParameters,recoTrackCovMatrix,momAtVtx,vert,charge,bs);
65 
66  if (chi2<chi2cut) {
67  //NOTE: tPC and tPCH use the same index for the same object
68  outputCollection.insert(tC[tindex],
69  std::make_pair(tPCH[tpindex],
70  -chi2));//-chi2 because the Association Map is ordered using std::greater
71  }
72  }
73  }
74  outputCollection.post_insert();
75  return outputCollection;
76 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
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:74
void post_insert()
post insert action
int j
Definition: DBlmapReader.cc:9
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:99
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:77
virtual reco::RecoToGenCollection TrackGenAssociatorByChi2Impl::associateRecoToGen ( const edm::Handle< edm::View< reco::Track > > &  tCH,
const edm::Handle< reco::GenParticleCollection > &  tPCH 
) const
inlineoverridevirtual

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

Implements reco::TrackToGenParticleAssociatorBaseImpl.

Definition at line 65 of file TrackGenAssociatorByChi2Impl.h.

References associateRecoToGen(), edm::HandleBase::id(), j, edm::RefVector< C, T, F >::push_back(), and edm::RefToBaseVector< T >::push_back().

66  {
68  for (unsigned int j=0; j<tCH->size();j++)
69  tc.push_back(edm::RefToBase<reco::Track>(tCH,j));
70 
72  for (unsigned int j=0; j<tPCH->size();j++)
74 
75  return associateRecoToGen(tc,tpc);
76  }
ProductID id() const
Definition: HandleBase.cc:15
int j
Definition: DBlmapReader.cc:9
reco::RecoToGenCollection associateRecoToGen(const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< reco::GenParticleCollection > &) const override
Association Sim To Reco with Collections (Gen Particle version)
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:62
double TrackGenAssociatorByChi2Impl::associateRecoToSim ( reco::TrackCollection::const_iterator  ,
TrackingParticleCollection::const_iterator  ,
const reco::BeamSpot  
) const
private

compare reco::TrackCollection and TrackingParticleCollection iterators: returns the chi2

double TrackGenAssociatorByChi2Impl::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 TrackGenAssociatorByChi2Impl.cc.

References track_associator::getChi2().

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

Member Data Documentation

double TrackGenAssociatorByChi2Impl::chi2cut
private

Definition at line 111 of file TrackGenAssociatorByChi2Impl.h.

bool TrackGenAssociatorByChi2Impl::onlyDiagonal
private

Definition at line 112 of file TrackGenAssociatorByChi2Impl.h.

const reco::BeamSpot* TrackGenAssociatorByChi2Impl::theBeamSpot
private

Definition at line 110 of file TrackGenAssociatorByChi2Impl.h.

const MagneticField* TrackGenAssociatorByChi2Impl::theMF
private

Definition at line 109 of file TrackGenAssociatorByChi2Impl.h.