CMS 3D CMS Logo

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

#include <TrackAssociatorByPositionImpl.h>

Inheritance diagram for TrackAssociatorByPositionImpl:
reco::TrackToTrackingParticleAssociatorBaseImpl

Public Types

enum  Method { Method::chi2, Method::dist, Method::momdr, Method::posdr }
 
typedef std::vector< SimHitTPPairSimHitTPAssociationList
 
typedef std::pair< TrackingParticleRef, TrackPSimHitRefSimHitTPPair
 

Public Member Functions

reco::RecoToSimCollection associateRecoToSim (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &) 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
 compare reco to sim the handle of reco::Track and TrackingParticle collections More...
 
 TrackAssociatorByPositionImpl (const TrackingGeometry *geo, const Propagator *prop, const SimHitTPAssociationList *assocList, double qMinCut, double qCut, double positionMinimumDistance, Method method, bool minIfNoMatch, bool considerAllSimHits)
 
- 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
 compare reco to sim the handle of reco::Track and TrackingParticle collections More...
 
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
 compare reco to sim the handle of reco::Track and TrackingParticle collections More...
 
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

FreeTrajectoryState getState (const reco::Track &) const
 
TrajectoryStateOnSurface getState (const TrackingParticleRef &, const SimHitTPAssociationList &simHitsTPAssoc) const
 
double quality (const TrajectoryStateOnSurface &, const TrajectoryStateOnSurface &) const
 

Private Attributes

bool theConsiderAllSimHits
 
const TrackingGeometrytheGeometry
 
Method theMethod
 
bool theMinIfNoMatch
 
double thePositionMinimumDistance
 
const PropagatorthePropagator
 
double theQCut
 
double theQminCut
 
const SimHitTPAssociationListtheSimHitsTPAssoc
 

Detailed Description

Class that performs the association of reco::Tracks and TrackingParticles based on position in muon detector

Author
vlimant

Definition at line 30 of file TrackAssociatorByPositionImpl.h.

Member Typedef Documentation

Definition at line 34 of file TrackAssociatorByPositionImpl.h.

Definition at line 33 of file TrackAssociatorByPositionImpl.h.

Member Enumeration Documentation

Enumerator
chi2 
dist 
momdr 
posdr 

Definition at line 35 of file TrackAssociatorByPositionImpl.h.

35 { chi2, dist, momdr, posdr};

Constructor & Destructor Documentation

TrackAssociatorByPositionImpl::TrackAssociatorByPositionImpl ( const TrackingGeometry geo,
const Propagator prop,
const SimHitTPAssociationList assocList,
double  qMinCut,
double  qCut,
double  positionMinimumDistance,
Method  method,
bool  minIfNoMatch,
bool  considerAllSimHits 
)
inline

Definition at line 37 of file TrackAssociatorByPositionImpl.h.

References associateRecoToSim(), associateSimToReco(), and quality().

45  :
46  theGeometry(geo),
47  thePropagator(prop),
48  theSimHitsTPAssoc(assocList),
49  theQminCut(qMinCut),
50  theQCut(qCut),
53  theMinIfNoMatch(minIfNoMatch),
54  theConsiderAllSimHits(considerAllSimHits) {}
const SimHitTPAssociationList * theSimHitsTPAssoc

Member Function Documentation

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

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

Implements reco::TrackToTrackingParticleAssociatorBaseImpl.

Definition at line 113 of file TrackAssociatorByPositionImpl.cc.

References edm::AssociationMap< Tag >::insert(), TrajectoryStateOnSurface::isValid(), edm::AssociationMap< Tag >::post_insert(), jets_cff::quality, edm::RefToBaseVector< T >::size(), edm::RefVector< C, T, F >::size(), and TrajectoryStateOnSurface::surface().

Referenced by TrackAssociatorByPositionImpl().

