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

#include <TrackAssociatorByChi2.h>

Inheritance diagram for TrackAssociatorByChi2:
TrackAssociatorBase

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 edm::Event *event=0, const edm::EventSetup *setup=0) const override
 Association Reco To Sim with Collections. More...
 
virtual reco::RecoToSimCollection associateRecoToSim (edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event=0, const edm::EventSetup *setup=0) 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 edm::Event *event=0, const edm::EventSetup *setup=0) const override
 Association Sim To Reco with Collections. More...
 
virtual reco::SimToRecoCollection associateSimToReco (edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event=0, const edm::EventSetup *setup=0) const override
 compare reco to sim the handle of reco::Track and TrackingParticle collections More...
 
 TrackAssociatorByChi2 (const edm::ESHandle< MagneticField > mF, const edm::ParameterSet &conf)
 Constructor with PSet. More...
 
 TrackAssociatorByChi2 (const edm::ESHandle< MagneticField > mF, double chi2Cut, bool onlyDiag, const edm::InputTag &beamspotSrc)
 Constructor with magnetic field, double, bool and InputTag. More...
 
 ~TrackAssociatorByChi2 ()
 Destructor. More...
 
- Public Member Functions inherited from TrackAssociatorBase
virtual
reco::RecoToSimCollectionSeed 
associateRecoToSim (edm::Handle< edm::View< TrajectorySeed > > &, edm::Handle< TrackingParticleCollection > &, const edm::Event *event, const edm::EventSetup *setup) const
 
virtual
reco::RecoToSimCollectionTCandidate 
associateRecoToSim (edm::Handle< TrackCandidateCollection > &, edm::Handle< TrackingParticleCollection > &, const edm::Event *event, const edm::EventSetup *setup) const
 
virtual
reco::SimToRecoCollectionSeed 
associateSimToReco (edm::Handle< edm::View< TrajectorySeed > > &, edm::Handle< TrackingParticleCollection > &, const edm::Event *event, const edm::EventSetup *setup) const
 
virtual
reco::SimToRecoCollectionTCandidate 
associateSimToReco (edm::Handle< TrackCandidateCollection > &, edm::Handle< TrackingParticleCollection > &, const edm::Event *event, const edm::EventSetup *setup) const
 
 TrackAssociatorBase ()
 Constructor. More...
 
virtual ~TrackAssociatorBase ()
 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

edm::InputTag bsSrc
 
double chi2cut
 
bool onlyDiagonal
 
edm::ESHandle< 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 28 of file TrackAssociatorByChi2.h.

Member Typedef Documentation

typedef std::map<double, SimTrack> TrackAssociatorByChi2::Chi2SimMap

Definition at line 31 of file TrackAssociatorByChi2.h.

Definition at line 32 of file TrackAssociatorByChi2.h.

Definition at line 33 of file TrackAssociatorByChi2.h.

Constructor & Destructor Documentation

TrackAssociatorByChi2::TrackAssociatorByChi2 ( const edm::ESHandle< MagneticField mF,
const edm::ParameterSet conf 
)
inline

Constructor with PSet.

Definition at line 36 of file TrackAssociatorByChi2.h.

References onlyDiagonal, and theMF.

36  :
37  chi2cut(conf.getParameter<double>("chi2cut")),
38  onlyDiagonal(conf.getParameter<bool>("onlyDiagonal")),
39  bsSrc(conf.getParameter<edm::InputTag>("beamSpot")) {
40  theMF=mF;
41  if (onlyDiagonal)
42  edm::LogInfo("TrackAssociator") << " ---- Using Off Diagonal Covariance Terms = 0 ---- " << "\n";
43  else
44  edm::LogInfo("TrackAssociator") << " ---- Using Off Diagonal Covariance Terms != 0 ---- " << "\n";
45  }
T getParameter(std::string const &) const
edm::ESHandle< MagneticField > theMF
TrackAssociatorByChi2::TrackAssociatorByChi2 ( const edm::ESHandle< MagneticField mF,
double  chi2Cut,
bool  onlyDiag,
const edm::InputTag beamspotSrc 
)
inline

Constructor with magnetic field, double, bool and InputTag.

Definition at line 48 of file TrackAssociatorByChi2.h.

References bsSrc, chi2cut, onlyDiagonal, and theMF.

48  {
49  chi2cut=chi2Cut;
50  onlyDiagonal=onlyDiag;
51  theMF=mF;
52  bsSrc = beamspotSrc;
53  }
edm::ESHandle< MagneticField > theMF
TrackAssociatorByChi2::~TrackAssociatorByChi2 ( )
inline

Destructor.

Definition at line 56 of file TrackAssociatorByChi2.h.

56 {}

Member Function Documentation

