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
TrackAssociatorByChi2Impl Class Reference

#include <TrackAssociatorByChi2Impl.h>

Inheritance diagram for TrackAssociatorByChi2Impl:
reco::TrackToTrackingParticleAssociatorBaseImpl

Public Types

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

Public Member Functions

virtual reco::RecoToSimCollection associateRecoToSim (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &) const override
 Association Reco To Sim with Collections. More...
 
virtual 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...
 
virtual reco::SimToRecoCollection associateSimToReco (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &) const override
 Association Sim To Reco with Collections. More...
 
virtual 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::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::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 associateRecoToSim (reco::TrackCollection::const_iterator, TrackingParticleCollection::const_iterator, const reco::BeamSpot &) const
 compare reco::TrackCollection and TrackingParticleCollection iterators: returns the chi2 More...
 
double compareTracksParam (reco::TrackCollection::const_iterator, edm::SimTrackContainer::const_iterator, const math::XYZTLorentzVectorD &, const GlobalVector &, const reco::TrackBase::CovarianceMatrix &, const reco::BeamSpot &) const
 compare reco::TrackCollection and edm::SimTrackContainer iterators: returns the chi2 More...
 
RecoToSimPairAssociation compareTracksParam (const reco::TrackCollection &, const edm::SimTrackContainer &, const edm::SimVertexContainer &, const reco::BeamSpot &) const
 compare collections reco to sim 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...
 
std::pair< bool,
reco::TrackBase::ParameterVector
parametersAtClosestApproach (const Basic3DVector< double > &, const Basic3DVector< double > &, float, const reco::BeamSpot &) const
 propagate the track parameters of TrackinParticle from production vertex to the point of closest approach to the beam line. 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 35 of file TrackAssociatorByChi2Impl.h.

Member Typedef Documentation

Definition at line 38 of file TrackAssociatorByChi2Impl.h.

Definition at line 39 of file TrackAssociatorByChi2Impl.h.

Definition at line 40 of file TrackAssociatorByChi2Impl.h.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 57 of file TrackAssociatorByChi2Impl.h.

59  :
60  theMF(&mF),
61  theBeamSpot(&bs),
62  chi2cut(chi2Cut),
63  onlyDiagonal(onlyDiag) {
64  }
const reco::BeamSpot * theBeamSpot

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

References edm::RefToBaseVector< T >::begin(), beam_dqm_sourceclient-live_cfg::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().