114  {
115  RecoToSimCollection outputCollection;
116  //for each reco track find a matching tracking particle
117  std::pair<unsigned int,unsigned int> minPair;
118  const double dQmin_default=1542543;
119  double dQmin=dQmin_default;
120 
121  //cdj edm::Handle<SimHitTPAssociationProducer::SimHitTPAssociationList> simHitsTPAssoc;
122 
123  //warning: make sure the TP collection used in the map is the same used in the associator!
124  //e->getByLabel(_simHitTpMapTag,simHitsTPAssoc);
125 
126  for (unsigned int Ti=0; Ti!=tCH.size();++Ti){
127  //initial state (initial OR inner OR outter)
128  FreeTrajectoryState iState = getState(*(tCH)[Ti]);
129 
130  bool atLeastOne=false;
131  // for each tracking particle, find a state position and the plane to propagate the track to.
132  for (unsigned int TPi=0;TPi!=tPCH.size();++TPi) {
133  //get a state in the muon system
134  TrajectoryStateOnSurface simReferenceState = getState((tPCH)[TPi],*theSimHitsTPAssoc);
135  if (!simReferenceState.isValid()) continue;
136 
137  //propagate the TRACK to the surface
138  TrajectoryStateOnSurface trackReferenceState = thePropagator->propagate(iState,simReferenceState.surface());
139  if (!trackReferenceState.isValid()) continue;
140 
141  //comparison
142  double dQ= quality(trackReferenceState,simReferenceState);
143  if (dQ < theQCut){
144  atLeastOne=true;
145  outputCollection.insert(tCH[Ti],
146  std::make_pair(tPCH[TPi],-dQ));//association map with quality, is order greater-first
147  edm::LogVerbatim("TrackAssociatorByPositionImpl")<<"track number: "<<Ti
148  <<" associated with dQ: "<<dQ
149  <<" to TrackingParticle number: " <<TPi;}
150  if (dQ < dQmin){
151  dQmin=dQ;
152  minPair = std::make_pair(Ti,TPi);}
153  }//loop over tracking particles
154  if (theMinIfNoMatch && !atLeastOne && dQmin!=dQmin_default){
155  outputCollection.insert(tCH[minPair.first],
156  std::make_pair(tPCH[minPair.second],-dQmin));}
157  }//loop over tracks
158  outputCollection.post_insert();
159  return outputCollection;
160 }
double quality(const TrajectoryStateOnSurface &, const TrajectoryStateOnSurface &) const
FreeTrajectoryState getState(const reco::Track &) const
const SimHitTPAssociationList * theSimHitsTPAssoc
const SurfaceType & surface() const
void post_insert()
post insert action
size_type size() const
void insert(const key_type &k, const data_type &v)
insert an association
TrajectoryStateOnSurface propagate(STA const &state, SUR const &surface) const
Definition: Propagator.h:53
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
SimToRecoCollection TrackAssociatorByPositionImpl::associateSimToReco ( const edm::RefToBaseVector< reco::Track > &  tCH,
const edm::RefVector< TrackingParticleCollection > &  tPCH 
) const
overridevirtual

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

Implements reco::TrackToTrackingParticleAssociatorBaseImpl.

Definition at line 164 of file TrackAssociatorByPositionImpl.cc.

References edm::AssociationMap< Tag >::insert(), TrajectoryStateOnSurface::isValid(), edm::AssociationMap< Tag >::post_insert(), jets_cff::quality, edm::RefToBaseVector< T >::size(), edm::RefVector< C, T, F >::size(), and TrajectoryStateOnSurface::surface().

Referenced by TrackAssociatorByPositionImpl().

165  {
166  SimToRecoCollection outputCollection;
167  //for each tracking particle, find matching tracks.
168 
169  std::pair<unsigned int,unsigned int> minPair;
170  const double dQmin_default=1542543;
171  double dQmin=dQmin_default;
172 
173  //edm::Handle<SimHitTPAssociationProducer::SimHitTPAssociationList> simHitsTPAssoc;
174 
175  //warning: make sure the TP collection used in the map is the same used in the associator!
176  //e->getByLabel(_simHitTpMapTag,simHitsTPAssoc);
177 
178  for (unsigned int TPi=0;TPi!=tPCH.size();++TPi){
179  //get a state in the muon system
180  TrajectoryStateOnSurface simReferenceState= getState((tPCH)[TPi],*theSimHitsTPAssoc);
181 
182  if (!simReferenceState.isValid()) continue;
183  bool atLeastOne=false;
184  // propagate every track from any state (initial, inner, outter) to the surface
185  // and make the position test
186  for (unsigned int Ti=0; Ti!=tCH.size();++Ti){
187  //initial state
188  FreeTrajectoryState iState = getState(*(tCH)[Ti]);
189 
190  //propagation to surface
191  TrajectoryStateOnSurface trackReferenceState = thePropagator->propagate(iState,simReferenceState.surface());
192  if (!trackReferenceState.isValid()) continue;
193 
194  //comparison
195  double dQ= quality(trackReferenceState, simReferenceState);
196  if (dQ < theQCut){
197  atLeastOne=true;
198  outputCollection.insert(tPCH[TPi],
199  std::make_pair(tCH[Ti],-dQ));//association map with quality, is order greater-first
200  edm::LogVerbatim("TrackAssociatorByPositionImpl")<<"TrackingParticle number: "<<TPi
201  <<" associated with dQ: "<<dQ
202  <<" to track number: "<<Ti;}
203  if (dQ <dQmin){
204  dQmin=dQ;
205  minPair = std::make_pair(TPi,Ti);}
206  }//loop over tracks
207  if (theMinIfNoMatch && !atLeastOne && dQmin!=dQmin_default){
208  outputCollection.insert(tPCH[minPair.first],
209  std::make_pair(tCH[minPair.second],-dQmin));}
210  }//loop over tracking particles
211 
212  outputCollection.post_insert();
213  return outputCollection;
214 }
double quality(const TrajectoryStateOnSurface &, const TrajectoryStateOnSurface &) const
FreeTrajectoryState getState(const reco::Track &) const
const SimHitTPAssociationList * theSimHitsTPAssoc
const SurfaceType & surface() const
void post_insert()
post insert action
size_type size() const
void insert(const key_type &k, const data_type &v)
insert an association
TrajectoryStateOnSurface propagate(STA const &state, SUR const &surface) const
Definition: Propagator.h:53
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
FreeTrajectoryState TrackAssociatorByPositionImpl::getState ( const reco::Track track) const
private

