CMS 3D CMS Logo

TrackAssociatorByPosition.cc

Go to the documentation of this file.
00001 #include "SimTracker/TrackAssociation/interface/TrackAssociatorByPosition.h"
00002 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00003 
00004 #include <TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h>
00005 #include <Geometry/CommonDetUnit/interface/GeomDet.h>
00006 
00007 #include "DataFormats/Math/interface/deltaR.h"
00008 //#include "PhysicsTools/Utilities/interface/DeltaR.h"
00009 
00010 using namespace edm;
00011 using namespace reco;
00012 
00013 TrajectoryStateOnSurface TrackAssociatorByPosition::getState(const TrackingParticle & st)const{
00014   TrackingParticle* simtrack = const_cast<TrackingParticle*>(&st);
00015   //loop over PSimHits
00016   const PSimHit * psimhit=0;
00017   const BoundPlane * plane=0;
00018   double dLim=thePositionMinimumDistance;
00019 
00020   //    look for the further most hit beyond a certain limit
00021   LogDebug("TrackAssociatorByPosition")<<(int)(simtrack->trackerPSimHit_end()-simtrack->trackerPSimHit_begin())<<" PSimHits.";
00022   for (std::vector<PSimHit> ::const_iterator psit=simtrack->trackerPSimHit_begin();psit!=simtrack->trackerPSimHit_end();++psit){
00023     //get the detid
00024     DetId dd(psit->detUnitId());
00025     LogDebug("TrackAssociatorByPosition")<<psit-simtrack->trackerPSimHit_begin()
00026                                          <<"] PSimHit on: "<<dd.rawId();
00027     //get the surface from the global geometry
00028     const GeomDet * gd=theGeometry->idToDet(dd);
00029     if (!gd){edm::LogError("TrackAssociatorByPosition")<<"no geomdet for: "<<dd.rawId()<<". will fail.";
00030       return TrajectoryStateOnSurface();}
00031     double d=gd->surface().toGlobal(psit->localPosition()).mag();
00032     if (d>dLim ){
00033       dLim=d;
00034       psimhit=&(*psit);
00035       plane=&gd->surface();}
00036   }
00037 
00038 
00039   if (psimhit && plane){
00040     //build a trajectorystate on this surface    
00041     SurfaceSide surfaceside = atCenterOfSurface;
00042     GlobalPoint initialPoint=plane->toGlobal(psimhit->localPosition());
00043     GlobalVector initialMomentum=plane->toGlobal(psimhit->momentumAtEntry());
00044     int initialCharge =  (psimhit->particleType()>0) ? -1:1;
00045     CartesianTrajectoryError initialCartesianErrors(HepSymMatrix(6,0)); //no error at initial state  
00046     const GlobalTrajectoryParameters initialParameters(initialPoint,initialMomentum,initialCharge,thePropagator->magneticField());
00047     return TrajectoryStateOnSurface(initialParameters,initialCartesianErrors,*plane,surfaceside);}
00048   else{
00049     //    edm::LogError("TrackAssociatorByPosition")<<"no corresponding PSimHit for a tracking particle. will fail.";
00050     return TrajectoryStateOnSurface();}
00051 }
00052 
00053 FreeTrajectoryState TrackAssociatorByPosition::getState(const reco::Track & track)const{
00054   static TrajectoryStateTransform transformer;
00055   //may be you want to do more than that if track does not go to IP
00056   return transformer. initialFreeState(track,thePropagator->magneticField());
00057 }
00058 
00059 double TrackAssociatorByPosition::quality(const TrajectoryStateOnSurface & tr, const TrajectoryStateOnSurface & sim) const {
00060   switch(theMethod){
00061   case 0:
00062     {
00063       AlgebraicVector5 v(tr.localParameters().vector() - sim.localParameters().vector());
00064       AlgebraicSymMatrix55 m(tr.localError().matrix());
00065       int ierr = ! m.Invert();
00066       if (ierr!=0) edm::LogInfo("TrackAssociatorByPosition")<<"error inverting the error matrix:\n"<<m;
00067       double est = ROOT::Math::Similarity(v,m);
00068       return est;
00069       break;
00070     }
00071   case 1:
00072     {
00073       return (tr.globalPosition() - sim.globalPosition()).mag();
00074       break;
00075     }
00076   case 2:
00077     {
00078       return  (deltaR<double>(tr.globalDirection().eta(),tr.globalDirection().phi(),
00079                               sim.globalDirection().eta(),sim.globalDirection().phi()));
00080       break;
00081     }
00082   case 3:
00083     {
00084       return (deltaR<double>(tr.globalPosition().eta(),tr.globalPosition().phi(),
00085                              sim.globalPosition().eta(),sim.globalPosition().phi()));
00086       break;
00087     }
00088   }
00089   //should never be reached
00090   edm::LogError("TrackAssociatorByPosition")<<"option: "<<theMethod<<" has not been recognized. association has no meaning.";
00091   return -1;
00092 }
00093 
00094 
00095 RecoToSimCollection TrackAssociatorByPosition::associateRecoToSim(edm::RefToBaseVector<reco::Track>& tCH, 
00096                                                                   edm::RefVector<TrackingParticleCollection>& tPCH,
00097                                                                   const edm::Event * e ) const{
00098   RecoToSimCollection  outputCollection;
00099   //for each reco track find a matching tracking particle
00100   std::pair<uint,uint> minPair;
00101   const double dQmin_default=1542543;
00102   double dQmin=dQmin_default;
00103   for (uint Ti=0; Ti!=tCH.size();++Ti){
00104     //initial state (initial OR inner OR outter)
00105     FreeTrajectoryState iState = getState(*(tCH)[Ti]);
00106 
00107     bool atLeastOne=false;
00108     //    for each tracking particle, find a state position and the plane to propagate the track to.
00109     for (uint TPi=0;TPi!=tPCH.size();++TPi) {
00110       //get a state in the muon system 
00111       TrajectoryStateOnSurface simReferenceState = getState(*(tPCH)[TPi]);
00112       if (!simReferenceState.isValid()) continue;
00113 
00114       //propagate the TRACK to the surface
00115       TrajectoryStateOnSurface trackReferenceState = thePropagator->propagate(iState,simReferenceState.surface());
00116       if (!trackReferenceState.isValid()) continue; 
00117       
00118       //comparison
00119       double dQ= quality(trackReferenceState,simReferenceState);
00120       if (dQ < theQCut){
00121         atLeastOne=true;
00122         outputCollection.insert(tCH[Ti],
00123                                 std::make_pair(edm::Ref<TrackingParticleCollection>(tPCH,TPi),-dQ));//association map with quality, is order greater-first
00124         edm::LogVerbatim("TrackAssociatorByPosition")<<"track number: "<<Ti
00125                                                      <<" associated with dQ: "<<dQ
00126                                                      <<" to TrackingParticle number: " <<TPi;}
00127       if (dQ < dQmin){
00128         dQmin=dQ;
00129         minPair = std::make_pair(Ti,TPi);}
00130     }//loop over tracking particles
00131     if (theMinIfNoMatch && !atLeastOne && dQmin!=dQmin_default){
00132       outputCollection.insert(tCH[minPair.first],
00133                               std::make_pair(edm::Ref<TrackingParticleCollection>(tPCH,minPair.second),-dQmin));}
00134   }//loop over tracks
00135   outputCollection.post_insert();
00136   return outputCollection;
00137 }
00138 
00139 
00140 
00141 SimToRecoCollection TrackAssociatorByPosition::associateSimToReco(edm::RefToBaseVector<reco::Track>& tCH, 
00142                                                                   edm::RefVector<TrackingParticleCollection>& tPCH,
00143                                                                   const edm::Event * e ) const {
00144   SimToRecoCollection  outputCollection;
00145   //for each tracking particle, find matching tracks.
00146 
00147   std::pair<uint,uint> minPair;
00148   const double dQmin_default=1542543;
00149   double dQmin=dQmin_default;
00150   for (uint TPi=0;TPi!=tPCH.size();++TPi){
00151     //get a state in the muon system
00152     TrajectoryStateOnSurface simReferenceState= getState(*(tPCH)[TPi]);
00153       
00154     if (!simReferenceState.isValid()) continue; 
00155     bool atLeastOne=false;
00156     //  propagate every track from any state (initial, inner, outter) to the surface 
00157     //  and make the position test
00158     for (uint Ti=0; Ti!=tCH.size();++Ti){
00159       //initial state
00160       FreeTrajectoryState iState = getState(*(tCH)[Ti]);
00161         
00162       //propagation to surface
00163       TrajectoryStateOnSurface trackReferenceState = thePropagator->propagate(iState,simReferenceState.surface());
00164       if (!trackReferenceState.isValid()) continue;
00165         
00166       //comparison
00167       double dQ= quality(trackReferenceState, simReferenceState);
00168       if (dQ < theQCut){
00169         atLeastOne=true;
00170         outputCollection.insert(edm::Ref<TrackingParticleCollection>(tPCH,TPi),
00171                                 std::make_pair(tCH[Ti],-dQ));//association map with quality, is order greater-first
00172         edm::LogVerbatim("TrackAssociatorByPosition")<<"TrackingParticle number: "<<TPi
00173                                                      <<" associated with dQ: "<<dQ
00174                                                      <<" to track number: "<<Ti;}
00175       if (dQ <dQmin){
00176         dQmin=dQ;
00177         minPair = std::make_pair(TPi,Ti);}
00178     }//loop over tracks
00179     if (theMinIfNoMatch && !atLeastOne && dQmin!=dQmin_default){
00180       outputCollection.insert(edm::Ref<TrackingParticleCollection>(tPCH,minPair.first),
00181                               std::make_pair(tCH[minPair.second],-dQmin));}
00182   }//loop over tracking particles
00183   
00184   outputCollection.post_insert();
00185   return outputCollection;
00186 }

Generated on Tue Jun 9 17:47:55 2009 for CMSSW by  doxygen 1.5.4