122  {
123 
124  const reco::BeamSpot& bs = *theBeamSpot;
125 
126  RecoToSimCollection outputCollection;
127 
128  //dereference the edm::Refs only once
129  std::vector<TrackingParticle const*> tPC;
130  tPC.reserve(tPCH.size());
131  for(auto const& ref: tPCH) {
132  tPC.push_back(&(*ref));
133  }
134 
135  int tindex=0;
136  for (RefToBaseVector<reco::Track>::const_iterator rt=tC.begin(); rt!=tC.end(); rt++, tindex++){
137 
138  LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION===========" << "\n"
139  << "rec::Track #"<<tindex<<" with pt=" << (*rt)->pt() << "\n"
140  << "===========================================" << "\n";
141 
142  TrackBase::ParameterVector rParameters = (*rt)->parameters();
143 
144  TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
145  if (onlyDiagonal){
146  for (unsigned int i=0;i<5;i++){
147  for (unsigned int j=0;j<5;j++){
148  if (i!=j) recoTrackCovMatrix(i,j)=0;
149  }
150  }
151  }
152 
153  recoTrackCovMatrix.Invert();
154 
155  int tpindex =0;
156  for (auto tp=tPC.begin(); tp!=tPC.end(); tp++, ++tpindex){
157 
158  //skip tps with a very small pt
159  //if (sqrt((*tp)->momentum().perp2())<0.5) continue;
160  int charge = (*tp)->charge();
161  if (charge==0) continue;
162  Basic3DVector<double> momAtVtx((*tp)->momentum().x(),(*tp)->momentum().y(),(*tp)->momentum().z());
163  Basic3DVector<double> vert=(Basic3DVector<double>) (*tp)->vertex();
164 
165  double chi2 = getChi2(rParameters,recoTrackCovMatrix,momAtVtx,vert,charge,bs);
166 
167  if (chi2<chi2cut) {
168  outputCollection.insert(tC[tindex],
169  std::make_pair(tPCH[tpindex],
170  -chi2));//-chi2 because the Association Map is ordered using std::greater
171  }
172  }
173  }
174  outputCollection.post_insert();
175  return outputCollection;
176 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
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: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::RecoToSimCollection TrackAssociatorByChi2Impl::associateRecoToSim ( const edm::Handle< edm::View< reco::Track > > &  tCH,
const edm::Handle< TrackingParticleCollection > &  tPCH 
) const
inlineoverridevirtual

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

Reimplemented from reco::TrackToTrackingParticleAssociatorBaseImpl.

Definition at line 78 of file TrackAssociatorByChi2Impl.h.

79  {
80  return TrackToTrackingParticleAssociatorBaseImpl::associateRecoToSim(tCH,tPCH);
81  }
double TrackAssociatorByChi2Impl::associateRecoToSim ( reco::TrackCollection::const_iterator  ,
TrackingParticleCollection::const_iterator  ,
const reco::BeamSpot  
) const
private

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

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

References edm::RefToBaseVector< T >::begin(), edm::RefVector< C, T, F >::begin(), beam_dqm_sourceclient-live_cfg::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().

180  {
181  const reco::BeamSpot& bs = *theBeamSpot;
182 
183  SimToRecoCollection outputCollection;
184 
185  int tpindex =0;
186  for (auto tp=tPCH.begin(); tp!=tPCH.end(); tp++, ++tpindex){
187 
188  //skip tps with a very small pt
189  //if (sqrt(tp->momentum().perp2())<0.5) continue;
190  int charge = (*tp)->charge();
191  if (charge==0) continue;
192 
193  LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION===========" << "\n"
194  << "TrackingParticle #"<<tpindex<<" with pt=" << sqrt((*tp)->momentum().perp2()) << "\n"
195  << "===========================================" << "\n";
196 
197  Basic3DVector<double> momAtVtx((*tp)->momentum().x(),(*tp)->momentum().y(),(*tp)->momentum().z());
198  Basic3DVector<double> vert((*tp)->vertex().x(),(*tp)->vertex().y(),(*tp)->vertex().z());
199 
200  int tindex=0;
201  for (RefToBaseVector<reco::Track>::const_iterator rt=tC.begin(); rt!=tC.end(); rt++, tindex++){
202 
203  TrackBase::ParameterVector rParameters = (*rt)->parameters();
204  TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
205  if (onlyDiagonal) {
206  for (unsigned int i=0;i<5;i++){
207  for (unsigned int j=0;j<5;j++){
208  if (i!=j) recoTrackCovMatrix(i,j)=0;
209  }
210  }
211  }
212  recoTrackCovMatrix.Invert();
213 
214  double chi2 = getChi2(rParameters,recoTrackCovMatrix,momAtVtx,vert,charge,bs);
215 
216  if (chi2<chi2cut) {
217  outputCollection.insert(*tp,
218  std::make_pair(tC[tindex],
219  -chi2));//-chi2 because the Association Map is ordered using std::greater
220  }
221  }
222  }
223  outputCollection.post_insert();
224  return outputCollection;
225 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
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: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::SimToRecoCollection TrackAssociatorByChi2Impl::associateSimToReco ( const edm::Handle< edm::View< reco::Track > > &  tCH,
const edm::Handle< TrackingParticleCollection > &  tPCH 
) const
inlineoverridevirtual

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

Reimplemented from reco::TrackToTrackingParticleAssociatorBaseImpl.

Definition at line 85 of file TrackAssociatorByChi2Impl.h.

86  {
87  return TrackToTrackingParticleAssociatorBaseImpl::associateSimToReco(tCH,tPCH);
88  }
double TrackAssociatorByChi2Impl::compareTracksParam ( reco::TrackCollection::const_iterator  ,
edm::SimTrackContainer::const_iterator  ,
const math::XYZTLorentzVectorD ,
const GlobalVector ,
const reco::TrackBase::CovarianceMatrix ,
const reco::BeamSpot  
) const
private

compare reco::TrackCollection and edm::SimTrackContainer iterators: returns the chi2

TrackAssociatorByChi2Impl::RecoToSimPairAssociation TrackAssociatorByChi2Impl::compareTracksParam ( const reco::TrackCollection rtColl,
const edm::SimTrackContainer stColl,
const edm::SimVertexContainer svColl,
const reco::BeamSpot bs 
) const
private

compare collections reco to sim

Definition at line 40 of file TrackAssociatorByChi2Impl.cc.

References beam_dqm_sourceclient-live_cfg::chi2, reco::deltaPhi(), f, i, and j.

43  {
44 
45  RecoToSimPairAssociation outputVec;
46 
47  for (TrackCollection::const_iterator track=rtColl.begin(); track!=rtColl.end(); track++){
48  Chi2SimMap outMap;
49 
50  TrackBase::ParameterVector rParameters = track->parameters();
51 
52  TrackBase::CovarianceMatrix recoTrackCovMatrix = track->covariance();
53  if (onlyDiagonal){
54  for (unsigned int i=0;i<5;i++){
55  for (unsigned int j=0;j<5;j++){
56  if (i!=j) recoTrackCovMatrix(i,j)=0;
57  }
58  }
59  }
60  recoTrackCovMatrix.Invert();
61 
62  for (SimTrackContainer::const_iterator st=stColl.begin(); st!=stColl.end(); st++){
63 
64  Basic3DVector<double> momAtVtx(st->momentum().x(),st->momentum().y(),st->momentum().z());
65  Basic3DVector<double> vert = (Basic3DVector<double>) svColl[st->vertIndex()].position();
66 
67  std::pair<bool,reco::TrackBase::ParameterVector> params = parametersAtClosestApproach(vert, momAtVtx, st->charge(), bs);
68  if (params.first){
69  TrackBase::ParameterVector sParameters = params.second;
70 
71  TrackBase::ParameterVector diffParameters = rParameters - sParameters;
72  diffParameters[2] = reco::deltaPhi(diffParameters[2],0.f);
73  double chi2 = ROOT::Math::Dot(diffParameters * recoTrackCovMatrix, diffParameters);
74  chi2/=5;
75  if (chi2<chi2cut) outMap[chi2]=*st;
76  }
77  }
78  outputVec.push_back(RecoToSimPair(*track,outMap));
79  }
80  return outputVec;
81 }
int i
Definition: DBlmapReader.cc:9
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:74
std::vector< RecoToSimPair > RecoToSimPairAssociation
std::pair< reco::Track, Chi2SimMap > RecoToSimPair
int j
Definition: DBlmapReader.cc:9
double f[11][100]
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
std::pair< bool, reco::TrackBase::ParameterVector > parametersAtClosestApproach(const Basic3DVector< double > &, const Basic3DVector< double > &, float, const reco::BeamSpot &) const
propagate the track parameters of TrackinParticle from production vertex to the point of closest appr...
std::map< double, SimTrack > Chi2SimMap
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:77
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 83 of file TrackAssociatorByChi2Impl.cc.

References track_associator::getChi2().

88  {
89  return track_associator::getChi2(rParameters, recoTrackCovMatrix,momAtVtx, vert, charge, *theMF, bs);
90 }
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
pair< bool, TrackBase::ParameterVector > TrackAssociatorByChi2Impl::parametersAtClosestApproach ( const Basic3DVector< double > &  vertex,
const Basic3DVector< double > &  momAtVtx,
float  charge,
const reco::BeamSpot bs 
) const
private

propagate the track parameters of TrackinParticle from production vertex to the point of closest approach to the beam line.

Definition at line 114 of file TrackAssociatorByChi2Impl.cc.

References reco::trackingParametersAtClosestApproachToBeamSpot().

117  {
119 }
std::pair< bool, reco::TrackBase::ParameterVector > trackingParametersAtClosestApproachToBeamSpot(const Basic3DVector< double > &vertex, const Basic3DVector< double > &momAtVtx, float charge, const MagneticField &magField, const BeamSpot &bs)

Member Data Documentation

double TrackAssociatorByChi2Impl::chi2cut
private

Definition at line 128 of file TrackAssociatorByChi2Impl.h.

bool TrackAssociatorByChi2Impl::onlyDiagonal
private

Definition at line 129 of file TrackAssociatorByChi2Impl.h.

const reco::BeamSpot* TrackAssociatorByChi2Impl::theBeamSpot
private

Definition at line 127 of file TrackAssociatorByChi2Impl.h.

const MagneticField* TrackAssociatorByChi2Impl::theMF
private

Definition at line 126 of file TrackAssociatorByChi2Impl.h.