Definition at line 72 of file TrackAssociatorByPositionImpl.cc.

References trajectoryStateTransform::initialFreeState().

72  {
73  //may be you want to do more than that if track does not go to IP
75 }
virtual const MagneticField * magneticField() const =0
FreeTrajectoryState initialFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
TrajectoryStateOnSurface TrackAssociatorByPositionImpl::getState ( const TrackingParticleRef st,
const SimHitTPAssociationList simHitsTPAssoc 
) const
private

Definition at line 14 of file TrackAssociatorByPositionImpl.cc.

References SurfaceSideDefinition::atCenterOfSurface, KineDebug3::count(), edmIntegrityCheck::d, createTree::dd, end, LogDebug, mag(), GeomDet::surface(), Surface::toGlobal(), and DetId::Tracker.

14  {
15 
16  using SimHitTPPair = std::pair<TrackingParticleRef, TrackPSimHitRef>;
17  SimHitTPPair clusterTPpairWithDummyTP(st,TrackPSimHitRef());//SimHit is dummy:
18  // sorting only the cluster is needed
19  auto range = std::equal_range(simHitsTPAssoc.begin(), simHitsTPAssoc.end(),
20  clusterTPpairWithDummyTP, [](const SimHitTPPair& iLHS, const SimHitTPPair& iRHS) -> bool {
21  return iLHS.first.key()> iRHS.first.key();
22  });
23 
24  // TrackingParticle* simtrack = const_cast<TrackingParticle*>(&st);
25  //loop over PSimHits
26  const PSimHit * psimhit=nullptr;
27  const BoundPlane * plane=nullptr;
28  double dLim=thePositionMinimumDistance;
29 
30  // look for the further most hit beyond a certain limit
31  auto start=range.first;
32  auto end=range.second;
33  LogDebug("TrackAssociatorByPositionImpl")<<range.second-range.first<<" PSimHits.";
34 
35  unsigned int count=0;
36  for (auto ip=start;ip!=end;++ip){
37 
38  TrackPSimHitRef psit = ip->second;
39 
40  //get the detid
41  DetId dd(psit->detUnitId());
42 
43  if (!theConsiderAllSimHits && dd.det()!=DetId::Tracker) continue;
44 
45  LogDebug("TrackAssociatorByPositionImpl")<<count++<<"] PSimHit on: "<<dd.rawId();
46  //get the surface from the global geometry
47  const GeomDet * gd=theGeometry->idToDet(dd);
48  if (!gd){edm::LogError("TrackAssociatorByPositionImpl")<<"no geomdet for: "<<dd.rawId()<<". will skip.";
49  continue;}
50  double d=gd->surface().toGlobal(psit->localPosition()).mag();
51  if (d>dLim ){
52  dLim=d;
53  psimhit=&(*psit);
54  plane=&gd->surface();}
55  }
56 
57 
58  if (psimhit && plane){
59  //build a trajectorystate on this surface
61  GlobalPoint initialPoint=plane->toGlobal(psimhit->localPosition());
62  GlobalVector initialMomentum=plane->toGlobal(psimhit->momentumAtEntry());
63  int initialCharge = (psimhit->particleType()>0) ? -1:1;
64  CartesianTrajectoryError initialCartesianErrors; //no error at initial state
65  const GlobalTrajectoryParameters initialParameters(initialPoint,initialMomentum,initialCharge,thePropagator->magneticField());
66  return TrajectoryStateOnSurface(initialParameters,initialCartesianErrors,*plane,surfaceside);}
67  else{
68  // edm::LogError("TrackAssociatorByPositionImpl")<<"no corresponding PSimHit for a tracking particle. will fail.";
69  return TrajectoryStateOnSurface();}
70 }
#define LogDebug(id)
Definition: start.py:1
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:106
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
virtual const GeomDet * idToDet(DetId) const =0
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
virtual const MagneticField * magneticField() const =0
std::pair< TrackingParticleRef, TrackPSimHitRef > SimHitTPPair
#define end
Definition: vmac.h:39
Definition: DetId.h:18
edm::Ref< edm::PSimHitContainer > TrackPSimHitRef
double TrackAssociatorByPositionImpl::quality ( const TrajectoryStateOnSurface tr,
const TrajectoryStateOnSurface sim 
) const
private

