CMS 3D CMS Logo

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

#include <TrackAssociatorByPosition.h>

Inheritance diagram for TrackAssociatorByPosition:
TrackAssociatorBase

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
 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
 compare reco to sim the handle of reco::Track and TrackingParticle collections More...
 
double quality (const TrajectoryStateOnSurface &, const TrajectoryStateOnSurface &) const
 
 TrackAssociatorByPosition (const edm::ParameterSet &iConfig, const TrackingGeometry *geo, const Propagator *prop)
 Constructor with propagator and PSet. More...
 
 ~TrackAssociatorByPosition ()
 Destructor. More...
 
- Public Member Functions inherited from TrackAssociatorBase
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 More...
 
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::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 More...
 
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

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

Private Attributes

edm::InputTag _simHitTpMapTag
 
bool theConsiderAllSimHits
 
const TrackingGeometrytheGeometry
 
unsigned int theMethod
 
bool theMinIfNoMatch
 
double thePositionMinimumDistance
 
const PropagatorthePropagator
 
double theQCut
 
double theQminCut
 

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 TrackAssociatorByPosition.h.

Constructor & Destructor Documentation

TrackAssociatorByPosition::TrackAssociatorByPosition ( const edm::ParameterSet iConfig,
const TrackingGeometry geo,
const Propagator prop 
)
inline

Constructor with propagator and PSet.

Definition at line 35 of file TrackAssociatorByPosition.h.

References _simHitTpMapTag, edm::ParameterSet::getParameter(), AlCaHLTBitMon_QueryRunRegistry::string, theConsiderAllSimHits, theGeometry, theMethod, theMinIfNoMatch, thePositionMinimumDistance, thePropagator, theQCut, and theQminCut.

37  {
38  theGeometry = geo;
39  thePropagator = prop;
40  theMinIfNoMatch = iConfig.getParameter<bool>("MinIfNoMatch");
41  theQminCut = iConfig.getParameter<double>("QminCut");
42  theQCut = iConfig.getParameter<double>("QCut");
43  thePositionMinimumDistance = iConfig.getParameter<double>("positionMinimumDistance");
44  std::string meth= iConfig.getParameter<std::string>("method");
45  if (meth=="chi2"){ theMethod =0; }
46  else if (meth=="dist"){theMethod =1;}
47  else if (meth=="momdr"){theMethod = 2;}
48  else if (meth=="posdr"){theMethod = 3;}
49  else{
50  edm::LogError("TrackAssociatorByPosition")<<meth<<" mothed not recognized. Use dr or chi2."; }
51 
52  theConsiderAllSimHits = iConfig.getParameter<bool>("ConsiderAllSimHits");
53  _simHitTpMapTag = iConfig.getParameter<edm::InputTag>("simHitTpMapTag");
54  };
T getParameter(std::string const &) const
const TrackingGeometry * theGeometry
TrackAssociatorByPosition::~TrackAssociatorByPosition ( )
inline

Destructor.

Definition at line 58 of file TrackAssociatorByPosition.h.

58  {
59  };

Member Function Documentation

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

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

Implements TrackAssociatorBase.

Definition at line 110 of file TrackAssociatorByPosition.cc.

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

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

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

Implements TrackAssociatorBase.

Definition at line 163 of file TrackAssociatorByPosition.cc.

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

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

Definition at line 69 of file TrackAssociatorByPosition.cc.

References trajectoryStateTransform::initialFreeState().

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

Definition at line 14 of file TrackAssociatorByPosition.cc.

References SurfaceSideDefinition::atCenterOfSurface, prof2calltree::count, createTree::dd, end, LogDebug, mag(), SimHitTPAssociationProducer::simHitTPAssociationListGreater(), dqm_diff::start, GeomDet::surface(), Surface::toGlobal(), and DetId::Tracker.

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

Definition at line 74 of file TrackAssociatorByPosition.cc.

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

74  {
75  switch(theMethod){
76  case 0:
77  {
80  int ierr = ! m.Invert();
81  if (ierr!=0) edm::LogInfo("TrackAssociatorByPosition")<<"error inverting the error matrix:\n"<<m;
82  double est = ROOT::Math::Similarity(v,m);
83  return est;
84  break;
85  }
86  case 1:
87  {
88  return (tr.globalPosition() - sim.globalPosition()).mag();
89  break;
90  }
91  case 2:
92  {
93  return (deltaR<double>(tr.globalDirection().eta(),tr.globalDirection().phi(),
94  sim.globalDirection().eta(),sim.globalDirection().phi()));
95  break;
96  }
97  case 3:
98  {
99  return (deltaR<double>(tr.globalPosition().eta(),tr.globalPosition().phi(),
100  sim.globalPosition().eta(),sim.globalPosition().phi()));
101  break;
102  }
103  }
104  //should never be reached
105  edm::LogError("TrackAssociatorByPosition")<<"option: "<<theMethod<<" has not been recognized. association has no meaning.";
106  return -1;
107 }
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

edm::InputTag TrackAssociatorByPosition::_simHitTpMapTag
private

Definition at line 91 of file TrackAssociatorByPosition.h.

Referenced by TrackAssociatorByPosition().

bool TrackAssociatorByPosition::theConsiderAllSimHits
private

Definition at line 87 of file TrackAssociatorByPosition.h.

Referenced by TrackAssociatorByPosition().

const TrackingGeometry* TrackAssociatorByPosition::theGeometry
private

Definition at line 80 of file TrackAssociatorByPosition.h.

Referenced by TrackAssociatorByPosition().

unsigned int TrackAssociatorByPosition::theMethod
private

Definition at line 82 of file TrackAssociatorByPosition.h.

Referenced by TrackAssociatorByPosition().

bool TrackAssociatorByPosition::theMinIfNoMatch
private

Definition at line 85 of file TrackAssociatorByPosition.h.

Referenced by TrackAssociatorByPosition().

double TrackAssociatorByPosition::thePositionMinimumDistance
private

Definition at line 86 of file TrackAssociatorByPosition.h.

Referenced by TrackAssociatorByPosition().

const Propagator* TrackAssociatorByPosition::thePropagator
private

Definition at line 81 of file TrackAssociatorByPosition.h.

Referenced by TrackAssociatorByPosition().

double TrackAssociatorByPosition::theQCut
private

Definition at line 84 of file TrackAssociatorByPosition.h.

Referenced by TrackAssociatorByPosition().

double TrackAssociatorByPosition::theQminCut
private

Definition at line 83 of file TrackAssociatorByPosition.h.

Referenced by TrackAssociatorByPosition().