RecoToSimCollection TrackAssociatorByChi2::associateRecoToSim ( const edm::RefToBaseVector< reco::Track > &  tC,
const edm::RefVector< TrackingParticleCollection > &  tPCH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const
overridevirtual

Association Reco To Sim with Collections.

Implements TrackAssociatorBase.

Definition at line 121 of file TrackAssociatorByChi2.cc.

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

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

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

Reimplemented from TrackAssociatorBase.

Definition at line 74 of file TrackAssociatorByChi2.h.

References TrackAssociatorBase::associateRecoToSim(), event(), and HcalObjRepresent::setup().

77  {
78  return TrackAssociatorBase::associateRecoToSim(tCH,tPCH,event,setup);
79  }
virtual reco::RecoToSimCollection associateRecoToSim(edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event, const edm::EventSetup *setup) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
double TrackAssociatorByChi2::associateRecoToSim ( reco::TrackCollection::const_iterator  ,
TrackingParticleCollection::const_iterator  ,
const reco::BeamSpot  
) const
private

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

SimToRecoCollection TrackAssociatorByChi2::associateSimToReco ( const edm::RefToBaseVector< reco::Track > &  tC,
const edm::RefVector< TrackingParticleCollection > &  tPCH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const
overridevirtual

Association Sim To Reco with Collections.

Implements TrackAssociatorBase.

Definition at line 183 of file TrackAssociatorByChi2.cc.

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

186  {
187  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
188  e->getByLabel(bsSrc,recoBeamSpotHandle);
189  reco::BeamSpot bs = *recoBeamSpotHandle;
190 
191  SimToRecoCollection outputCollection;
192 
193  //dereference the edm::Ref only once
194  std::vector<const TrackingParticle*> tPC;
195  tPC.reserve(tPCH.size());
196  for(auto const& ref: tPCH) {
197  tPC.push_back(&(*ref));
198  }
199 
200  int tpindex =0;
201  for (auto tp=tPC.begin(); tp!=tPC.end(); tp++, ++tpindex){
202 
203  //skip tps with a very small pt
204  //if (sqrt((*tp)->momentum().perp2())<0.5) continue;
205  int charge = (*tp)->charge();
206  if (charge==0) continue;
207 
208  LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION===========" << "\n"
209  << "TrackingParticle #"<<tpindex<<" with pt=" << sqrt((*tp)->momentum().perp2()) << "\n"
210  << "===========================================" << "\n";
211 
212  Basic3DVector<double> momAtVtx((*tp)->momentum().x(),(*tp)->momentum().y(),(*tp)->momentum().z());
213  Basic3DVector<double> vert((*tp)->vertex().x(),(*tp)->vertex().y(),(*tp)->vertex().z());
214 
215  int tindex=0;
216  for (RefToBaseVector<reco::Track>::const_iterator rt=tC.begin(); rt!=tC.end(); rt++, tindex++){
217 
218  TrackBase::ParameterVector rParameters = (*rt)->parameters();
219  TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
220  if (onlyDiagonal) {
221  for (unsigned int i=0;i<5;i++){
222  for (unsigned int j=0;j<5;j++){
223  if (i!=j) recoTrackCovMatrix(i,j)=0;
224  }
225  }
226  }
227  recoTrackCovMatrix.Invert();
228 
229  double chi2 = getChi2(rParameters,recoTrackCovMatrix,momAtVtx,vert,charge,bs);
230 
231  if (chi2<chi2cut) {
232  //NOTE: tPCH and tPC have the same index for the same object
233  outputCollection.insert(tPCH[tpindex],
234  std::make_pair(tC[tindex],
235  -chi2));//-chi2 because the Association Map is ordered using std::greater
236  }
237  }
238  }
239  outputCollection.post_insert();
240  return outputCollection;
241 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
const_iterator end() const
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:73
void post_insert()
post insert action
T sqrt(T t)
Definition: SSEVec.h:48
int j
Definition: DBlmapReader.cc:9
void insert(const key_type &k, const data_type &v)
insert an association
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 begin() const
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:76
virtual reco::SimToRecoCollection TrackAssociatorByChi2::associateSimToReco ( edm::Handle< edm::View< reco::Track > > &  tCH,
edm::Handle< TrackingParticleCollection > &  tPCH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const
inlineoverridevirtual

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

Reimplemented from TrackAssociatorBase.

Definition at line 83 of file TrackAssociatorByChi2.h.

References TrackAssociatorBase::associateSimToReco(), event(), and HcalObjRepresent::setup().

86  {
87  return TrackAssociatorBase::associateSimToReco(tCH,tPCH,event,setup);
88  }
virtual reco::SimToRecoCollection associateSimToReco(edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event, const edm::EventSetup *setup) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
double TrackAssociatorByChi2::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

TrackAssociatorByChi2::RecoToSimPairAssociation TrackAssociatorByChi2::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 TrackAssociatorByChi2.cc.

References 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
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::pair< reco::Track, Chi2SimMap > RecoToSimPair
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:73
int j
Definition: DBlmapReader.cc:9
double f[11][100]
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:12
std::map< double, SimTrack > Chi2SimMap
std::vector< RecoToSimPair > RecoToSimPairAssociation
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:76
double TrackAssociatorByChi2::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 TrackAssociatorByChi2.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
edm::ESHandle< MagneticField > theMF
pair< bool, TrackBase::ParameterVector > TrackAssociatorByChi2::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 TrackAssociatorByChi2.cc.

References reco::trackingParametersAtClosestApproachToBeamSpot().

117  {
118  return reco::trackingParametersAtClosestApproachToBeamSpot(vertex,momAtVtx,charge, *theMF, bs);
119 }
edm::ESHandle< MagneticField > theMF
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

edm::InputTag TrackAssociatorByChi2::bsSrc
private

Definition at line 130 of file TrackAssociatorByChi2.h.

Referenced by TrackAssociatorByChi2().

double TrackAssociatorByChi2::chi2cut
private

Definition at line 128 of file TrackAssociatorByChi2.h.

Referenced by TrackAssociatorByChi2().

bool TrackAssociatorByChi2::onlyDiagonal
private

Definition at line 129 of file TrackAssociatorByChi2.h.

Referenced by TrackAssociatorByChi2().

edm::ESHandle<MagneticField> TrackAssociatorByChi2::theMF
private

Definition at line 127 of file TrackAssociatorByChi2.h.

Referenced by TrackAssociatorByChi2().