Definition at line 77 of file TrackAssociatorByPositionImpl.cc.

References vertices_cff::chi2, PV3DBase< T, PVType, FrameType >::eta(), TrajectoryStateOnSurface::globalDirection(), TrajectoryStateOnSurface::globalPosition(), TrajectoryStateOnSurface::localError(), TrajectoryStateOnSurface::localParameters(), funct::m, mag(), LocalTrajectoryError::matrix(), PV3DBase< T, PVType, FrameType >::phi(), findQualityFiles::v, and LocalTrajectoryParameters::vector().

Referenced by TrackAssociatorByPositionImpl().

77  {
78  switch(theMethod){
79  case Method::chi2:
80  {
83  int ierr = ! m.Invert();
84  if (ierr!=0) edm::LogInfo("TrackAssociatorByPositionImpl")<<"error inverting the error matrix:\n"<<m;
85  double est = ROOT::Math::Similarity(v,m);
86  return est;
87  break;
88  }
89  case Method::dist:
90  {
91  return (tr.globalPosition() - sim.globalPosition()).mag();
92  break;
93  }
94  case Method::momdr:
95  {
96  return (deltaR<double>(tr.globalDirection().eta(),tr.globalDirection().phi(),
97  sim.globalDirection().eta(),sim.globalDirection().phi()));
98  break;
99  }
100  case Method::posdr:
101  {
102  return (deltaR<double>(tr.globalPosition().eta(),tr.globalPosition().phi(),
103  sim.globalPosition().eta(),sim.globalPosition().phi()));
104  break;
105  }
106  }
107  //should never be reached
108  edm::LogError("TrackAssociatorByPositionImpl")<<"option: "<<static_cast<int>(theMethod)<<" has not been recognized. association has no meaning.";
109  return -1;
110 }
const LocalTrajectoryParameters & localParameters() const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
GlobalPoint globalPosition() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
AlgebraicVector5 vector() const
const AlgebraicSymMatrix55 & matrix() const
const LocalTrajectoryError & localError() const
ROOT::Math::SVector< double, 5 > AlgebraicVector5
T eta() const
Definition: PV3DBase.h:76
GlobalVector globalDirection() const

Member Data Documentation

bool TrackAssociatorByPositionImpl::theConsiderAllSimHits
private

Definition at line 79 of file TrackAssociatorByPositionImpl.h.

const TrackingGeometry* TrackAssociatorByPositionImpl::theGeometry
private

Definition at line 71 of file TrackAssociatorByPositionImpl.h.

Method TrackAssociatorByPositionImpl::theMethod
private

Definition at line 77 of file TrackAssociatorByPositionImpl.h.

bool TrackAssociatorByPositionImpl::theMinIfNoMatch
private

Definition at line 78 of file TrackAssociatorByPositionImpl.h.

double TrackAssociatorByPositionImpl::thePositionMinimumDistance
private

Definition at line 76 of file TrackAssociatorByPositionImpl.h.

const Propagator* TrackAssociatorByPositionImpl::thePropagator
private

Definition at line 72 of file TrackAssociatorByPositionImpl.h.

double TrackAssociatorByPositionImpl::theQCut
private

Definition at line 75 of file TrackAssociatorByPositionImpl.h.

double TrackAssociatorByPositionImpl::theQminCut
private

Definition at line 74 of file TrackAssociatorByPositionImpl.h.

const SimHitTPAssociationList* TrackAssociatorByPositionImpl::theSimHitsTPAssoc
private

Definition at line 73 of file TrackAssociatorByPositionImpl.h.