![]() |
![]() |
#include <TrackAssociatorByPosition.h>
Class that performs the association of reco::Tracks and TrackingParticles based on position in muon detector
Definition at line 32 of file TrackAssociatorByPosition.h.
TrackAssociatorByPosition::TrackAssociatorByPosition | ( | const edm::ParameterSet & | iConfig, |
const TrackingGeometry * | geo, | ||
const Propagator * | prop | ||
) | [inline] |
Constructor with propagator and PSet.
Definition at line 37 of file TrackAssociatorByPosition.h.
References _simHitTpMapTag, edm::ParameterSet::getParameter(), AlCaHLTBitMon_QueryRunRegistry::string, theConsiderAllSimHits, theGeometry, theMethod, theMinIfNoMatch, thePositionMinimumDistance, thePropagator, theQCut, and theQminCut.
{ theGeometry = geo; thePropagator = prop; theMinIfNoMatch = iConfig.getParameter<bool>("MinIfNoMatch"); theQminCut = iConfig.getParameter<double>("QminCut"); theQCut = iConfig.getParameter<double>("QCut"); thePositionMinimumDistance = iConfig.getParameter<double>("positionMinimumDistance"); std::string meth= iConfig.getParameter<std::string>("method"); if (meth=="chi2"){ theMethod =0; } else if (meth=="dist"){theMethod =1;} else if (meth=="momdr"){theMethod = 2;} else if (meth=="posdr"){theMethod = 3;} else{ edm::LogError("TrackAssociatorByPosition")<<meth<<" mothed not recognized. Use dr or chi2."; } theConsiderAllSimHits = iConfig.getParameter<bool>("ConsiderAllSimHits"); _simHitTpMapTag = iConfig.getParameter<edm::InputTag>("simHitTpMapTag"); };
TrackAssociatorByPosition::~TrackAssociatorByPosition | ( | ) | [inline] |
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 [virtual] |
compare reco to sim the handle of reco::Track and TrackingParticle collections
Implements TrackAssociatorBase.
Definition at line 109 of file TrackAssociatorByPosition.cc.
References edm::Event::getByLabel(), edm::AssociationMap< Tag >::insert(), TrajectoryStateOnSurface::isValid(), edm::AssociationMap< Tag >::post_insert(), edm::RefToBaseVector< T >::size(), edm::RefVector< C, T, F >::size(), and TrajectoryStateOnSurface::surface().
{ RecoToSimCollection outputCollection; //for each reco track find a matching tracking particle std::pair<unsigned int,unsigned int> minPair; const double dQmin_default=1542543; double dQmin=dQmin_default; //warning: make sure the TP collection used in the map is the same used in the associator! e->getByLabel(_simHitTpMapTag,simHitsTPAssoc); for (unsigned int Ti=0; Ti!=tCH.size();++Ti){ //initial state (initial OR inner OR outter) FreeTrajectoryState iState = getState(*(tCH)[Ti]); bool atLeastOne=false; // for each tracking particle, find a state position and the plane to propagate the track to. for (unsigned int TPi=0;TPi!=tPCH.size();++TPi) { //get a state in the muon system TrajectoryStateOnSurface simReferenceState = getState((tPCH)[TPi]); if (!simReferenceState.isValid()) continue; //propagate the TRACK to the surface TrajectoryStateOnSurface trackReferenceState = thePropagator->propagate(iState,simReferenceState.surface()); if (!trackReferenceState.isValid()) continue; //comparison double dQ= quality(trackReferenceState,simReferenceState); if (dQ < theQCut){ atLeastOne=true; outputCollection.insert(tCH[Ti], std::make_pair(edm::Ref<TrackingParticleCollection>(tPCH,TPi),-dQ));//association map with quality, is order greater-first edm::LogVerbatim("TrackAssociatorByPosition")<<"track number: "<<Ti <<" associated with dQ: "<<dQ <<" to TrackingParticle number: " <<TPi;} if (dQ < dQmin){ dQmin=dQ; minPair = std::make_pair(Ti,TPi);} }//loop over tracking particles if (theMinIfNoMatch && !atLeastOne && dQmin!=dQmin_default){ outputCollection.insert(tCH[minPair.first], std::make_pair(edm::Ref<TrackingParticleCollection>(tPCH,minPair.second),-dQmin));} }//loop over tracks outputCollection.post_insert(); return outputCollection; }
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 [virtual] |
compare reco to sim the handle of reco::Track and TrackingParticle collections
Implements TrackAssociatorBase.
Definition at line 160 of file TrackAssociatorByPosition.cc.
References edm::Event::getByLabel(), edm::AssociationMap< Tag >::insert(), TrajectoryStateOnSurface::isValid(), edm::AssociationMap< Tag >::post_insert(), edm::RefToBaseVector< T >::size(), edm::RefVector< C, T, F >::size(), and TrajectoryStateOnSurface::surface().
{ SimToRecoCollection outputCollection; //for each tracking particle, find matching tracks. std::pair<unsigned int,unsigned int> minPair; const double dQmin_default=1542543; double dQmin=dQmin_default; //warning: make sure the TP collection used in the map is the same used in the associator! e->getByLabel(_simHitTpMapTag,simHitsTPAssoc); for (unsigned int TPi=0;TPi!=tPCH.size();++TPi){ //get a state in the muon system TrajectoryStateOnSurface simReferenceState= getState((tPCH)[TPi]); if (!simReferenceState.isValid()) continue; bool atLeastOne=false; // propagate every track from any state (initial, inner, outter) to the surface // and make the position test for (unsigned int Ti=0; Ti!=tCH.size();++Ti){ //initial state FreeTrajectoryState iState = getState(*(tCH)[Ti]); //propagation to surface TrajectoryStateOnSurface trackReferenceState = thePropagator->propagate(iState,simReferenceState.surface()); if (!trackReferenceState.isValid()) continue; //comparison double dQ= quality(trackReferenceState, simReferenceState); if (dQ < theQCut){ atLeastOne=true; outputCollection.insert(edm::Ref<TrackingParticleCollection>(tPCH,TPi), std::make_pair(tCH[Ti],-dQ));//association map with quality, is order greater-first edm::LogVerbatim("TrackAssociatorByPosition")<<"TrackingParticle number: "<<TPi <<" associated with dQ: "<<dQ <<" to track number: "<<Ti;} if (dQ <dQmin){ dQmin=dQ; minPair = std::make_pair(TPi,Ti);} }//loop over tracks if (theMinIfNoMatch && !atLeastOne && dQmin!=dQmin_default){ outputCollection.insert(edm::Ref<TrackingParticleCollection>(tPCH,minPair.first), std::make_pair(tCH[minPair.second],-dQmin));} }//loop over tracking particles outputCollection.post_insert(); return outputCollection; }
FreeTrajectoryState TrackAssociatorByPosition::getState | ( | const reco::Track & | track | ) | const [private] |
Definition at line 68 of file TrackAssociatorByPosition.cc.
References trajectoryStateTransform::initialFreeState().
{ //may be you want to do more than that if track does not go to IP return trajectoryStateTransform::initialFreeState(track,thePropagator->magneticField()); }
TrajectoryStateOnSurface TrackAssociatorByPosition::getState | ( | const TrackingParticleRef | st | ) | const [private] |
Definition at line 13 of file TrackAssociatorByPosition.cc.
References SurfaceSideDefinition::atCenterOfSurface, prof2calltree::count, createTree::dd, end, LogDebug, mag(), SimHitTPAssociationProducer::simHitTPAssociationListGreater(), dqm_diff::start, GeomDet::surface(), Surface::toGlobal(), and align::Tracker.
{ std::pair<TrackingParticleRef, TrackPSimHitRef> clusterTPpairWithDummyTP(st,TrackPSimHitRef());//SimHit is dummy: for simHitTPAssociationListGreater // sorting only the cluster is needed auto range = std::equal_range(simHitsTPAssoc->begin(), simHitsTPAssoc->end(), clusterTPpairWithDummyTP, SimHitTPAssociationProducer::simHitTPAssociationListGreater); // TrackingParticle* simtrack = const_cast<TrackingParticle*>(&st); //loop over PSimHits const PSimHit * psimhit=0; const BoundPlane * plane=0; double dLim=thePositionMinimumDistance; // look for the further most hit beyond a certain limit auto start=range.first; auto end=range.second; LogDebug("TrackAssociatorByPosition")<<range.second-range.first<<" PSimHits."; unsigned int count=0; for (auto ip=start;ip!=end;++ip){ TrackPSimHitRef psit = ip->second; //get the detid DetId dd(psit->detUnitId()); if (!theConsiderAllSimHits && dd.det()!=DetId::Tracker) continue; LogDebug("TrackAssociatorByPosition")<<count++<<"] PSimHit on: "<<dd.rawId(); //get the surface from the global geometry const GeomDet * gd=theGeometry->idToDet(dd); if (!gd){edm::LogError("TrackAssociatorByPosition")<<"no geomdet for: "<<dd.rawId()<<". will skip."; continue;} double d=gd->surface().toGlobal(psit->localPosition()).mag(); if (d>dLim ){ dLim=d; psimhit=&(*psit); plane=&gd->surface();} } if (psimhit && plane){ //build a trajectorystate on this surface SurfaceSideDefinition::SurfaceSide surfaceside = SurfaceSideDefinition::atCenterOfSurface; GlobalPoint initialPoint=plane->toGlobal(psimhit->localPosition()); GlobalVector initialMomentum=plane->toGlobal(psimhit->momentumAtEntry()); int initialCharge = (psimhit->particleType()>0) ? -1:1; CartesianTrajectoryError initialCartesianErrors; //no error at initial state const GlobalTrajectoryParameters initialParameters(initialPoint,initialMomentum,initialCharge,thePropagator->magneticField()); return TrajectoryStateOnSurface(initialParameters,initialCartesianErrors,*plane,surfaceside);} else{ // edm::LogError("TrackAssociatorByPosition")<<"no corresponding PSimHit for a tracking particle. will fail."; return TrajectoryStateOnSurface();} }
double TrackAssociatorByPosition::quality | ( | const TrajectoryStateOnSurface & | tr, |
const TrajectoryStateOnSurface & | sim | ||
) | const |
Definition at line 73 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().
{ switch(theMethod){ case 0: { AlgebraicVector5 v(tr.localParameters().vector() - sim.localParameters().vector()); AlgebraicSymMatrix55 m(tr.localError().matrix()); int ierr = ! m.Invert(); if (ierr!=0) edm::LogInfo("TrackAssociatorByPosition")<<"error inverting the error matrix:\n"<<m; double est = ROOT::Math::Similarity(v,m); return est; break; } case 1: { return (tr.globalPosition() - sim.globalPosition()).mag(); break; } case 2: { return (deltaR<double>(tr.globalDirection().eta(),tr.globalDirection().phi(), sim.globalDirection().eta(),sim.globalDirection().phi())); break; } case 3: { return (deltaR<double>(tr.globalPosition().eta(),tr.globalPosition().phi(), sim.globalPosition().eta(),sim.globalPosition().phi())); break; } } //should never be reached edm::LogError("TrackAssociatorByPosition")<<"option: "<<theMethod<<" has not been recognized. association has no meaning."; return -1; }
Definition at line 92 of file TrackAssociatorByPosition.h.
Referenced by TrackAssociatorByPosition().
edm::Handle<SimHitTPAssociationProducer::SimHitTPAssociationList> TrackAssociatorByPosition::simHitsTPAssoc [mutable, private] |
Definition at line 91 of file TrackAssociatorByPosition.h.
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().