CMS 3D CMS Logo

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